atomes 1.2.1
atomes: an atomic scale modeling tool box
Loading...
Searching...
No Matches
read_c3d.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-2025 by CNRS and University of Strasbourg */
15
22/*
23* This file: 'read_c3d.c'
24*
25* Contains:
26*
27
28 - The functions to read Chem3d atomic coordinates
29
30*
31* List of functions:
32
33 int c3d_get_atom_coordinates ();
34 int open_c3d_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
49extern void check_for_species (double v, int ato);
50
57{
58 int i, j, k;
59 double v;
60 gchar * lia[5] = {"a", "b", "c", "d", "e"};
61 this_reader -> nspec = 0;
62 active_project -> steps = this_reader -> steps;
63 active_project -> natomes = this_reader -> natomes;
65 this_reader -> z = allocdouble (1);
66 this_reader -> nsps = allocint (1);
67#ifdef OPENMP
68 int v_dummy;
69 int res;
70 int numth = omp_get_max_threads ();
71 gboolean doatoms = FALSE;
72 gchar * saved_line;
73 if (this_reader -> steps < numth)
74 {
75 if (numth >= 2*(this_reader -> steps-1))
76 {
77 doatoms = TRUE;
78 }
79 else
80 {
81 numth = this_reader -> steps;
82 }
83 }
84
85 if (doatoms)
86 {
87 // OpenMP on atoms
88 res = 0;
89 for (i=0; i<this_reader -> steps; i++)
90 {
91 k = i*(this_reader -> natomes + 1) + 1;
92 #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)
93 for (j=0; j<this_reader -> natomes; j++)
94 {
95 if (res == 2) goto enda;
96 this_line = g_strdup_printf ("%s", coord_line[k+j]);
97 this_word = strtok_r (this_line, " ", & saved_line);
98 if (! this_word)
99 {
100 format_error (i+1, j+1, lia[0], k+j);
101 res = 2;
102 goto enda;
103 }
105 v_dummy = 0;
106 if (! v)
107 {
108 #pragma omp critical
109 v_dummy = set_v_dummy (this_word);
110 }
111 if (v || v_dummy)
112 {
113 if (! i)
114 {
115 v = v + v_dummy * 0.1;
116 #pragma omp critical
117 check_for_species (v, j);
118 }
119 this_word = strtok_r (NULL, " ", & saved_line);
120 if (! this_word)
121 {
122 format_error (i+1, j+1, lia[1], k+j);
123 res = 2;
124 goto enda;
125 }
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].x = 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].y = string_to_double ((gpointer)this_word);
142 this_word = strtok_r (NULL, " ", & saved_line);
143 if (! this_word)
144 {
145 format_error (i+1, j+1, lia[4], k+j);
146 res = 2;
147 goto enda;
148 }
149 active_project -> atoms[i][j].z = string_to_double ((gpointer)this_word);
150 }
151 else
152 {
153 format_error (i+1, j+1, lia[0], k+j);
154 res = 2;
155 goto enda;
156 }
157 g_free (this_line);
158 enda:;
159 }
160 if (res == 2) break;
161 }
162 }
163 else
164 {
165 res = 0;
166 #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)
167 for (i=0; i<this_reader -> steps; i++)
168 {
169 if (res == 2) goto ends;
170 k = i*(this_reader -> natomes + 1) + 1;
171 for (j=0; j<this_reader -> natomes; j++)
172 {
173 this_line = g_strdup_printf ("%s", coord_line[k+j]);
174 this_word = strtok_r (this_line, " ", & saved_line);
175 if (! this_word)
176 {
177 format_error (i+1, j+1, lia[0], k+j);
178 res = 2;
179 goto ends;
180 }
182 v_dummy = 0;
183 if (! v)
184 {
185 #pragma omp critical
186 v_dummy = set_v_dummy (this_word);
187 }
188 if (v || v_dummy)
189 {
190 if (! i)
191 {
192 v = v + v_dummy * 0.1;
193 check_for_species (v, j);
194 }
195 this_word = strtok_r (NULL, " ", & saved_line);
196 if (! this_word)
197 {
198 format_error (i+1, j+1, lia[1], k+j);
199 res = 2;
200 goto ends;
201 }
202 this_word = strtok_r (NULL, " ", & saved_line);
203 if (! this_word)
204 {
205 format_error (i+1, j+1, lia[2], k+j);
206 res = 2;
207 goto ends;
208 }
209 active_project -> atoms[i][j].x = string_to_double ((gpointer)this_word);
210 this_word = strtok_r (NULL, " ", & saved_line);
211 if (! this_word)
212 {
213 format_error (i+1, j+1, lia[3], k+j);
214 res = 2;
215 goto ends;
216 }
217 active_project -> atoms[i][j].y = string_to_double ((gpointer)this_word);
218 this_word = strtok_r (NULL, " ", & saved_line);
219 if (! this_word)
220 {
221 format_error (i+1, j+1, lia[4], k+j);
222 res = 2;
223 goto ends;
224 }
225 active_project -> atoms[i][j].z = string_to_double ((gpointer)this_word);
226 }
227 else
228 {
229 format_error (i+1, j+1, lia[0], k+j);
230 res = 2;
231 goto ends;
232 }
233 g_free (this_line);
234 }
235 ends:;
236 }
237 }
238 g_free (coord_line);
239 if (res == 2) return 2;
240#else
241 line_node * tmp_line;
242 tail = head;
243 k = 0;
244 for (i=0; i<active_project -> steps; i++)
245 {
246 tmp_line = tail;
247 tail = tail -> next;
248 g_free (tmp_line);
249 k ++;
250 for (j=0; j<active_project -> natomes; j++)
251 {
252 this_line = g_strdup_printf ("%s", tail -> line);
253 this_word = strtok (this_line, " ");
254 if (! this_word)
255 {
256 format_error (i+1, j+1, lia[0], k+j);
257 return 2;
258 }
260 if (v)
261 {
262 if (! i) check_for_species (v, j);
263 this_word = strtok (NULL, " ");
264 if (! this_word)
265 {
266 format_error (i+1, j+1, lia[1], k+j);
267 return 2;
268 }
269 this_word = strtok (NULL, " ");
270 if (! this_word)
271 {
272 format_error (i+1, j+1, lia[2], k+j);
273 return 2;
274 }
275 active_project -> atoms[i][j].x = string_to_double ((gpointer)this_word);
276 this_word = strtok (NULL, " ");
277 if (! this_word)
278 {
279 format_error (i+1, j+1, lia[3], k+j);
280 return 2;
281 }
282 active_project -> atoms[i][j].y = string_to_double ((gpointer)this_word);
283 this_word = strtok (NULL, " ");
284 if (! this_word)
285 {
286 format_error (i+1, j+1, lia[4], k+j);
287 return 2;
288 }
289 active_project -> atoms[i][j].z = string_to_double ((gpointer)this_word);
290 }
291 else
292 {
293 format_error (i+1, j+1, lia[0], k+j);
294 return 2;
295 }
296 tmp_line = tail;
297 tail = tail -> next;
298 g_free (tmp_line);
299 k ++;
300 }
301 }
302#endif
303 for (i=1; i<active_project -> steps; i++)
304 {
305 for (j=0; j<active_project -> natomes; j++)
306 {
307 active_project -> atoms[i][j].sp = active_project -> atoms[0][j].sp;
308 }
309 }
310 return 0;
311}
312
320int open_c3d_file (int linec)
321{
322 int res;
323#ifdef OPENMP
324 this_line = g_strdup_printf ("%s", coord_line[0]);
325 this_word = strtok (this_line, " ");
326 if (! this_word)
327 {
328 add_reader_info ("Wrong file format - cannot find the number of atoms !", 0);
329 add_reader_info ("Wrong file format - first line is corrupted !", 0);
330 res = 2;
331 goto end;
332 }
333 this_reader -> natomes = (int)string_to_double ((gpointer)this_word);
334 reader_info ("c3d", "Number of atoms", this_reader -> natomes);
335 g_free (this_line);
336 if (linec%(this_reader -> natomes + 1) != 0)
337 {
338 res = 2;
339 }
340 else
341 {
342 this_reader -> steps = linec / (this_reader -> natomes + 1);
343 reader_info ("c3d", "Number of steps", this_reader -> steps);
344 res = (this_reader -> steps > 1) ? 2 : c3d_get_atom_coordinates ();
345 }
346#else
347 this_line = g_strdup_printf ("%s", head -> line);
348 this_word = strtok (this_line, " ");
349 if (! this_word)
350 {
351 add_reader_info ("Wrong file format - cannot find the number of atoms !", 0);
352 add_reader_info ("Wrong file format - the first line is corrupted !", 0);
353 res = 2;
354 goto end;
355 }
356 this_reader -> natomes = (int)string_to_double ((gpointer)this_word);
357 reader_info ("chem3d", "Number of atoms", this_reader -> natomes);
358 g_free (this_line);
359 if (linec%(this_reader -> natomes + 1) != 0)
360 {
361 res = 2;
362 }
363 else
364 {
365 this_reader -> steps = linec / (this_reader -> natomes + 1);
366 reader_info ("chem3d", "Number of steps", this_reader -> steps);
367 res = (this_reader -> steps > 1) ? 2 : c3d_get_atom_coordinates ();
368 }
369#endif
370 end:
371 return res;
372}
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:163
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:61
Function declarations for reading atomes project file Function declarations for saving atomes proje...
int c3d_get_atom_coordinates()
get the atomic coordinates from the C3D file
Definition read_c3d.c:56
int open_c3d_file(int linec)
open C3D file
Definition read_c3d.c:320
void check_for_species(double v, int ato)
Fill the species for each atom and the associated data.
Definition read_coord.c:220
line_node * head
Definition read_coord.c:77
gchar ** coord_line
Definition read_coord.c:74
coord_file * this_reader
Definition read_coord.c:73
void reader_info(gchar *type, gchar *sinf, int val)
display reader information
Definition read_coord.c:125
line_node * tail
Definition read_coord.c:78
void add_reader_info(gchar *info, int mid)
append information message to the reader information
Definition read_coord.c:88
char * this_word
Definition read_coord.c:76
void format_error(int stp, int ato, gchar *mot, int line)
Message to display an error message.
Definition read_coord.c:140
gchar * this_line
Definition read_coord.c:75
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:166
Functions declaration to read atomic coordinates.
GtkWidget * res[2]
Definition w_encode.c:212