atomes 1.2.1
atomes: an atomic scale modeling tool box
Loading...
Searching...
No Matches
read_isaacs.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
23/*
24* This file: 'read_isaacs.c'
25*
26* Contains:
27*
28
29 - The functions to read an ISAACS XML file
30 - The functions to write an ISAACS XML file
31
32*
33* List of functions:
34
35 int XmlwriterFilename (const char *uri);
36 int write_xml (const char * filetosave);
37 int get_spec_from_data (xmlChar * data);
38 int setprop (xmlNodePtr pnode);
39 int testopening (char * tdata, char * tfichier);
40 int setchemistry (xmlNodePtr xsnode);
41 int setbox (xmlNodePtr boxnode);
42 int setpbc (xmlNodePtr pbcnode);
43 int setcutoffs (xmlNodePtr cutnode);
44 int settime(xmlNodePtr timenode);
45 int check_xml (const char * filetocheck);
46
47 size_t strfind (char * ida);
48
49 gboolean file_exists(const char * filename);
50
51 gchar * open_xml (const char * filetoread);
52
53 xmlNodePtr findnode (xmlNodePtr startnode, char * nname);
54
55*/
56
57#include <libxml/encoding.h>
58#include <libxml/xmlwriter.h>
59#include <libxml/xmlreader.h>
60#include <libxml/parser.h>
61
62#include "global.h"
63#include "interface.h"
64#include "callbacks.h"
65
66extern int open_coordinate_file (int id);
67
68#define NFORMATS 7
69
70char * reg_types[NFORMATS] = {"XYZ file",
71 "Chem3D file",
72 "CPMD trajectory",
73 "VASP trajectory",
74 "multiple XYZ file",
75 "multiple Chem3D file",
76 "PDB file"};
77
85size_t strfind (char * ida)
86{
87 size_t a, b, c;
88
89 a = strlen(ida);
90 b = 0;
91 for ( c=0 ; c < a ; c++)
92 {
93 if (ida[c] != ' ')
94 {
95 b ++;
96 }
97 }
98 return b;
99}
100
108int XmlwriterFilename (const char * uri)
109{
110 int rc;
111 int i, j;
112 xmlTextWriterPtr writer;
113
114 char isaacinfo[16] = " I.S.A.A.C.S. v";
115 char xmlinfo[11]=" XML file ";
116 const xmlChar intro[50]="";
117 char * val;
118 char * ncut;
119 size_t lgt;
120
121 /* Create a new XmlWriter for uri, with no compression. */
122 writer = xmlNewTextWriterFilename(uri, 0);
123 if (writer == NULL) return 0;
124 rc = xmlTextWriterSetIndent(writer, 1);
125 if (rc < 0) return 0;
126 /* Start the document with the xml default for the version,
127 * encoding MY_ENCODING and the default for the standalone
128 * declaration. */
129 rc = xmlTextWriterStartDocument(writer, NULL, MY_ENCODING, NULL);
130 if (rc < 0) return 0;
131
132
133 strcpy((char *)intro, isaacinfo);
134 strcat((char *)intro, "1.0");
135 strcat((char *)intro, xmlinfo);
136 rc = xmlTextWriterWriteComment(writer, intro);
137 if (rc < 0) return 0;
138 rc = xmlTextWriterStartElement(writer, BAD_CAST "isaacs-xml");
139 if (rc < 0) return 0;
140
141 /* Start the "data" element. Since this is the first
142 * element, this will be the root element of the document.
143 <data>
144 <type> </type>
145 <file> <file>
146 </data>*/
147
148 rc = xmlTextWriterWriteComment(writer, (const xmlChar *)" Data format and file containing the configuration(s) ");
149 if (rc < 0) return 0;
150 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"data");
151 if (rc < 0) return 0;
152 int reg = 0;
153 if (active_project -> tfile < 2)
154 {
155 reg = (active_project -> steps > 1) ? 4 : 0;
156 }
157 else if (active_project -> tfile == 2)
158 {
159 reg = (active_project -> steps > 1) ? 5 : 1;
160 }
161 else if (active_project -> tfile < 5)
162 {
163 reg = 2;
164 }
165 else if (active_project -> tfile < 6)
166 {
167 reg = 3;
168 }
169 else if (active_project -> tfile < 9)
170 {
171 reg = 6;
172 }
173 else
174 {
175 return 0;
176 }
177
178 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"type", "%s", reg_types[reg]);
179 if (rc < 0) return 0;
180 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"file", "%s", active_project -> coordfile);
181 if (rc < 0) return 0;
182 rc = xmlTextWriterEndElement(writer);
183 if (rc < 0) return 0;
184
185 /* Start the "chemistry" element.
186 <chemistry>
187 <atoms> </atoms>
188 <species number=" ">
189 <active_project -> label id="0"> </active_project -> label>
190 ...
191 </species>
192 <element symbol="">
193 <name> </name>
194 <z> </z>
195 <mass> </mass>
196 <rad> </rad>
197 <radius> </radius>
198 <nscatt> </nscatt>
199 <xscatt> </xscatt>
200 </element>
201 </chemistry> */
202
203 rc = xmlTextWriterWriteComment(writer, (const xmlChar *)" Chemistry information ");
204 if (rc < 0) return 0;
205 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"chemistry");
206 if (rc < 0) return 0;
207 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"atoms", "%d", active_project -> natomes);
208 if (rc < 0) return 0;
209 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"species");
210 if (rc < 0) return 0;
211 val = g_strdup_printf("%d",active_project -> nspec);
212 rc = xmlTextWriterWriteAttribute(writer, BAD_CAST (const xmlChar *)"number", BAD_CAST val);
213 g_free (val);
214 if (rc < 0) return 0;
215 for ( i=0 ; i<active_project -> nspec ; i++)
216 {
217 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"label");
218 if (rc < 0) return 0;
219 val = g_strdup_printf("%d",i);
220 rc = xmlTextWriterWriteAttribute(writer, BAD_CAST (const xmlChar *)"id", BAD_CAST val);
221 g_free (val);
222 if (rc < 0) return 0;
223 rc = xmlTextWriterWriteFormatString (writer, "%s", active_chem -> label[i]);
224 if (rc < 0) return 0;
225 rc = xmlTextWriterEndElement(writer);
226 if (rc < 0) return 0;
227 }
228 rc = xmlTextWriterEndElement(writer);
229 if (rc < 0) return 0;
230 for ( i=0 ; i<active_project -> nspec ; i++)
231 {
232 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"element");
233 if (rc < 0) return 0;
234 rc = xmlTextWriterWriteAttribute(writer, BAD_CAST (const xmlChar *)"symbol", BAD_CAST active_chem -> label[i]);
235 if (rc < 0) return 0;
236 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"name", "%s ", active_chem -> element[i]);
237 if (rc < 0) return 0;
238 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"z", "%d", (int)active_chem -> chem_prop[CHEM_Z][i]);
239 if (rc < 0) return 0;
240 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"mass", "%f", active_chem -> chem_prop[CHEM_M][i]);
241 if (rc < 0) return 0;
242 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"rad", "%d", 0);
243 if (rc < 0) return 0;
244 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"radius", "%f", active_chem -> chem_prop[CHEM_R][i]);
245 if (rc < 0) return 0;
246 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"nscatt", "%f", active_chem -> chem_prop[CHEM_N][i]);
247 if (rc < 0) return 0;
248 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"xscatt", "%f", active_chem -> chem_prop[CHEM_Z][i]);
249 if (rc < 0) return 0;
250 rc = xmlTextWriterEndElement(writer);
251 if (rc < 0) return 0;
252 }
253 rc = xmlTextWriterEndElement(writer);
254 if (rc < 0) return 0;
255
256 /* Start the "box" element.
257 <box>
258 <edges>
259 <a></a>
260 <b></b>
261 <c></c>
262 </edges>
263 <angles> </angles>
264 <vectors>
265 <a.x></a.x>
266 <a.y></a.y>
267 <a.z></a.z>
268 <b.x></b.x>
269 <b.y></b.y>
270 <b.z></b.z>
271 <c.x></c.x>
272 <c.y></c.y>
273 <c.z></c.z>
274 </vectors>
275 </box> */
276
277 rc = xmlTextWriterWriteComment(writer, (const xmlChar *)" Box information ");
278 if (rc < 0) return 0;
279 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"box");
280 if (rc < 0) return 0;
281 if (active_cell -> ltype == 1)
282 {
283 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"edges");
284 if (rc < 0) return 0;
285 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"a", "%f", active_box -> param[0][0]);
286 if (rc < 0) return 0;
287 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"b", "%f", active_box -> param[0][1]);
288 if (rc < 0) return 0;
289 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"c", "%f", active_box -> param[0][2]);
290 if (rc < 0) return 0;
291 rc = xmlTextWriterEndElement(writer);
292 if (rc < 0) return 0;
293 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"angles");
294 if (rc < 0) return 0;
295 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"alpha", "%f", active_box -> param[1][0]);
296 if (rc < 0) return 0;
297 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"beta", "%f", active_box -> param[1][1]);
298 if (rc < 0) return 0;
299 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"gamma", "%f", active_box -> param[1][2]);
300 if (rc < 0) return 0;
301 rc = xmlTextWriterEndElement(writer);
302 if (rc < 0) return 0;
303 }
304 else
305 {
306 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"edges");
307 if (rc < 0) return 0;
308 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"a", "%f", 0.0);
309 if (rc < 0) return 0;
310 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"b", "%f", 0.0);
311 if (rc < 0) return 0;
312 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"c", "%f", 0.0);
313 if (rc < 0) return 0;
314 rc = xmlTextWriterEndElement(writer);
315 if (rc < 0) return 0;
316 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"angles");
317 if (rc < 0) return 0;
318 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"alpha", "%f", 0.0);
319 if (rc < 0) return 0;
320 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"beta", "%f", 0.0);
321 if (rc < 0) return 0;
322 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"gamma", "%f", 0.0);
323 if (rc < 0) return 0;
324 rc = xmlTextWriterEndElement(writer);
325 if (rc < 0) return 0;
326 }
327 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"vectors");
328 if (rc < 0) return 0;
329 if (active_cell -> ltype == 2)
330 {
331 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"a.x", "%f", active_box -> vect[0][0]);
332 if (rc < 0) return 0;
333 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"a.y", "%f", active_box -> vect[0][1]);
334 if (rc < 0) return 0;
335 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"a.z", "%f", active_box -> vect[0][2]);
336 if (rc < 0) return 0;
337 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"b.x", "%f", active_box -> vect[1][0]);
338 if (rc < 0) return 0;
339 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"b.y", "%f", active_box -> vect[1][1]);
340 if (rc < 0) return 0;
341 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"b.z", "%f", active_box -> vect[1][2]);
342 if (rc < 0) return 0;
343 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"c.x", "%f", active_box -> vect[2][0]);
344 if (rc < 0) return 0;
345 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"c.y", "%f", active_box -> vect[2][1]);
346 if (rc < 0) return 0;
347 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"c.z", "%f", active_box -> vect[2][2]);
348 if (rc < 0) return 0;
349 }
350 else
351 {
352 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"a.x", "%f", 0.0);
353 if (rc < 0) return 0;
354 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"a.y", "%f", 0.0);
355 if (rc < 0) return 0;
356 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"a.z", "%f", 0.0);
357 if (rc < 0) return 0;
358 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"b.x", "%f", 0.0);
359 if (rc < 0) return 0;
360 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"b.y", "%f", 0.0);
361 if (rc < 0) return 0;
362 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"b.z", "%f", 0.0);
363 if (rc < 0) return 0;
364 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"c.x", "%f", 0.0);
365 if (rc < 0) return 0;
366 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"c.y", "%f", 0.0);
367 if (rc < 0) return 0;
368 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"c.z", "%f", 0.0);
369 if (rc < 0) return 0;
370 }
371 rc = xmlTextWriterEndElement(writer);
372 if (rc < 0) return 0;
373 rc = xmlTextWriterEndElement(writer);
374 if (rc < 0) return 0;
375 /* Start the "pbc" element.
376 <pbc>
377 <apply></apply>
378 <fractional></fractional>
379 <fractype></fractype>
380 </pbc>
381 */
382
383 rc = xmlTextWriterWriteComment(writer, (const xmlChar *)" PBC information ");
384 if (rc < 0) return 0;
385 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"pbc");
386 if (rc < 0) return 0;
387 if (active_cell -> pbc)
388 {
389 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"apply", "%s", "TRUE");
390 }
391 else
392 {
393 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"apply", "%s", "FALSE");
394 }
395 if (rc < 0) return 0;
396 if (active_cell -> frac)
397 {
398 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"fractional", "%s", "TRUE");
399 }
400 else
401 {
402 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"fractional", "%s", "FALSE");
403 }
404 if (rc < 0) return 0;
405 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"fractype", "%d", active_cell -> frac);
406 if (rc < 0) return 0;
407 rc = xmlTextWriterEndElement(writer);
408 if (rc < 0) return 0;
409
410 /* Start the "cutoffs" element.
411 <cutoffs>
412 <total></total>
413 <partials>
414 <active_chem -> label[spi]-active_chem -> label[spj]></active_chem -> label[spi]-active_chem -> label[spj]>
415 </partials>
416 </cutoffs>
417 */
418
419 rc = xmlTextWriterWriteComment(writer, (const xmlChar *)" Cutoffs information ");
420 if (rc < 0) return 0;
421 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"cutoffs");
422 if (rc < 0) return 0;
423 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"total", "%f", active_chem -> grtotcutoff);
424 if (rc < 0) return 0;
425 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"partials");
426 if (rc < 0) return 0;
427 for ( i=0 ; i<active_project -> nspec ; i++ )
428 {
429 for ( j=0 ; j<active_project -> nspec ; j++ )
430 {
431 lgt = 10;
432 lgt += strfind (active_chem -> label[i]) + strlen("-") + strfind(active_chem -> label[j]);
433 ncut = g_malloc0 (lgt*sizeof*ncut);
434 strncpy (ncut, active_chem -> label[i], strfind(active_chem -> label[i]));
435 strcat (ncut, "-");
436 strncat (ncut, active_chem -> label[j], strfind(active_chem -> label[j]));
437 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST ncut, "%f", active_chem -> cutoffs[i][j]);
438 if (rc < 0) return 0;
439 g_free (ncut);
440 ncut=NULL;
441 }
442 }
443 rc = xmlTextWriterEndElement(writer);
444 if (rc < 0) return 0;
445 rc = xmlTextWriterEndElement(writer);
446 if (rc < 0) return 0;
447
448 /* Start the "time series" element.
449 <time-series>
450 <dt></dt>
451 <unit></unit>
452 <ndt></ndt>
453 </time-series>
454 */
455
456 if (active_project -> steps > 1)
457 {
458 rc = xmlTextWriterWriteComment(writer, (const xmlChar *)" Time series ");
459 if (rc < 0) return 0;
460 rc = xmlTextWriterStartElement(writer, BAD_CAST (const xmlChar *)"time-series");
461 if (rc < 0) return 0;
462 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"dt", "%f", 0.0);
463 if (rc < 0) return 0;
464 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"unit", "t [%s]", untime[active_project -> tunit]);
465 if (rc < 0) return 0;
466 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"ndt", "%d", 1);
467 if (rc < 0) return 0;
468 rc = xmlTextWriterEndElement(writer);
469 if (rc < 0) return 0;
470 }
471
472 /* Start the "project" element.
473 <project></project>
474 */
475
476 rc = xmlTextWriterWriteComment(writer, (const xmlChar *)" Apply project ");
477 if (rc < 0) return 0;
478 if (active_project -> run)
479 {
480 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"project", "%s", "TRUE");
481 }
482 else
483 {
484 rc = xmlTextWriterWriteFormatElement(writer, BAD_CAST (const xmlChar *)"project", "%s", "FALSE");
485 }
486 if (rc < 0) return 0;
487// Closing "</isaacs-xml>"
488 rc = xmlTextWriterEndElement(writer);
489 if (rc < 0) return 0;
490
491 rc = xmlTextWriterEndDocument(writer);
492 if (rc < 0) return 0;
493
494 xmlFreeTextWriter(writer);
495 return 1;
496}
497
505int write_xml (const char * filetosave)
506{
507 /* first, the file version */
508 int res=XmlwriterFilename(filetosave);
509 /*
510 * Cleanup function for the XML library.
511 */
512 xmlCleanupParser();
513 /*
514 * this is to debug memory for regression tests
515 */
516 xmlMemoryDump();
517 return res;
518}
519
527gboolean file_exists(const char * filename)
528{
529 FILE * file;
530 if ((file = fopen(filename, "r")))
531 {
532 fclose(file);
533 return TRUE;
534 }
535 return FALSE;
536}
537
546xmlNodePtr findnode (xmlNodePtr startnode, char * nname)
547{
548 xmlNodePtr tmp;
549
550 tmp = startnode;
551 while (g_strcmp0 ((char *)(tmp -> name), nname) != 0)
552 {
553 if (tmp -> next == NULL)
554 {
555 return NULL;
556 }
557 tmp = tmp -> next;
558 }
559 return tmp;
560}
561
569int get_spec_from_data (xmlChar * data)
570{
571 int i;
572 for ( i=0 ; i<active_project -> nspec ; i++)
573 {
574 if (g_strcmp0 (exact_name((char *)data), active_chem -> label[i]) == 0)
575 {
576 return i;
577 }
578 }
579 return -1;
580}
581
589int setprop (xmlNodePtr pnode)
590{
591 int i, res;
592 xmlAttrPtr pspec;
593 xmlNodePtr idn, ptnode, ptn;
594 xmlChar * data;
595
596 ptnode=pnode;
597 for ( i=0 ; i<active_project -> nspec ; i++)
598 {
599 ptnode = findnode(ptnode, "element");
600 if (ptnode == NULL)
601 {
602 res=6;
603 goto pend;
604 }
605 pspec = ptnode -> properties;
606 idn = pspec -> children;
607 data = xmlNodeGetContent(idn);
608 res = get_spec_from_data(data);
609 if (res < 0)
610 {
611 res=6;
612 goto pend;
613 }
614 ptnode = ptnode -> children;
615 ptn = findnode(ptnode, "name");
616 if (ptn == NULL)
617 {
618 res=6;
619 goto pend;
620 }
621 data = xmlNodeGetContent(ptn);
622 if (g_strcmp0 (exact_name((char *)data), active_chem -> element[res]) != 0)
623 {
624 res=6;
625 goto pend;
626 }
627 ptn = findnode(ptnode, "z");
628 if (ptn == NULL)
629 {
630 res=6;
631 goto pend;
632 }
633 data= xmlNodeGetContent(ptn);
634 if ((int)active_chem -> chem_prop[CHEM_Z][res] != (int)string_to_double ((gpointer)data))
635 {
636 res=6;
637 goto pend;
638 }
639 ptn = findnode(ptnode, "mass");
640 if (ptn == NULL)
641 {
642 res=6;
643 goto pend;
644 }
645 data= xmlNodeGetContent(ptn);
646 active_chem -> chem_prop[CHEM_M][res] = string_to_double ((gpointer)data);
647 ptn = findnode(ptnode, "rad");
648 if (ptn == NULL)
649 {
650 res=6;
651 goto pend;
652 }
653 data= xmlNodeGetContent(ptn);
654 // val = string_to_double ((gpointer)data);
655 ptn = findnode(ptnode, "radius");
656 if (ptn == NULL)
657 {
658 res=6;
659 goto pend;
660 }
661 data= xmlNodeGetContent(ptn);
662 active_chem -> chem_prop[CHEM_R][res] = string_to_double ((gpointer)data);
663 ptn = findnode(ptnode, "nscatt");
664 if (ptn == NULL)
665 {
666 res=6;
667 goto pend;
668 }
669 data= xmlNodeGetContent(ptn);
670 active_chem -> chem_prop[CHEM_N][res] = string_to_double ((gpointer)data);
671 ptn = findnode(ptnode, "xscatt");
672 if (ptn == NULL)
673 {
674 res=6;
675 goto pend;
676 }
677 data= xmlNodeGetContent(ptn);
678 active_chem -> chem_prop[CHEM_X][res] = string_to_double ((gpointer)data);
679 ptnode = ptnode -> parent;
680 ptnode = ptnode -> next;
681 }
682 res = 0;
683pend:
684 return res;
685}
686
695int testopening (char * tdata, char * tfichier)
696{
697 int i, j;
698 char * err;
699
700 j = -1;
701 for ( i=0 ; i<NFORMATS ; i ++ )
702 {
703 if (g_strcmp0 (reg_types[i], tdata) == 0)
704 {
705 j = i;
706 break;
707 }
708 }
709 if (j > -1)
710 {
711 if (file_exists (tfichier))
712 {
713 switch (j)
714 {
715 case 0:
716 active_project -> tfile = 0;
717 break;
718 case 1:
719 active_project -> tfile = 2;
720 break;
721 case 2:
722 active_project -> tfile = 3;
723 break;
724 case 3:
725 active_project -> tfile = 5;
726 break;
727 case 4:
728 active_project -> tfile = 0;
729 break;
730 case 5:
731 active_project -> tfile = 2;
732 break;
733 case 6:
734 active_project -> tfile = 7;
735 break;
736 }
737 active_project -> coordfile = tfichier;
738 if (open_coordinate_file (active_project -> tfile) == 0)
739 {
740 return 0;
741 }
742 else
743 {
744 return 11;
745 }
746 }
747 else
748 {
749 j=4;
750 err=g_strdup_printf("The %s\n %s\n from the XML file does not exist\n", tdata, tfichier);
751 show_error (err, 0, MainWindow);
752 g_free (err);
753 return j;
754 }
755 }
756 else
757 {
758 show_error ("Unknown data format in XML file\n", 0, MainWindow);
759 j=4;
760 return j;
761 }
762}
763
771int setchemistry (xmlNodePtr xsnode)
772{
773 xmlNodePtr idn, cs, xnode;
774 xmlAttrPtr xspec;
775 xmlChar * data;
776 int ats, res, i;
777 double val;
778
779 xnode = xsnode -> children;
780 cs = findnode(xnode, "atoms");
781 if (cs == NULL)
782 {
783 res=5;
784 goto xend;
785 }
786 data = xmlNodeGetContent(cs);
787 val = string_to_double ((gpointer)data);
788 ats = val;
789 if (ats != active_project -> natomes)
790 {
791 show_warning ("The number of atoms in the XML file\n"
792 "is not the same that the number of atoms\n"
793 "in the file containing the coordinates.\n"
794 "Other information from the XML file\n"
795 "will be ignored\n", MainWindow);
796 res=5;
797 goto xend;
798 }
799 else
800 {
801 cs = findnode(xnode, "species");
802 if (cs == NULL)
803 {
804 res=5;
805 goto xend;
806 }
807 xspec = cs -> properties;
808 idn = xspec -> children;
809 data = xmlNodeGetContent(idn);
810 val = string_to_double ((gpointer)data);
811 ats = val;
812 if (ats != active_project -> nspec)
813 {
814 show_warning ("The number of chemical species in the XML file\n"
815 "is not the same that the number of chemical species\n"
816 "in the file that contains the atomic coordinates.\n", MainWindow);
817 res = 5;
818 goto xend;
819 }
820 else
821 {
822 cs = cs -> children;
823 ats = 0;
824 for ( i=0 ; i<active_project -> nspec ; i++)
825 {
826 cs = findnode(cs, "label");
827 if (cs == NULL)
828 {
829 res=5;
830 goto xend;
831 }
832 data = xmlNodeGetContent(cs);
833 res = get_spec_from_data (data);
834 if (res < 0)
835 {
836 res=5;
837 goto xend;
838 }
839 else
840 {
841 ats ++;
842 }
843 }
844 if (ats != active_project -> nspec)
845 {
846 res=5;
847 goto xend;
848 }
849 xnode = xsnode -> children;
850 res = setprop (xnode);
851 }
852 }
853
854xend:
855 return res;
856}
857
865int setbox (xmlNodePtr boxnode)
866{
867 int box;
868 xmlNodePtr bnode, ba, bb;
869 xmlChar * data;
870 double val;
871
872 bnode = boxnode -> children;
873 ba=findnode(bnode, "edges");
874 if (ba == NULL)
875 {
876 box=7;
877 goto bend;
878 }
879 ba = ba -> children;
880 bb=findnode(ba, "a");
881 if (bb == NULL)
882 {
883 box=7;
884 goto bend;
885 }
886 data = xmlNodeGetContent(bb);
887 val = string_to_double ((gpointer)data);
888 active_box -> param[0][0]= val;
889 bb=findnode(ba, "b");
890 if (bb == NULL)
891 {
892 box=7;
893 goto bend;
894 }
895 data = xmlNodeGetContent(bb);
896 val = string_to_double ((gpointer)data);
897 active_box -> param[0][1]= val;
898 bb=findnode(ba, "c");
899 if (bb == NULL)
900 {
901 box=7;
902 goto bend;
903 }
904 data = xmlNodeGetContent(bb);
905 val = string_to_double ((gpointer)data);
906 active_box -> param[0][2]= val;
907 ba=findnode(bnode, "angles");
908 if (ba == NULL)
909 {
910 box=7;
911 goto bend;
912 }
913 ba = ba -> children;
914 bb=findnode(ba, "alpha");
915 if (bb == NULL)
916 {
917 box=7;
918 goto bend;
919 }
920 data = xmlNodeGetContent(bb);
921 val = string_to_double ((gpointer)data);
922 active_box -> param[1][0] = val;
923 bb=findnode(ba, "beta");
924 if (bb == NULL)
925 {
926 box=7;
927 goto bend;
928 }
929 data = xmlNodeGetContent(bb);
930 val = string_to_double ((gpointer)data);
931 active_box -> param[1][1]= val;
932 bb=findnode(ba, "gamma");
933 if (bb == NULL)
934 {
935 box=7;
936 goto bend;
937 }
938 if (active_box -> param[0][0] != 0.0 && active_box -> param[0][1] != 0.0 && active_box -> param[0][2] != 0.0
939 && active_box -> param[1][0] != 0.0 && active_box -> param[1][1] != 0.0 && active_box -> param[1][2] != 0.0)
940 {
941 active_cell -> ltype = 1;
942 }
943 else
944 {
945 active_cell -> ltype = 0;
946 }
947 data = xmlNodeGetContent(bb);
948 val = string_to_double ((gpointer)data);
949 active_box -> param[1][2]= val;
950 ba=findnode(bnode, "vectors");
951 if (ba == NULL)
952 {
953 box=7;
954 goto bend;
955 }
956 ba = ba -> children;
957 bb=findnode(ba, "a.x");
958 if (bb == NULL)
959 {
960 box=7;
961 goto bend;
962 }
963 data = xmlNodeGetContent(bb);
964 val = string_to_double ((gpointer)data);
965 active_box -> vect[0][0]= val;
966 bb=findnode(ba, "a.y");
967 if (bb == NULL)
968 {
969 box=7;
970 goto bend;
971 }
972 data = xmlNodeGetContent(bb);
973 val = string_to_double ((gpointer)data);
974 active_box -> vect[0][1]= val;
975 bb=findnode(ba, "a.z");
976 if (bb == NULL)
977 {
978 box=7;
979 goto bend;
980 }
981 data = xmlNodeGetContent(bb);
982 val = string_to_double ((gpointer)data);
983 active_box -> vect[0][2]= val;
984 bb=findnode(ba, "b.x");
985 if (bb == NULL)
986 {
987 box=7;
988 goto bend;
989 }
990 data = xmlNodeGetContent(bb);
991 val = string_to_double ((gpointer)data);
992 active_box -> vect[1][0]= val;
993 bb=findnode(ba, "b.y");
994 if (bb == NULL)
995 {
996 box=7;
997 goto bend;
998 }
999 data = xmlNodeGetContent(bb);
1000 val = string_to_double ((gpointer)data);
1001 active_box -> vect[1][1]= val;
1002 bb=findnode(ba, "b.z");
1003 if (bb == NULL)
1004 {
1005 box=7;
1006 goto bend;
1007 }
1008 data = xmlNodeGetContent(bb);
1009 val = string_to_double ((gpointer)data);
1010 active_box -> vect[1][2]= val;
1011 bb=findnode(ba, "c.x");
1012 if (bb == NULL)
1013 {
1014 box=7;
1015 goto bend;
1016 }
1017 data = xmlNodeGetContent(bb);
1018 val = string_to_double ((gpointer)data);
1019 active_box -> vect[2][0]= val;
1020 bb=findnode(ba, "c.y");
1021 if (bb == NULL)
1022 {
1023 box=7;
1024 goto bend;
1025 }
1026 data = xmlNodeGetContent(bb);
1027 val = string_to_double ((gpointer)data);
1028 active_box -> vect[2][1]= val;
1029 bb=findnode(ba, "c.z");
1030 if (bb == NULL)
1031 {
1032 box=7;
1033 goto bend;
1034 }
1035 data = xmlNodeGetContent(bb);
1036 val = string_to_double ((gpointer)data);
1037 active_box -> vect[2][2]= val;
1038
1039 if (active_box -> vect[0][0] != 0.0 && active_box -> vect[0][1] != 0.0 && active_box -> vect[0][2] != 0.0
1040 && active_box -> vect[1][0] != 0.0 && active_box -> vect[1][1] != 0.0 && active_box -> vect[1][2] != 0.0
1041 && active_box -> vect[2][0] != 0.0 && active_box -> vect[2][1] != 0.0 && active_box -> vect[2][2] != 0.0)
1042 {
1043 active_cell -> ltype = 2;
1044 }
1045 box = 0;
1046
1047bend:
1048 return box;
1049}
1050
1058int setpbc (xmlNodePtr pbcnode)
1059{
1060 int pbc;//, j;
1061 xmlNodePtr pnode, bnode;
1062 xmlChar * data;
1063 //double val;
1064
1065 pnode = pbcnode -> children;
1066 bnode=findnode(pnode, "apply");
1067 if (bnode == NULL)
1068 {
1069 pbc=8;
1070 goto pend;
1071 }
1072 data = xmlNodeGetContent(bnode);
1073 active_cell -> pbc = 0;
1074 if (g_strcmp0 ((char *)data, (char *)"TRUE") == 0)
1075 {
1076 active_cell -> pbc = 1;
1077 }
1078 bnode=findnode(pnode, "fractional");
1079 if (bnode == NULL)
1080 {
1081 pbc=8;
1082 goto pend;
1083 }
1084 data = xmlNodeGetContent(bnode);
1085 active_cell -> frac = 0;
1086 if (g_strcmp0 ((char *)data, (char *)"TRUE") == 0)
1087 {
1088 bnode=findnode(pnode, "fractype");
1089 if (bnode == NULL)
1090 {
1091 pbc=8;
1092 goto pend;
1093 }
1094 data = xmlNodeGetContent(bnode);
1095 //val=string_to_double ((gpointer)data);
1096 //j = val;
1097 active_cell -> frac = 1;
1098 }
1099 pbc=0;
1100
1101pend:
1102 return pbc;
1103}
1104
1112int setcutoffs (xmlNodePtr cutnode)
1113{
1114 int cut;
1115 int i, j;
1116 xmlNodePtr cnode, cn;
1117 xmlChar * data;
1118 double val;
1119 char * ncut;
1120 size_t lgt;
1121
1122 cnode = cutnode -> children;
1123 cn = findnode(cnode, "total");
1124 if (cn == NULL)
1125 {
1126 cut=9;
1127 goto cend;
1128 }
1129 data = xmlNodeGetContent(cn);
1130 val=string_to_double ((gpointer)data);
1131 active_chem -> grtotcutoff = val;
1132 cn = findnode(cnode, "partials");
1133 if (cn == NULL)
1134 {
1135 cut=9;
1136 goto cend;
1137 }
1138 cnode = cn -> children;
1139 for ( i=0 ; i<active_project -> nspec ; i++ )
1140 {
1141 for ( j=0 ; j<active_project -> nspec ; j++ )
1142 {
1143 lgt=10;
1144 lgt+=strfind(active_chem -> label[i]) + strlen("-") + strfind(active_chem -> label[j]);
1145 ncut = g_malloc0 (lgt*sizeof*ncut);
1146 strncpy (ncut, active_chem -> label[i], strfind(active_chem -> label[i]));
1147 strcat (ncut, "-");
1148 strncat (ncut, active_chem -> label[j], strfind(active_chem -> label[j]));
1149 cn = findnode(cnode, ncut);
1150 if (cn == NULL)
1151 {
1152 cut=9;
1153 goto cend;
1154 }
1155 data = xmlNodeGetContent(cn);
1156 val=string_to_double ((gpointer)data);
1157 active_chem -> cutoffs[i][j] = val;
1158 g_free (ncut);
1159 ncut=NULL;
1160 }
1161 }
1162 cut=0;
1163
1164cend:
1165 return cut;
1166}
1167
1175int settime(xmlNodePtr timenode)
1176{
1177 int i, j;
1178 //double val, delta_t;
1179 //int tunit, ndtbs;
1180 int tps;
1181 xmlChar * data;
1182 xmlNodePtr tnode, tn;
1183
1184 tnode = timenode -> children;
1185 tn = findnode(tnode, "dt");
1186 if (tn == NULL)
1187 {
1188 tps=10;
1189 goto tend;
1190 }
1191 data = xmlNodeGetContent(tn);
1192 //val = string_to_double ((gpointer)data);
1193 //delta_t = val;
1194 tn = findnode(tnode, "unit");
1195 if (tn == NULL)
1196 {
1197 tps=10;
1198 goto tend;
1199 }
1200 data = xmlNodeGetContent(tn);
1201 j=-1;
1202 for ( i=0 ; i<5 ; i++)
1203 {
1204 if (g_strcmp0 ((char *)data, g_strdup_printf ("t [%s]", untime[i])) == 0) j=i;
1205 }
1206 if (j == -1)
1207 {
1208 tps=10;
1209 goto tend;
1210 }
1211 if (j < 5)
1212 {
1213 //tunit = j+1;
1214 }
1215 tn = findnode(tnode, "ndt");
1216 if (tn != NULL)
1217 {
1218 data = xmlNodeGetContent(tn);
1219 //val = string_to_double ((gpointer)data);
1220 //ndtbs = val;
1221 }
1222 tps=0;
1223
1224tend:
1225 return tps;
1226}
1227
1235int check_xml (const char * filetocheck)
1236{
1237 int res;
1238 xmlDoc * doc;
1239 xmlTextReaderPtr reader;
1240 const xmlChar xisaacs[11]="isaacs-xml";
1241 xmlChar * cdata, * cfichier;
1242 xmlNodePtr racine, node;
1243 /*
1244 * build an xmlReader for that file
1245 */
1246 reader = xmlReaderForFile(filetocheck, NULL, 0);
1247 if (reader != NULL)
1248 {
1249 res = 0;
1250 doc = xmlParseFile(filetocheck);
1251 if (doc == NULL)
1252 {
1253 res=1;
1254 goto end;
1255 }
1256 else
1257 {
1258 racine = xmlDocGetRootElement(doc);
1259 if (g_strcmp0 ((char *)(racine->name), (char *)xisaacs) != 0)
1260 {
1261 res=2;
1262 goto end;
1263 }
1264 else
1265 {
1266 node = racine -> children;
1267 node = findnode(node, "data");
1268 if (node == NULL)
1269 {
1270 res=3;
1271 goto end;
1272 }
1273 node = node -> children;
1274 node = findnode(node, "type");
1275 if (node == NULL)
1276 {
1277 res=3;
1278 goto end;
1279 }
1280 cdata = xmlNodeGetContent(node);
1281 node = findnode(node, "file");
1282 if (node == NULL)
1283 {
1284 res=3;
1285 goto end;
1286 }
1287 cfichier = xmlNodeGetContent(node);
1288 res= testopening((char *)cdata, (char *)cfichier);
1289 if (res != 0)
1290 {
1291 goto end;
1292 }
1293 else
1294 {
1295 node = racine -> children;
1296 node = findnode(node, "chemistry");
1297 if (node == NULL)
1298 {
1299 res=3;
1300 goto end;
1301 }
1302 res= setchemistry(node);
1303 if (res != 0)
1304 {
1305 goto end;
1306 }
1307 node = racine -> children;
1308 node = findnode(node, "pbc");
1309 if (node == NULL)
1310 {
1311 res=3;
1312 goto end;
1313 }
1314 res=setpbc (node);
1315 if (res != 0)
1316 {
1317 goto end;
1318 }
1319 node = racine -> children;
1320 node = findnode(node, "box");
1321 if (node == NULL)
1322 {
1323 res=3;
1324 goto end;
1325 }
1326 active_cell -> box = g_malloc0(sizeof*active_cell -> box);
1327 active_box = & active_cell -> box[0];
1328 res=setbox (node);
1329 if (res != 0)
1330 {
1331 goto end;
1332 }
1333 node = racine -> children;
1334 node = findnode(node, "cutoffs");
1335 if (node == NULL)
1336 {
1337 res=3;
1338 goto end;
1339 }
1340 active_chem -> cutoffs = allocddouble (active_project -> nspec, active_project -> nspec);
1341 res=setcutoffs (node);
1342 if (res != 0)
1343 {
1344 goto end;
1345 }
1346 if (active_project -> steps > 1)
1347 {
1348 node = racine -> children;
1349 node = findnode(node, "time-series");
1350 if (node == NULL)
1351 {
1352 res=3;
1353 goto end;
1354 }
1355 res=settime (node);
1356 if (res != 0)
1357 {
1358 goto end;
1359 }
1360 }
1361 node = racine -> children;
1362 node = findnode(node, "project");
1363 if (node == NULL)
1364 {
1365 res=3;
1366 goto end;
1367 }
1368 cdata=xmlNodeGetContent(node);
1369 if (g_strcmp0 ((char *)cdata, (char *)"TRUE") == 0)
1370 {
1371 active_project -> run = 1;
1372 }
1373 else
1374 {
1375 active_project -> run = 0;
1376 }
1377 }
1378 }
1379 }
1380 }
1381 else
1382 {
1383 res=1;
1384 }
1385
1386end:
1387 /*
1388 * Cleanup function for the XML library.
1389 */
1390 if (reader != NULL) xmlFreeTextReader(reader);
1391 /*
1392 *Free the global variables that may
1393 *have been allocated by the parser.
1394 */
1395 xmlCleanupParser();
1396
1397 return res;
1398}
1399
1407gchar * open_xml (const char * filetoread)
1408{
1409 int oxml;
1410 oxml = check_xml(filetoread);
1411 switch(oxml)
1412 {
1413 case 1:
1414 return ("Impossible to open the Isaacs project file (XML) file\n");
1415 break;
1416 case 2:
1417 return ("The file is not a valid Isaacs project file (XML) file\n");
1418 break;
1419 case 3:
1420 return ("The Isaacs project file (XML) file is incomplete\n");
1421 break;
1422 case 4:
1423 return ("Problem(s) in the &lt;data&gt; section of the Isaacs project file (XML)\n");
1424 break;
1425 case 5:
1426 return ("Problem(s) in the &lt;chemistry&gt; section of the Isaacs project file (XML)\n");
1427 break;
1428 case 6:
1429 return ("Problem(s) in the &lt;element&gt; section of the Isaacs project file (XML)\n");
1430 break;
1431 case 7:
1432 return ("Problem(s) in the &lt;box&gt; section of the Isaacs project file (XML)\n");
1433 break;
1434 case 8:
1435 return ("Problem(s) in the &lt;pbc&gt; section of the Isaacs project file (XML)\n");
1436 break;
1437 case 9:
1438 return ("Problem(s) in the &lt;cutoffs&gt; section of the Isaacs project file (XML)\n");
1439 break;
1440 case 10:
1441 return ("Problem(s) in the &lt;time-series&gt; section of the Isaacs project file (XML)\n");
1442 break;
1443 default:
1444 return NULL;
1445 break;
1446 }
1447}
Callback declarations for main window.
gchar * param[2]
void label(cairo_t *cr, double val, int axe, int p, project *this_proj)
draw axis label
Definition labels.c:56
float val
Definition dlp_init.c:117
gchar * filetoread
char * untime[5]
Definition global.c:150
double ** allocddouble(int xal, int yal)
allocate a double ** pointer
Definition global.c:475
GtkWidget * MainWindow
Definition global.c:201
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...
cell_info * active_cell
Definition project.c:50
#define CHEM_N
Definition global.h:302
chemical_data * active_chem
Definition project.c:48
#define CHEM_R
Definition global.h:301
box_info * active_box
Definition project.c:51
#define MY_ENCODING
Definition global.h:91
#define CHEM_M
Definition global.h:300
#define CHEM_X
Definition global.h:303
#define CHEM_Z
Definition global.h:299
project * active_project
Definition project.c:47
struct box box
box layout data structure
Definition glwin.h:349
void show_warning(char *warning, GtkWidget *win)
show warning
Definition interface.c:260
void show_error(char *error, int val, GtkWidget *win)
show error message
Definition interface.c:293
gchar * exact_name(gchar *name)
short cut to print string without spaces
Definition interface.c:370
Messaging function declarations.
int open_coordinate_file(int id)
try to open coordinate file, type is based of id
Definition callbacks.c:1266
int setprop(xmlNodePtr pnode)
read chemical properties from XML node
int setbox(xmlNodePtr boxnode)
read box properties from node
xmlNodePtr findnode(xmlNodePtr startnode, char *nname)
find XML node
int settime(xmlNodePtr timenode)
read MD information from node
gboolean file_exists(const char *filename)
file exists ?
int testopening(char *tdata, char *tfichier)
test atomic coordinates file opening
int write_xml(const char *filetosave)
write XML file
char * reg_types[NFORMATS]
Definition read_isaacs.c:70
int check_xml(const char *filetocheck)
check the opening of ISAACS XML file
size_t strfind(char *ida)
size of a string without spaces
Definition read_isaacs.c:85
#define NFORMATS
Definition read_isaacs.c:68
int XmlwriterFilename(const char *uri)
write ISAACS XML file
int get_spec_from_data(xmlChar *data)
get atomic species from data
int setchemistry(xmlNodePtr xsnode)
read chemistry data from node
int setcutoffs(xmlNodePtr cutnode)
read bond cutoffs from node
gchar * open_xml(const char *filetoread)
Open ISAACS XML file.
int setpbc(xmlNodePtr pbcnode)
read the PBC information from node
Definition glwin.h:351
int b
Definition tab-1.c:95
int c
Definition tab-1.c:95
int a
Definition tab-1.c:95
GtkWidget * res[2]
Definition w_encode.c:212
int element
Definition w_periodic.c:61