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