atomes 1.1.15
atomes: an atomic scale modeling tool box
Loading...
Searching...
No Matches
read_xyz.c
Go to the documentation of this file.
1/* This file is part of the 'atomes' software
2
3'atomes' is free software: you can redistribute it and/or modify it under the terms
4of the GNU Affero General Public License as published by the Free Software Foundation,
5either version 3 of the License, or (at your option) any later version.
6
7'atomes' is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
8without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
9See the GNU General Public License for more details.
10
11You should have received a copy of the GNU Affero General Public License along with 'atomes'.
12If not, see <https://www.gnu.org/licenses/>
13
14Copyright (C) 2022-2024 by CNRS and University of Strasbourg */
15
22/*
23* This file: 'read_xyz.c'
24*
25* Contains:
26*
27
28 - The functions to read XYZ atomic coordinates
29
30*
31* List of functions:
32
33 int xyz_get_atom_coordinates ();
34 int open_xyz_file (int linec);
35
36*/
37
38#include "global.h"
39#include "glview.h"
40#include "callbacks.h"
41#include "interface.h"
42#include "project.h"
43#include "bind.h"
44#include "readers.h"
45#ifdef OPENMP
46# include <omp.h>
47#endif
48
55{
56 int i, j, k;
57 double v;
58 gchar * lia[4] = {"a", "b", "c", "d"};
59 this_reader -> nspec = 0;
60 active_project -> steps = this_reader -> steps;
61 active_project -> natomes = this_reader -> natomes;
63 this_reader -> z = allocdouble (1);
64 this_reader -> nsps = allocint (1);
65#ifdef OPENMP
66 int v_dummy;
67 int res;
68 int numth = omp_get_max_threads ();
69 gboolean doatoms = FALSE;
70 gchar * saved_line;
71 if (this_reader -> steps < numth)
72 {
73 if (numth >= 2*(this_reader -> steps-1))
74 {
75 doatoms = TRUE;
76 }
77 else
78 {
79 numth = this_reader -> steps;
80 }
81 }
82
83 if (doatoms)
84 {
85 // OpenMP on atoms
86 res = 0;
87 for (i=0; i<this_reader -> steps; i++)
88 {
89 k = i*(this_reader -> natomes + 2) + 2;
90 #pragma omp parallel for num_threads(numth) private(j,v,v_dummy,this_line,saved_line,this_word) shared(i,k,lia,coord_line,this_reader,active_project,res)
91 for (j=0; j<this_reader -> natomes; j++)
92 {
93 if (res == 2) goto enda;
94 this_line = g_strdup_printf ("%s", coord_line[k+j]);
95 saved_line = g_strdup_printf ("%s", this_line);
96 this_word = strtok_r (this_line, " ", & saved_line);
97 if (! this_word)
98 {
99 format_error (i+1, j+1, lia[0], k+j);
100 res = 2;
101 goto enda;
102 }
104 v_dummy = 0;
105 if (! v)
106 {
107 #pragma omp critical
108 v_dummy = set_v_dummy (this_word);
109 }
110 if (v || v_dummy)
111 {
112 if (! i)
113 {
114 v = v + v_dummy * 0.1;
115 #pragma omp critical
116 check_for_species (v, j);
117 }
118 this_word = strtok_r (NULL, " ", & saved_line);
119 if (! this_word)
120 {
121 format_error (i+1, j+1, lia[1], k+j);
122 res = 2;
123 goto enda;
124 }
125 active_project -> atoms[i][j].x = string_to_double ((gpointer)this_word);
126 this_word = strtok_r (NULL, " ", & saved_line);
127 if (! this_word)
128 {
129 format_error (i+1, j+1, lia[2], k+j);
130 res = 2;
131 goto enda;
132 }
133 active_project -> atoms[i][j].y = string_to_double ((gpointer)this_word);
134 this_word = strtok_r (NULL, " ", & saved_line);
135 if (! this_word)
136 {
137 format_error (i+1, j+1, lia[3], k+j);
138 res = 2;
139 goto enda;
140 }
141 active_project -> atoms[i][j].z = string_to_double ((gpointer)this_word);
142 }
143 else
144 {
145 format_error (i+1, j+1, lia[0], k+j);
146 res = 2;
147 goto enda;
148 }
149 g_free (this_line);
150 enda:;
151 }
152 if (res == 2) break;
153 }
154 }
155 else
156 {
157 res = 0;
158 #pragma omp parallel for num_threads(numth) private(i,j,k,v,v_dummy,this_line,saved_line,this_word) shared(lia,coord_line,this_reader,active_project,res)
159 for (i=0; i<this_reader -> steps; i++)
160 {
161 if (res == 2) goto ends;
162 k = i*(this_reader -> natomes + 2) + 2;
163 for (j=0; j<this_reader -> natomes; j++)
164 {
165 this_line = g_strdup_printf ("%s", coord_line[k+j]);
166 saved_line = g_strdup_printf ("%s", this_line);
167 this_word = strtok_r (this_line, " ", & saved_line);
168 if (! this_word)
169 {
170 format_error (i+1, j+1, lia[0], k+j);
171 res = 2;
172 goto ends;
173 }
175 v_dummy = 0;
176 if (! v)
177 {
178 #pragma omp critical
179 v_dummy = set_v_dummy (this_word);
180 }
181 if (v || v_dummy)
182 {
183 if (! i)
184 {
185 v = v + v_dummy * 0.1;
186 check_for_species (v, j);
187 }
188 this_word = strtok_r (NULL, " ", & saved_line);
189 if (! this_word)
190 {
191 format_error (i+1, j+1, lia[1], k+j);
192 res = 2;
193 goto ends;
194 }
195 active_project -> atoms[i][j].x = string_to_double ((gpointer)this_word);
196 this_word = strtok_r (NULL, " ", & saved_line);
197 if (! this_word)
198 {
199 format_error (i+1, j+1, lia[1], k+j);
200 res = 2;
201 goto ends;
202 }
203 active_project -> atoms[i][j].y = string_to_double ((gpointer)this_word);
204 this_word = strtok_r (NULL, " ", & saved_line);
205 if (! this_word)
206 {
207 format_error (i+1, j+1, lia[2], k+j);
208 res = 2;
209 goto ends;
210 }
211 active_project -> atoms[i][j].z = string_to_double ((gpointer)this_word);
212 }
213 else
214 {
215 format_error (i+1, j+1, lia[3], k+j);
216 res = 2;
217 goto ends;
218 }
219 g_free (this_line);
220 }
221 ends:;
222 }
223 }
224 g_free (coord_line);
225 if (res == 2) return 2;
226#else
227 line_node * tmp_line;
228 tail = head;
229 k = 0;
230 for (i=0; i<active_project -> steps; i++)
231 {
232 for (j=0; j<2; j++)
233 {
234 tmp_line = tail;
235 tail = tail -> next;
236 g_free (tmp_line);
237 k ++;
238 }
239 for (j=0; j<active_project -> natomes; j++)
240 {
241 this_line = g_strdup_printf ("%s", tail -> line);
242 this_word = strtok (this_line, " ");
243 if (! this_word)
244 {
245 format_error (i+1, j+1, lia[0], k+j);
246 return 2;
247 }
249 if (v)
250 {
251 if (! i) check_for_species (v, j);
252 this_word = strtok (NULL, " ");
253 if (! this_word)
254 {
255 format_error (i+1, j+1, lia[1], k+j);
256 return 2;
257 }
258 active_project -> atoms[i][j].x = string_to_double ((gpointer)this_word);
259 this_word = strtok (NULL, " ");
260 if (! this_word)
261 {
262 format_error (i+1, j+1, lia[2], k+j);
263 return 2;
264 }
265 active_project -> atoms[i][j].y = string_to_double ((gpointer)this_word);
266 this_word = strtok (NULL, " ");
267 if (! this_word)
268 {
269 format_error (i+1, j+1, lia[3], k+j);
270 return 2;
271 }
272 active_project -> atoms[i][j].z = string_to_double ((gpointer)this_word);
273 }
274 else
275 {
276 format_error (i+1, j+1, lia[0], k+j);
277 return 2;
278 }
279 tmp_line = tail;
280 tail = tail -> next;
281 g_free (tmp_line);
282 k ++;
283 }
284 }
285#endif
286 for (i=1; i<active_project -> steps; i++)
287 {
288 for (j=0; j<active_project -> natomes; j++)
289 {
290 active_project -> atoms[i][j].sp = active_project -> atoms[0][j].sp;
291 }
292 }
293 return 0;
294}
295
303int open_xyz_file (int linec)
304{
305 int res;
306#ifdef OPENMP
307 this_line = g_strdup_printf ("%s", coord_line[0]);
308 this_word = strtok (this_line, " ");
309 if (! this_word)
310 {
311 add_reader_info ("Wrong file format - cannot find the number of atoms !", 0);
312 add_reader_info ("Wrong file format - first line is corrupted !", 0);
313 res = 2;
314 goto end;
315 }
316 this_reader -> natomes = (int)string_to_double ((gpointer)this_word);
317 reader_info ("xyz", "Number of atoms", this_reader -> natomes);
318 g_free (this_line);
319 if (linec%(this_reader -> natomes + 2) != 0)
320 {
321 res = 2;
322 }
323 else
324 {
325 this_reader -> steps = linec / (this_reader -> natomes + 2);
326 reader_info ("xyz", "Number of steps", this_reader -> steps);
328 }
329#else
330 this_line = g_strdup_printf ("%s", head -> line);
331 this_word = strtok (this_line, " ");
332 if (! this_word)
333 {
334 add_reader_info ("Wrong file format - cannot find the number of atoms !", 0);
335 add_reader_info ("Wrong file format - the first line is corrupted !", 0);
336 res = 2;
337 goto end;
338 }
339 this_reader -> natomes = (int)string_to_double ((gpointer)this_word);
340 reader_info ("xyz", "Number of atoms", this_reader -> natomes);
341 g_free (this_line);
342 if (linec%(this_reader -> natomes + 2) != 0)
343 {
344 res = 2;
345 }
346 else
347 {
348 this_reader -> steps = linec / (this_reader -> natomes + 2);
349 reader_info ("xyz", "Number of steps", this_reader -> steps);
351 }
352#endif
353 end:
354 return res;
355}
Binding to the Fortran90 subroutines.
double get_z_from_periodic_table(gchar *lab)
get Z from atom label
Definition w_library.c:304
Callback declarations for main window.
void allocatoms(project *this_proj)
allocate project data
Definition open_p.c:160
int atoms[NUM_STYLES][2]
double * allocdouble(int val)
allocate a double * pointer
Definition global.c:459
int * allocint(int val)
allocate an int * pointer
Definition global.c:314
double string_to_double(gpointer string)
convert string to double
Definition global.c:624
Global variable declarations Global convenience function declarations Global data structure defin...
project * active_project
Definition project.c:47
Variable declarations related to the OpenGL window Function declarations related to the OpenGL wind...
Messaging function declarations.
double z
Definition ogl_draw.c:57
Function declarations for reading atomes project file Function declarations for saving atomes proje...
void check_for_species(double v, int ato)
Fill the species for each atom and the associated data.
Definition read_coord.c:196
line_node * head
Definition read_coord.c:75
gchar ** coord_line
Definition read_coord.c:72
coord_file * this_reader
Definition read_coord.c:71
void reader_info(gchar *type, gchar *sinf, int val)
display reader information
Definition read_coord.c:101
line_node * tail
Definition read_coord.c:76
void add_reader_info(gchar *info, int mid)
append information message to the reader information
Definition read_coord.c:86
char * this_word
Definition read_coord.c:74
void format_error(int stp, int ato, gchar *mot, int line)
Message to display an error message.
Definition read_coord.c:116
gchar * this_line
Definition read_coord.c:73
int set_v_dummy(gchar *this_word)
check if dummy is used for unknown species, if not then ask what to do
Definition read_coord.c:142
int xyz_get_atom_coordinates()
get the atomic coordinates from the XYZ file
Definition read_xyz.c:54
int open_xyz_file(int linec)
open XYZ file
Definition read_xyz.c:303
Functions declaration to read atomic coordinates.
GtkWidget * res[2]
Definition w_encode.c:212