atomes 1.1.15
atomes: an atomic scale modeling tool box
Loading...
Searching...
No Matches
read_coord.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_coord.c'
24*
25* Contains:
26*
27
28 - The general functions to import atomic coordinates
29
30*
31* List of functions:
32
33 gboolean set_dummy_in_use (gchar * this_word);
34
35 int open_coord_file (gchar * filename, int fti);
36
37 void add_reader_info (gchar * info, int mid);
38 void reader_info (gchar * type, gchar * sinf, int val);
39 void format_error (int stp, int ato, gchar * mot, int line);
40 void check_for_species (double v, int ato);
41
42*/
43
44#include "global.h"
45#include "glview.h"
46#include "callbacks.h"
47#include "interface.h"
48#include "project.h"
49#include "bind.h"
50#include "cbuild_edit.h"
51#ifdef OPENMP
52# include <omp.h>
53#endif
54
55extern int open_xyz_file (int linec);
56extern int open_c3d_file (int linec);
57extern int open_pdb_file (int linec);
58extern int open_trj_file (int linec);
59extern int open_vas_file (int linec);
60extern int open_cif_file (int linec);
61extern int open_hist_file (int linec);
62extern void allocatoms (project * this_proj);
63extern chemical_data * alloc_chem_data (int spec);
64extern int build_crystal (gboolean visible, project * this_proj, gboolean to_wrap, gboolean show_clones, cell_info * cell, GtkWidget * widg);
65extern const gchar * dfi[2];
66
67extern atom_search * cif_search;
69
70FILE * coordf;
72gchar ** coord_line = NULL;
73gchar * this_line = NULL;
74char * this_word;
75line_node * head = NULL;
76line_node * tail = NULL;
77
86void add_reader_info (gchar * info, int mid)
87{
88 this_reader -> info = (this_reader -> info) ? g_strdup_printf ("%s\n%s", this_reader -> info, info) : g_strdup_printf ("%s", info);
89 if (! mid) this_reader -> mid = 0;
90}
91
101void reader_info (gchar * type, gchar * sinf, int val)
102{
103 g_print ("Reading coordinates [%s]: %s = %d\n", type, sinf, val);
104}
105
116void format_error (int stp, int ato, gchar * mot, int line)
117{
118 gchar * str;
119 if (ato < 0)
120 {
121 str = g_strdup_printf ("Wrong file format: error at step %d !\n"
122 "Wrong file format: record <b>%s</b> on line <b>%d</b> is corrupted !",
123 stp, mot, line);
124 }
125 else
126 {
127 str = g_strdup_printf ("Wrong file format: error at step %d, atom %d !\n"
128 "Wrong file format: record <b>%s</b> on line <b>%d</b> is corrupted !",
129 stp, ato, mot, line+1);
130 }
131 add_reader_info (str, 0);
132 g_free (str);
133}
134
143{
144 int i;
145 for (i=0; i<this_reader -> ndummy; i++)
146 {
147 if (g_strcmp0(this_reader -> dummy[i], this_word) == 0)
148 {
149 return i+1;
150 }
151 }
152 gchar ** dummy = NULL;
153 if (this_reader -> ndummy)
154 {
156 g_free (this_reader -> dummy);
157 }
158 this_reader -> ndummy ++;
159 this_reader -> dummy = g_malloc0(this_reader -> ndummy*sizeof*this_reader -> dummy);
160 if (dummy)
161 {
162 for (i=0; i<this_reader -> ndummy-1; i++)
163 {
164 this_reader -> dummy[i] = g_strdup_printf ("%s", dummy[i]);
165 }
166 }
167 this_reader -> dummy[this_reader -> ndummy-1] = g_strdup_printf ("%s", this_word);
168 // Dummy added, then do we use this dummy ?
169 gchar * str = g_strdup_printf ("Use dummy atom(s) for unknown %s species ?", this_word);
170 gboolean use_dummy = ask_yes_no ("Use dummy atom(s) ?", str, GTK_MESSAGE_QUESTION, MainWindow);
171 g_free (str);
172 if (use_dummy)
173 {
174 str = g_strdup_printf ("Using dummy atom(s) for unknown %s species", this_word);
175 add_reader_info (str, 1);
176 g_free (str);
177 return this_reader -> ndummy;
178 }
179 else
180 {
181 str = g_strdup_printf ("No dummy atom(s) for unknown %s species", this_word);
182 add_reader_info (str, 1);
183 g_free (str);
184 return 0;
185 }
186}
187
196void check_for_species (double v, int ato)
197{
198 int i;
199 gboolean add_spec;
200 add_spec = TRUE;
201 if (this_reader -> nspec)
202 {
203 for (i=0; i<this_reader -> nspec; i++)
204 {
205 if (this_reader -> z[i] == v)
206 {
207 if (this_reader -> cartesian)
208 {
209 active_project -> atoms[0][ato].sp = i;
210 }
211 else
212 {
213 this_reader -> lot[ato] = i;
214 }
215 this_reader -> nsps[i] ++;
216 add_spec = FALSE;
217 break;
218 }
219 }
220 }
221 if (add_spec)
222 {
223 if (this_reader -> nspec)
224 {
225 this_reader -> z = g_realloc (this_reader -> z, (this_reader -> nspec+1)*sizeof*this_reader -> z);
226 this_reader -> nsps = g_realloc (this_reader -> nsps, (this_reader -> nspec+1)*sizeof*this_reader -> nsps);
227 }
228 if (this_reader -> cartesian)
229 {
230 active_project -> atoms[0][ato].sp = this_reader -> nspec;
231 }
232 else
233 {
234 this_reader -> lot[ato] = this_reader -> nspec;
235 }
236 this_reader -> nsps[this_reader -> nspec] = 1;
237 this_reader -> z[this_reader -> nspec] = v;
238 this_reader -> nspec ++;
239 }
240}
241
250int open_coord_file (gchar * filename, int fti)
251{
252 int res = 0;
253#ifdef OPENMP
254 struct stat status;
255 res = stat (filename, & status);
256 if (res == -1)
257 {
258 add_reader_info ("Error - cannot get file statistics !", 0);
259 return 1;
260 }
261 int fsize = status.st_size;
262#endif
263 coordf = fopen (filename, dfi[0]);
264 if (! coordf)
265 {
266 add_reader_info ("Error - cannot open coordinates file !", 0);
267 return 1;
268 }
269 int i, j, k, l;
270#ifdef OPENMP
271 gchar * coord_content = g_malloc0(fsize*sizeof*coord_content);
272 fread (coord_content, fsize, 1, coordf);
273 fclose (coordf);
274 int linecount = 0;
275 for (j=0; j<fsize; j++) if (coord_content[j] == '\n') linecount ++;
276 coord_line = g_malloc0 (linecount*sizeof*coord_line);
277 coord_line[0] = & coord_content[0];
278 i = 1;
279 for (j=0; j<fsize; j++)
280 {
281 if (coord_content[j] == '\n')
282 {
283 coord_content[j] = '\0';
284 if (i < linecount)
285 {
286 coord_line[i] = & coord_content[j+1];
287 i ++;
288 }
289 }
290 }
291#else
292 gchar * buf = g_malloc0(LINE_SIZE*sizeof*buf);
293 head = NULL;
294 tail = NULL;
295 i = 0;
296 while (fgets(buf, LINE_SIZE, coordf))
297 {
298 if (head == NULL)
299 {
300 head = g_malloc0 (sizeof*head);
301 tail = g_malloc0 (sizeof*tail);
302 tail = head;
303 }
304 else
305 {
306 tail -> next = g_malloc0 (sizeof*tail -> next);
307 if (fti == 9 || fti == 10)
308 {
309 tail -> next -> prev = g_malloc0 (sizeof*tail -> next -> prev);
310 tail -> next -> prev = tail;
311 }
312 tail = tail -> next;
313 }
314 tail -> line = g_strdup_printf ("%s", buf);
315 tail -> line = substitute_string (tail -> line, "\n", "\0");
316 i ++;
317 }
318 g_free (buf);
319 fclose (coordf);
320#endif
321 if (i)
322 {
323 this_reader -> cartesian = TRUE;
324 if (fti < 2)
325 {
326 res = open_xyz_file (i);
327 }
328 else if (fti == 2)
329 {
330 res = open_c3d_file (i);
331 }
332 else if (fti < 5)
333 {
334 res = open_trj_file (i);
335 }
336 else if (fti < 7)
337 {
338 res = open_vas_file (i);
339 }
340 else if (fti > 6 && fti < 9)
341 {
342 res = open_pdb_file (i);
343 }
344 else if (fti == 9 || fti == 10)
345 {
346 if (fti == 10) cif_use_symmetry_positions = TRUE;
347 this_reader -> cartesian = FALSE;
348 res = open_cif_file (i);
349 }
350 else if (fti == 11)
351 {
352 res = open_hist_file (i);
353 }
354 }
355 else
356 {
357 res = 1;
358 }
359#ifndef OPENMP
360 if (tail) g_free (tail);
361#endif
362 if (! res)
363 {
364 if (fti == 9)
365 {
366 if (! this_reader -> cartesian)
367 {
368 // this_reader -> lattice.sp_group -> sid = 2;
369 // get_origin (this_reader -> lattice.sp_group);
371 {
372 i = build_crystal (FALSE, active_project, TRUE, FALSE, & this_reader -> lattice, MainWindow);
373 if (! i)
374 {
375 add_reader_info ("Error trying to build crystal using the CIF file parameters !\n"
376 "This usually comes from: \n"
377 "\t - incorrect space group description\n"
378 "\t - incomplete space group description\n"
379 "\t - missing space group setting\n"
380 "\t - incorrect space group setting\n", 0);
381 res = 3;
382 }
383 else if (i > 1)
384 {
385 add_reader_info ("Potential issue(s) when building crystal !\n"
386 "This usually comes from: \n"
387 "\t - incorrect space group description\n"
388 "\t - incomplete space group description\n"
389 "\t - missing space group setting\n"
390 "\t - incorrect space group setting\n", 1);
391 if (this_reader -> num_sym_pos)
392 {
393 add_reader_info ("\nAnother model will be built using included symmetry positions\n", 1);
395 }
396 }
397 }
398 }
399 }
400 for (i=0; i<active_project -> steps; i++)
401 {
402 for (j=0; j<active_project -> natomes; j++)
403 {
404 active_project -> atoms[i][j].id = j;
405 active_project -> atoms[i][j].show[0] = TRUE;
406 active_project -> atoms[i][j].show[1] = TRUE;
407 active_project -> atoms[i][j].label[0] = FALSE;
408 active_project -> atoms[i][j].label[1] = FALSE;
409 active_project -> atoms[i][j].pick[0] = FALSE;
410 active_project -> atoms[i][j].cloned = FALSE;
411 }
412 }
413 if (fti != 9 || this_reader -> cartesian)
414 {
415 active_project -> nspec = this_reader -> nspec;
418 k = l = 0;
419 reader_info (coord_files_ext[fti], "Number of species", active_project -> nspec);
420 for (i=0; i<active_project -> nspec; i++)
421 {
422 j = (int)this_reader -> z[i];
423 if (this_reader -> z[i] < 1.0)
424 {
425 active_chem -> label[i] = g_strdup_printf ("%s", this_reader -> dummy[l]);
426 active_chem -> element[i] = g_strdup_printf ("Dummy %s", this_reader -> dummy[l]);
427 active_chem -> chem_prop[CHEM_M][i] = 1.0;
428 active_chem -> chem_prop[CHEM_R][i] = 0.5;
429 l ++;
430 }
431 else
432 {
433 active_chem -> label[i] = g_strdup_printf ("%s", periodic_table_info[j].lab);
434 active_chem -> element[i] = g_strdup_printf ("%s", periodic_table_info[j].name);
435 active_chem -> chem_prop[CHEM_M][i] = set_mass_ (& j);
436 active_chem -> chem_prop[CHEM_R][i] = set_radius_ (& j, & k);
437 if (! active_chem -> chem_prop[CHEM_R][i])
438 {
439 gchar * str = g_strdup_printf ("For species %s, radius is equal to 0.0 !", active_chem -> label[i]);
440 add_reader_info (str, 1);
441 g_free (str);
442 }
443 active_chem -> chem_prop[CHEM_N][i] = set_neutron_ (& j);
444 active_chem -> chem_prop[CHEM_X][i] = active_chem -> chem_prop[CHEM_Z][i];
445 }
446 active_chem -> nsps[i] = this_reader -> nsps[i];
447 g_print ("Reading coordinates [%s]:\t %s, nsps[%d]= %d\n", coord_files_ext[fti], active_chem -> label[i], i+1, active_chem -> nsps[i]);
448 active_chem -> chem_prop[CHEM_Z][i] = this_reader -> z[i];
449 }
450 }
451 else
452 {
453 reader_info (coord_files_ext[fti], "Number of species", active_project -> nspec);
454 for (i=0; i<active_project -> nspec; i++)
455 {
456 g_print ("Reading coordinates [%s]:\t %s, nsps[%d]= %d\n", coord_files_ext[fti], active_chem -> label[i], i+1, active_chem -> nsps[i]);
457 }
458 }
459 }
460 if (! (fti == 9 && cif_use_symmetry_positions) || res)
461 {
462 if (cif_search)
463 {
464 g_free (cif_search);
465 cif_search = NULL;
466 }
467 if (cif_object)
468 {
469 g_free (cif_object);
470 cif_object = NULL;
471 }
473 }
474 return res;
475}
gchar * mot[2][2]
Definition popup.c:3293
Binding to the Fortran90 subroutines.
double set_radius_(int *, int *)
double set_neutron_(int *)
double set_mass_(int *)
gchar * substitute_string(gchar *init, gchar *o_motif, gchar *n_motif)
substitute all patterns in string
Definition w_library.c:372
Callback declarations for main window.
Function declarations for the crystal builder.
integer(kind=c_int) function chemistry()
Definition chemistry.F90:22
dummy_atom * dummy
Definition cpmd_atoms.c:73
void label(cairo_t *cr, double val, int axe, int p, project *this_proj)
draw axis label
Definition labels.c:56
int atoms[NUM_STYLES][2]
float val
Definition dlp_init.c:117
int activep
Definition global.c:159
GtkWidget * MainWindow
Definition global.c:201
gchar ** duplicate_strings(int num, gchar **old_val)
copy a list of strings
Definition global.c:544
gboolean cif_use_symmetry_positions
Definition global.c:189
Global variable declarations Global convenience function declarations Global data structure defin...
element_data periodic_table_info[]
Definition w_library.c:71
#define CHEM_N
Definition global.h:300
chemical_data * active_chem
Definition project.c:48
#define CHEM_R
Definition global.h:299
#define CHEM_M
Definition global.h:298
char * coord_files_ext[NCFORMATS+1]
Definition callbacks.c:100
#define LINE_SIZE
Definition global.h:486
#define CHEM_X
Definition global.h:301
#define CHEM_Z
Definition global.h:297
project * active_project
Definition project.c:47
Variable declarations related to the OpenGL window Function declarations related to the OpenGL wind...
gboolean ask_yes_no(gchar *title, gchar *text, int type, GtkWidget *widg)
ask yes or no for something: prepare dialog
Definition interface.c:356
Messaging function declarations.
integer(kind=c_int) function lattice(totl, lid, vectors, vmod, angles, lat, cfrac, apbc)
Definition lattice.F90:162
double z
Definition ogl_draw.c:57
Function declarations for reading atomes project file Function declarations for saving atomes proje...
void active_project_changed(int id)
change the active project
Definition update_p.c:175
int open_coord_file(gchar *filename, int fti)
open atomic coordinates file
Definition read_coord.c:250
int open_trj_file(int linec)
open CPMD file
Definition read_trj.c:251
FILE * coordf
Definition read_coord.c:70
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
int build_crystal(gboolean visible, project *this_proj, gboolean to_wrap, gboolean show_clones, cell_info *cell, GtkWidget *widg)
build crystal
void reader_info(gchar *type, gchar *sinf, int val)
display reader information
Definition read_coord.c:101
int open_c3d_file(int linec)
open C3D file
Definition read_c3d.c:322
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
chemical_data * alloc_chem_data(int spec)
allocate chemistry data
Definition open_p.c:186
int open_hist_file(int linec)
open DL-POLY history file
Definition read_hist.c:440
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
void allocatoms(project *this_proj)
allocate project data
Definition open_p.c:160
int open_cif_file(int linec)
open CIF file
Definition read_cif.c:1728
gchar * this_line
Definition read_coord.c:73
void check_for_species(double v, int ato)
Fill the species for each atom and the associated data.
Definition read_coord.c:196
const gchar * dfi[2]
Definition main.c:74
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
atom_search * cif_search
Definition read_cif.c:281
int open_xyz_file(int linec)
open XYZ file
Definition read_xyz.c:303
int open_vas_file(int linec)
open VASP file
Definition read_vas.c:244
int open_pdb_file(int linec)
open PDB file
Definition read_pdb.c:217
atomic_object * cif_object
Definition read_cif.c:282
int id
Definition global.h:941
int status
Definition w_advance.c:160
GtkWidget * res[2]
Definition w_encode.c:212
int element
Definition w_periodic.c:61
GtkWidget * lab
Definition workspace.c:73