atomes 1.1.15
atomes: an atomic scale modeling tool box
Loading...
Searching...
No Matches
dlp_print.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
25/*
26* This file: 'dlp_print.c'
27*
28* Contains:
29*
30
31 - The functions to handle the output of the DL-POLY FIELD file
32 - The functions to handle the output of the DL-POLY CONTROL file
33 - The functions to handle the output of the DL-POLY CONFIG file
34 - The functions to fill the structural element(s) tree models in the assistant
35
36*
37* List of functions:
38
39 int get_num_struct_to_print (field_molecule * fmol, int sid);
40 int get_pbc ();
41
42 gboolean print_this_imp_inv (imp_inv * inv, int di, int a, int b, int c, int d);
43 gboolean member_of_atom (field_atom* fat, int id);
44 gboolean print_ana ();
45
46 void print_field_prop (field_prop * pro, int st, field_molecule * mol);
47 void print_field_struct (field_struct * stru, field_molecule * mol);
48 void print_all_field_struct (field_molecule * mol, int str);
49 void print_dlp_improper_inversion (int di, GtkTextBuffer * buf, field_struct * dhii, int fi, GtkTreeStore * store, GtkTreeIter * iter);
50 void print_dlp_dihedral (int dih, GtkTextBuffer * buf, field_struct * dh, int fi, GtkTreeStore * store, GtkTreeIter * iter);
51 void print_dlp_angle (int ai, GtkTextBuffer * buf, field_struct * an, int fi, GtkTreeStore * store, GtkTreeIter * iter);
52 void print_dlp_bond (int bi, GtkTextBuffer * buf, field_struct * bd, int fi, GtkTreeStore * store, GtkTreeIter * iter);
53 void print_dlp_rigid (GtkTextBuffer * buf, field_rigid * rig);
54 void print_dlp_tet (GtkTextBuffer * buf, field_tethered * tet);
55 void print_dlp_pmf (GtkTextBuffer * buf, field_pmf * pmf);
56 void print_dlp_cons (GtkTextBuffer * buf, field_constraint * cons);
57 void print_dlp_shell (GtkTextBuffer * buf, field_molecule * fmol, field_shell * shell);
58 void print_dlp_atom (GtkTextBuffer * buf, int at, int numat);
59 void print_dlp_molecule (GtkTextBuffer * buf, field_molecule * fmol);
60 void print_dlp_body (GtkTextBuffer * buf, field_nth_body * body);
61 void print_dlp_tersoff_cross (GtkTextBuffer * buf, field_nth_body * body_a, field_nth_body * body_b);
62 void print_dlp_tersoff (GtkTextBuffer * buf, field_nth_body * body);
63 void print_dlp_field (GtkTextBuffer * buf);
64 void print_dlp_config (GtkTextBuffer * buf);
65 void print_int (GtkTextBuffer * buf, int data);
66 void print_control_int (GtkTextBuffer * buf, int data, gchar * info_a, gchar * info_b, gchar * key);
67 void print_float (GtkTextBuffer * buf, double data);
68 void print_control_float (GtkTextBuffer * buf, double data, gchar * info_a, gchar * info_b, gchar * key);
69 void print_sci (GtkTextBuffer * buf, double data);
70 void print_control_sci (GtkTextBuffer * buf, double data, gchar * info_a, gchar * info_b, gchar * key);
71 void print_string (GtkTextBuffer * buf, gchar * string);
72 void print_control_string (GtkTextBuffer * buf, gchar * string, gchar * info_a, gchar * info_b, gchar * key);
73 void print_control_key (GtkTextBuffer * buf, gchar * info, gchar * key);
74 void print_dlp_control (GtkTextBuffer * buf);
75
76*/
77
78#include "dlp_field.h"
79#include "interface.h"
80
81extern gboolean in_bond (int at, int bd[2]);
82extern int get_num_vdw_max ();
83extern gchar * get_body_element_name (field_nth_body * body, int aid, int nbd);
84
95{
96 int i, j, k, u, v;
97 field_atom* fat;
98 j = struct_id(st+7);
99#ifdef DEBUG
100 int w;
101 g_debug ("Prop - natomes= %d", j);
102#endif
103 for (i=0; i<j; i++)
104 {
105#ifdef DEBUG
106 g_debug (" at[%d]= %d", i, pro -> aid[i]);
107#endif
108 if (mol != NULL)
109 {
110 for (k=0; k< mol -> multi; k++)
111 {
112 u = mol -> atoms_id[pro -> aid[i]][k].a;
113 v = mol -> atoms_id[pro -> aid[i]][k].b;
114 fat = get_active_atom (mol -> id, u);
115 if (v > fat -> num)
116 {
117#ifdef DEBUG
118 g_debug ("********************** BIG BUG BIG BUG BIG BUG **********************");
119 g_debug (" TO CHEK:: multi= %d:: at.a= %d, at.b= %d", k, u, v);
120 g_debug ("********************** BIG BUG BIG BUG BIG BUG **********************");
121#endif
122 }
123 else
124 {
125#ifdef DEBUG
126 w = fat -> list[v];
127 g_debug (" multi= %d:: at.a= %d, at.b= %d, real_id= %d", k, u, v, w);
128#endif
129 }
130 }
131 }
132 }
133
134#ifdef DEBUG
135 g_debug ("Prop - key= %d", pro -> key);
136 j = fvalues[activef][st][pro -> key];
137 for (i=0; i<j; i++)
138 {
139 g_debug (" val[%d]= %f", i, pro -> val[i]);
140 }
141 g_debug ("Prop - show= %d", pro -> show);
142 g_debug ("Prop - use= %d", pro -> use);
143#endif
144}
145
155{
156#ifdef DEBUG
157 int i;
158 g_debug (" ");
159 g_debug ("Struct - st= %d", stru -> st);
160 g_debug ("Struct - id= %d", stru -> id);
161 g_debug ("Struct - num= %d", stru -> num);
162 g_debug ("Struct - av= %f", stru -> av);
163 g_debug ("Struct - natomes= %d", struct_id(stru -> st+7));
164 for (i=0; i<struct_id(stru -> st+7); i++)
165 {
166 g_debug (" at[%d]= %d", i, stru -> aid[i]);
167 }
168 g_debug ("Default property:: ");
169#endif
170 print_field_prop (stru -> def, stru -> st, NULL);
171 if (stru -> other != NULL)
172 {
173#ifdef DEBUG
174 g_debug ("Other property(ies):: ");
175#endif
176 field_prop * tmp_pr;
177 tmp_pr = stru -> other;
178 print_field_prop (tmp_pr, stru -> st, mol);
179 while (tmp_pr -> next != NULL)
180 {
181 print_field_prop (tmp_pr -> next, stru -> st, mol);
182 tmp_pr = tmp_pr -> next;
183 }
184 }
185}
186
196{
197 int i;
198 field_struct * tmp_s;
199#ifdef DEBUG
200 g_debug (" ");
201 g_debug ("IN MOL:: %d", mol -> id);
202 g_debug ("PRINTING STRUCT:: %d", str);
203 g_debug ("Total Num of Struct %d:: NUM= %d", str, mol -> nstruct[str]);
204#endif
205 tmp_s = mol -> first_struct[str];
206 for (i=0; i<mol -> nstruct[str]; i++)
207 {
208 print_field_struct (tmp_s, mol);
209 if (tmp_s -> next != NULL) tmp_s = tmp_s -> next;
210 }
211#ifdef DEBUG
212 g_debug ("END STRUCT :: %d", str);
213#endif
214}
215
216typedef struct imp_inv imp_inv;
218{
219 int a;
220 int b;
221 int c;
222 int d;
225};
226
239gboolean print_this_imp_inv (imp_inv * inv, int di, int a, int b, int c, int d)
240{
241 if (! inv) return TRUE;
242 while (inv)
243 {
244 if (di == 6)
245 {
246 if (inv -> a == a && inv -> d == d)
247 {
248 if (inv -> b == b && inv -> c == c) return FALSE;
249 if (inv -> b == c && inv -> c == b) return FALSE;
250 }
251 }
252 else
253 {
254 if (inv -> a == a && inv -> b == b) return FALSE;
255 /*
256 {
257 if (inv -> b == b && inv -> c == d && inv -> d == c) return FALSE;
258 if (inv -> b == c && inv -> c == b && inv -> d == d) return FALSE;
259 if (inv -> b == c && inv -> c == d && inv -> d == b) return FALSE;
260 if (inv -> b == d && inv -> c == b && inv -> d == c) return FALSE;
261 if (inv -> b == d && inv -> c == c && inv -> d == b) return FALSE;
262 }*/
263 }
264 inv = inv -> next;
265 }
266 return TRUE;
267}
268
277gboolean member_of_atom (field_atom* fat, int id)
278{
279 int i;
280 for (i=0; i<tmp_fat -> num; i++)
281 {
282 if (tmp_fat -> list[i] == id) return TRUE;
283 }
284 return FALSE;
285}
286
299void print_dlp_improper_inversion (int di, GtkTextBuffer * buf, field_struct * dhii, int fi, GtkTreeStore * store, GtkTreeIter * iter)
300{
301 int a, b, c, d, e, i, j, k, l, m, n, o, p, q, r, s, t, u, v;
302 gboolean show;
303 gchar * stra, * strb, * strc, * strd, * stre, * strf, * strg;
304 float w;
305 GtkTreeIter di_level;
306 int * ids = allocint(4);
307 imp_inv * first_imp_inv = NULL;
308 imp_inv * this_ii = NULL;
309
310 for (i=0; i<tmp_fat -> num; i++)
311 {
312 j = tmp_fat -> list[i];
313 if (tmp_proj -> atoms[0][j].coord[2] == fi)
314 {
315 if ((tmp_proj -> atoms[0][j].numv > 2 && di == 6) || (tmp_proj -> atoms[0][j].numv == 3 && di == 7))
316 {
317 a = ids[0] = tmp_fat -> list_id[i];
318 for (k=0; k<tmp_proj -> atoms[0][j].numv; k++)
319 {
320 l = tmp_proj -> atoms[0][j].vois[k];
321 if (tmp_proj -> atoms[0][l].faid == tmp_fbt -> id)
322 {
323 b = ids[1] = tmp_proj -> atoms[0][l].fid;
324 for (m=0; m<tmp_proj -> atoms[0][j].numv; m++)
325 {
326 if (m != k)
327 {
328 n = tmp_proj -> atoms[0][j].vois[m];
329 if (tmp_proj -> atoms[0][n].faid == tmp_fct -> id)
330 {
331 c = ids[2] = tmp_proj -> atoms[0][n].fid;
332 for (o=0; o<tmp_proj -> atoms[0][j].numv; o++)
333 {
334 if (o != k && o != m)
335 {
336 p = tmp_proj -> atoms[0][j].vois[o];
337 if (tmp_proj -> atoms[0][p].faid == tmp_fdt -> id)
338 {
339 d = ids[3] = tmp_proj -> atoms[0][p].fid;
340 if (print_this_imp_inv(first_imp_inv, di, a, b, c, d))
341 {
342 if (! first_imp_inv)
343 {
344 first_imp_inv = g_malloc0 (sizeof*first_imp_inv);
345 this_ii = first_imp_inv;
346 }
347 else
348 {
349 this_ii -> next = g_malloc0 (sizeof*this_ii);
350 this_ii -> next -> prev = this_ii;
351 this_ii = this_ii -> next;
352 }
353 this_ii -> a = a;
354 this_ii -> b = b;
355 this_ii -> c = c;
356 this_ii -> d = d;
357 tmp_fprop = get_active_prop_using_atoms (dhii -> other, 4, ids);
358 if (tmp_fprop == NULL)
359 {
360 tmp_fprop = dhii -> def;
361 if (buf == NULL) show = FALSE;
362 }
363 else if (buf == NULL)
364 {
365 show = tmp_fprop -> show;
366 }
367 if (buf != NULL && tmp_fprop -> use)
368 {
369 if (di == 6)
370 {
371 stra = g_strdup_printf ("%4s\t%d\t%d\t%d\t%d",fkeysw[activef][di+2][tmp_fprop -> key], b+1, a+1, c+1, d+1);
372 }
373 else
374 {
375 stra = g_strdup_printf ("%4s\t%d\t%d\t%d\t%d",fkeysw[activef][di+2][tmp_fprop -> key], a+1, b+1, c+1, d+1);
376 }
377 print_info (stra, NULL, buf);
378 g_free (stra);
379 for (e=0; e<fvalues[activef][di+1][tmp_fprop -> key]; e++)
380 {
381 stra = g_strdup_printf ("\t%15.10f", tmp_fprop -> val[e]);
382 print_info (stra, NULL, buf);
383 g_free (stra);
384 if (e == 2)
385 {
386 // Print 1-4 electrostatic interaction scale factor
387 stra = g_strdup_printf ("\t%15.10f", 0.0);
388 print_info (stra, NULL, buf);
389 g_free (stra);
390 // Print 1-4 van der Waals interaction scale factor
391 stra = g_strdup_printf ("\t%15.10f", 0.0);
392 print_info (stra, NULL, buf);
393 g_free (stra);
394 }
395 }
396 print_info ("\n", NULL, buf);
397 }
398 else if (buf == NULL)
399 {
400 stra = g_strdup_printf ("%d", a+1);
401 strb = g_strdup_printf ("%d", b+1);
402 strc = g_strdup_printf ("%d", c+1);
403 strd = g_strdup_printf ("%d", d+1);
404 w = 0.0;
405 for (e=0; e<tmp_fmol -> multi; e++)
406 {
407 q = tmp_fmol -> atoms_id[a][e].a;
408 r = tmp_fmol -> atoms_id[a][e].b;
409 s = get_active_atom (tmp_fmol -> id, q) -> list[r];
410 q = tmp_fmol -> atoms_id[b][e].a;
411 r = tmp_fmol -> atoms_id[b][e].b;
412 t = get_active_atom (tmp_fmol -> id, q) -> list[r];
413 q = tmp_fmol -> atoms_id[c][e].a;
414 r = tmp_fmol -> atoms_id[c][e].b;
415 u = get_active_atom (tmp_fmol -> id, q) -> list[r];
416 q = tmp_fmol -> atoms_id[d][e].a;
417 r = tmp_fmol -> atoms_id[d][e].b;
418 v = get_active_atom (tmp_fmol -> id, q) -> list[r];
419 if (di == 6)
420 {
421 w += dihedral_3d (& tmp_proj -> cell, 0,
422 & tmp_proj -> atoms[0][t],
423 & tmp_proj -> atoms[0][u],
424 & tmp_proj -> atoms[0][s],
425 & tmp_proj -> atoms[0][v]).angle;
426 }
427 else
428 {
429 w += inversion_3d (& tmp_proj -> cell, 0,
430 & tmp_proj -> atoms[0][s],
431 & tmp_proj -> atoms[0][t],
432 & tmp_proj -> atoms[0][u],
433 & tmp_proj -> atoms[0][v]).angle;
434 }
435 }
436
437 w /= tmp_fmol -> multi;
438 stre = g_strdup_printf ("%.3f", w);
439 strf = g_strdup_printf ("%s (%s)", fnames[activef][di+2][tmp_fprop -> key], fkeysw[activef][di+2][tmp_fprop -> key]);
440 strg = parameters_info (di+1, tmp_fprop -> key, fvars_dihedral[activef][tmp_fprop -> key], tmp_fprop -> val);
441 gtk_tree_store_append (store, & di_level, iter);
442 gtk_tree_store_set (store, & di_level, 0, 0,
443 1, stra,
444 2, strb,
445 3, strc,
446 4, strd,
447 5, 0,
448 6, stre,
449 7, show,
450 8, tmp_fprop -> use,
451 9, strf,
452 10, strg,
453 11, dhii -> id, -1);
454 g_free (stra);
455 g_free (strb);
456 g_free (strc);
457 g_free (strd);
458 g_free (stre);
459 g_free (strf);
460 g_free (strg);
461 }
462 }
463 }
464 }
465 }
466 }
467 }
468 }
469 }
470 }
471 }
472 }
473 }
474 if (first_imp_inv)
475 {
476 this_ii = first_imp_inv;
477 while (this_ii -> next)
478 {
479 this_ii = this_ii -> next;
480 g_free (this_ii -> prev);
481 }
482 g_free (this_ii);
483 }
484 g_free (ids);
485}
486
499void print_dlp_dihedral (int dih, GtkTextBuffer * buf, field_struct * dh, int fi, GtkTreeStore * store, GtkTreeIter * iter)
500{
501 int a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s;
502 gboolean show;
503 gchar * stra, * strb, * strc, * strd, * stre, * strf, * strg;
504 float v;
505 GtkTreeIter di_level;
506 gboolean same_atom = FALSE;
507 gboolean * already_done;
508 if (tmp_fat -> id == tmp_fdt -> id && tmp_fbt -> id && tmp_fct -> id)
509 {
510 same_atom = TRUE;
511 already_done = allocbool (tmp_fmol -> mol -> natoms);
512 }
513 int * ids = allocint(4);
514 for (i=0; i<tmp_fat -> num; i++)
515 {
516 j = tmp_fat -> list[i];
517 if (tmp_proj -> atoms[0][j].coord[2] == fi)
518 {
519 a = ids[0] = tmp_fat -> list_id[i];
520 if (same_atom) already_done[a] = TRUE;
521 for (k=0; k<tmp_proj -> atoms[0][j].numv; k++)
522 {
523 l = tmp_proj -> atoms[0][j].vois[k];
524 if (tmp_proj -> atoms[0][l].faid == tmp_fbt -> id)
525 {
526 b = ids[1] = tmp_proj -> atoms[0][l].fid;
527 for (m=0; m<tmp_proj -> atoms[0][l].numv; m++)
528 {
529 n = tmp_proj -> atoms[0][l].vois[m];
530 if (n != j && tmp_proj -> atoms[0][n].faid == tmp_fct -> id)
531 {
532 c = ids[2] = tmp_proj -> atoms[0][n].fid;
533 for (o=0; o<tmp_proj -> atoms[0][n].numv; o++)
534 {
535 p = tmp_proj -> atoms[0][n].vois[o];
536 d = ids[3] = tmp_proj -> atoms[0][p].fid;
537 if (p != j && p != l && tmp_proj -> atoms[0][p].faid == tmp_fdt -> id && (! same_atom || (same_atom && ! already_done[d])))
538 {
539 tmp_fprop = get_active_prop_using_atoms (dh -> other, 4, ids);
540 if (tmp_fprop == NULL)
541 {
542 tmp_fprop = dh -> def;
543 if (buf == NULL) show = FALSE;
544 }
545 else if (buf == NULL)
546 {
547 show = tmp_fprop -> show;
548 }
549 if (buf != NULL && tmp_fprop -> use)
550 {
551 stra = g_strdup_printf ("%4s\t%d\t%d\t%d\t%d",fkeysw[activef][dih+2][tmp_fprop -> key], a+1, b+1, c+1, d+1);
552 print_info (stra, NULL, buf);
553 g_free (stra);
554 for (q=0; q<fvalues[activef][dih+1][tmp_fprop -> key]; q++)
555 {
556 stra = g_strdup_printf ("\t%15.10f", tmp_fprop -> val[q]);
557 print_info (stra, NULL, buf);
558 g_free (stra);
559 if (q == 2)
560 {
561 // Print 1-4 electrostatic interaction scale factor
562 stra = g_strdup_printf ("\t%15.10f", 0.0);
563 print_info (stra, NULL, buf);
564 g_free (stra);
565 // Print 1-4 van der Waals interaction scale factor
566 stra = g_strdup_printf ("\t%15.10f", 0.0);
567 print_info (stra, NULL, buf);
568 g_free (stra);
569 }
570 }
571 print_info ("\n", NULL, buf);
572 }
573 else if (buf == NULL)
574 {
575 stra = g_strdup_printf ("%d", a+1);
576 strb = g_strdup_printf ("%d", b+1);
577 strc = g_strdup_printf ("%d", c+1);
578 strd = g_strdup_printf ("%d", d+1);
579 v = 0.0;
580 for (q=0; q<tmp_fmol -> multi; q++)
581 {
582 r = tmp_fmol -> atoms_id[a][q].a;
583 s = tmp_fmol -> atoms_id[a][q].b;
584 e = get_active_atom (tmp_fmol -> id, r) -> list[s];
585 r = tmp_fmol -> atoms_id[b][q].a;
586 s = tmp_fmol -> atoms_id[b][q].b;
587 f = get_active_atom (tmp_fmol -> id, r) -> list[s];
588 r = tmp_fmol -> atoms_id[c][q].a;
589 s = tmp_fmol -> atoms_id[c][q].b;
590 g = get_active_atom (tmp_fmol -> id, r) -> list[s];
591 r = tmp_fmol -> atoms_id[d][q].a;
592 s = tmp_fmol -> atoms_id[d][q].b;
593 h = get_active_atom (tmp_fmol -> id, r) -> list[s];
594 v += dihedral_3d (& tmp_proj -> cell, 0,
595 & tmp_proj -> atoms[0][e],
596 & tmp_proj -> atoms[0][f],
597 & tmp_proj -> atoms[0][g],
598 & tmp_proj -> atoms[0][h]).angle;
599 }
600 v /= tmp_fmol -> multi;
601 stre = g_strdup_printf ("%.3f", v);
602 strf = g_strdup_printf ("%s (%s)", fnames[activef][dih+2][tmp_fprop -> key], fkeysw[activef][dih+2][tmp_fprop -> key]);
603 strg = parameters_info (dih+1, tmp_fprop -> key, fvars_dihedral[activef][tmp_fprop -> key], tmp_fprop -> val);
604 gtk_tree_store_append (store, & di_level, iter);
605 gtk_tree_store_set (store, & di_level, 0, 0,
606 1, stra,
607 2, strb,
608 3, strc,
609 4, strd,
610 5, 0,
611 6, stre,
612 7, show,
613 8, tmp_fprop -> use,
614 9, strf,
615 10, strg,
616 11, dh -> id, -1);
617 g_free (stra);
618 g_free (strb);
619 g_free (strc);
620 g_free (strd);
621 g_free (stre);
622 g_free (strf);
623 g_free (strg);
624 }
625 }
626 }
627 }
628 }
629 }
630 }
631 }
632 }
633 g_free (ids);
634 if (same_atom) g_free (already_done);
635}
636
649void print_dlp_angle (int ai, GtkTextBuffer * buf, field_struct * an, int fi, GtkTreeStore * store, GtkTreeIter * iter)
650{
651 int a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, u;
652 gboolean show;
653 gchar * stra, * strb, * strc, * strd, * stre, * strf;
654 float v;
655 GtkTreeIter an_level;
656 int * ids = allocint(3);
657 gboolean same_atom = FALSE;
658 gboolean * already_done;
659 if (tmp_fat -> id == tmp_fct -> id)
660 {
661 same_atom = TRUE;
662 already_done = allocbool (tmp_fmol -> mol -> natoms);
663 }
664
665 for (i=0; i<tmp_fat -> num; i++)
666 {
667 j = tmp_fat -> list[i];
668 if (tmp_proj -> atoms[0][j].coord[2] == fi)
669 {
670 k = ids[0] = tmp_proj -> atoms[0][j].fid;
671 if (same_atom) already_done[k] = TRUE;
672 for (l=0; l<tmp_proj -> atoms[0][j].numv; l++)
673 {
674 m = tmp_proj -> atoms[0][j].vois[l];
675 if (tmp_proj -> atoms[0][m].faid == tmp_fbt -> id)
676 {
677 n = ids[1] = tmp_proj -> atoms[0][m].fid;
678 for (o=0; o<tmp_proj -> atoms[0][m].numv; o++)
679 {
680 p = tmp_proj -> atoms[0][m].vois[o];
681 q = ids[2] = tmp_proj -> atoms[0][p].fid;
682 if (p != j && tmp_proj -> atoms[0][p].faid == tmp_fct -> id && (! same_atom || (same_atom && ! already_done[q])))
683 {
684 tmp_fprop = get_active_prop_using_atoms (an -> other, 3, ids);
685 if (tmp_fprop == NULL)
686 {
687 tmp_fprop = an -> def;
688 if (buf == NULL) show = FALSE;
689 }
690 else if (buf == NULL)
691 {
692 show = tmp_fprop -> show;
693 }
694 if (buf != NULL && tmp_fprop -> use)
695 {
696 stra = g_strdup_printf ("%4s\t%d\t%d\t%d",fkeysw[activef][ai+2][tmp_fprop -> key], k+1, n+1, q+1);
697 print_info (stra, NULL, buf);
698 g_free (stra);
699 for (u=0; u<fvalues[activef][ai+1][tmp_fprop -> key]; u++)
700 {
701 stra = g_strdup_printf ("\t%15.10f", tmp_fprop -> val[u]);
702 print_info (stra, NULL, buf);
703 g_free (stra);
704 }
705 print_info ("\n", NULL, buf);
706 }
707 else if (buf == NULL)
708 {
709 stra = g_strdup_printf ("%d", k+1);
710 strb = g_strdup_printf ("%d", n+1);
711 strc = g_strdup_printf ("%d", q+1);
712 v = 0.0;
713 for (u=0; u<tmp_fmol -> multi; u++)
714 {
715 a = tmp_fmol -> atoms_id[k][u].a;
716 b = tmp_fmol -> atoms_id[k][u].b;
717 c = get_active_atom (tmp_fmol -> id, a) -> list[b];
718 d = tmp_fmol -> atoms_id[n][u].a;
719 e = tmp_fmol -> atoms_id[n][u].b;
720 f = get_active_atom (tmp_fmol -> id, d) -> list[e];
721 e = tmp_fmol -> atoms_id[q][u].a;
722 g = tmp_fmol -> atoms_id[q][u].b;
723 h = get_active_atom (tmp_fmol -> id, e) -> list[g];
724 v += angle_3d (& tmp_proj -> cell, 0,
725 & tmp_proj -> atoms[0][c],
726 & tmp_proj -> atoms[0][f],
727 & tmp_proj -> atoms[0][h]).angle;
728 }
729 v /= tmp_fmol -> multi;
730 strd = g_strdup_printf ("%.3f", v);
731 stre = g_strdup_printf ("%s (%s)", fnames[activef][ai+2][tmp_fprop -> key], fkeysw[activef][ai+2][tmp_fprop -> key]);
732 strf = parameters_info (ai+1, tmp_fprop -> key, fvars_angle[activef][tmp_fprop -> key], tmp_fprop -> val);
733 gtk_tree_store_append (store, & an_level, iter);
734 gtk_tree_store_set (store, & an_level, 0, 0,
735 1, stra,
736 2, strb,
737 3, strc,
738 4, 0,
739 5, strd,
740 6, show,
741 7, tmp_fprop -> use,
742 8, stre,
743 9, strf,
744 10, an -> id, -1);
745 g_free (stra);
746 g_free (strb);
747 g_free (strc);
748 g_free (strd);
749 g_free (stre);
750 g_free (strf);
751 }
752 }
753 }
754 }
755 }
756 }
757 }
758 g_free (ids);
759 if (same_atom) g_free (already_done);
760}
761
774void print_dlp_bond (int bi, GtkTextBuffer * buf, field_struct * bd, int fi, GtkTreeStore * store, GtkTreeIter * iter)
775{
776 int i, j, k, l, m, n, o, p, q, r, s, t, u;
777 gboolean show;
778 gchar * stra, * strb, * strc, * strd, * stre;
779 float v;
780 int * ids = allocint(2);
781 GtkTreeIter bd_level;
782 gboolean same_atom = FALSE;
783 gboolean * already_done;
784 if (tmp_fat -> id == tmp_fbt -> id)
785 {
786 same_atom = TRUE;
787 already_done = allocbool (tmp_fmol -> mol -> natoms);
788 }
789 for (i=0; i<tmp_fat -> num; i++)
790 {
791 j = tmp_fat -> list[i];
792 if (tmp_proj -> atoms[0][j].coord[2] == fi)
793 {
794 k = ids[0] = tmp_fat -> list_id[i];
795 if (same_atom) already_done[k] = TRUE;
796 for (l=0; l<tmp_proj -> atoms[0][j].numv; l++)
797 {
798 m = tmp_proj -> atoms[0][j].vois[l];
799 n = ids[1] = tmp_proj -> atoms[0][m].fid;
800 if (tmp_proj -> atoms[0][m].faid == tmp_fbt -> id && (! same_atom || (same_atom && ! already_done[n])))
801 {
802 tmp_fprop = get_active_prop_using_atoms (bd -> other, 2, ids);
803 if (tmp_fprop == NULL)
804 {
805 tmp_fprop = bd -> def;
806 if (buf == NULL) show = FALSE;
807 }
808 else if (buf == NULL)
809 {
810 show = tmp_fprop -> show;
811 }
812 if (buf != NULL && tmp_fprop -> use)
813 {
814 stra = g_strdup_printf ("%4s\t%d\t%d",fkeysw[activef][bi+2][tmp_fprop -> key], k+1, n+1);
815 print_info (stra, NULL, buf);
816 g_free (stra);
817 for (o=0; o<fvalues[activef][bi+1][tmp_fprop -> key]; o++)
818 {
819 stra = g_strdup_printf ("\t%15.10f", tmp_fprop -> val[o]);
820 print_info (stra, NULL, buf);
821 g_free (stra);
822 }
823 print_info ("\n", NULL, buf);
824 }
825 else if (buf == NULL)
826 {
827 stra = g_strdup_printf ("%d", k+1);
828 strb = g_strdup_printf ("%d", n+1);
829 v = 0.0;
830 for (o=0; o<tmp_fmol -> multi; o++)
831 {
832 p = tmp_fmol -> atoms_id[k][o].a;
833 q = tmp_fmol -> atoms_id[k][o].b;
834 r = get_active_atom (tmp_fmol -> id, p) -> list[q];
835 s = tmp_fmol -> atoms_id[n][o].a;
836 t = tmp_fmol -> atoms_id[n][o].b;
837 u = get_active_atom (tmp_fmol -> id, s) -> list[t];
838 v += distance_3d (& tmp_proj -> cell, 0, & tmp_proj -> atoms[0][r], & tmp_proj -> atoms[0][u]).length;
839 }
840 v /= tmp_fmol -> multi;
841 strc = g_strdup_printf ("%.3f", v);
842 strd = g_strdup_printf ("%s (%s)", fnames[activef][bi+2][tmp_fprop -> key], fkeysw[activef][bi+2][tmp_fprop -> key]);
843 stre = parameters_info (bi+1, tmp_fprop -> key, fvars_bond[activef][tmp_fprop -> key], tmp_fprop -> val);
844 gtk_tree_store_append (store, & bd_level, iter);
845 gtk_tree_store_set (store, & bd_level, 0, 0,
846 1, stra,
847 2, strb,
848 3, 0,
849 4, strc,
850 5, show,
851 6, tmp_fprop -> use,
852 7, strd,
853 8, stre,
854 9, bd -> id, -1);
855 g_free (stra);
856 g_free (strb);
857 g_free (strc);
858 g_free (strd);
859 g_free (stre);
860 }
861 }
862 }
863 }
864 }
865 g_free (ids);
866 if (same_atom) g_free (already_done);
867}
868
877void print_dlp_rigid (GtkTextBuffer * buf, field_rigid * rig)
878{
879 gchar * str;
880 int h, i, j, k, l, m, n;
881 str = g_strdup_printf ("%d\t", rig -> num);
882 j = 1;
883 n = rig -> num;
884 if (rig -> num > 15)
885 {
886 k = rig -> num - 15;
887 l = k / 16;
888 j += l + 1;
889 n = k - l*16;
890 }
891 h = 0;
892 for (i=0; i<j; i++)
893 {
894 k = (i) ? 1 : 0;
895 l = (j == 1 || (j > 1 && i == j-1)) ? n : 15+k;
896 for (m=0; m<l; m++)
897 {
898 str = g_strdup_printf ("%s\t%d", str, rig -> list[h]+1);
899 h ++;
900 }
901 str = g_strdup_printf ("%s\n", str);
902 }
903 print_info (str, NULL, buf);
904 g_free (str);
905}
906
915void print_dlp_tet (GtkTextBuffer * buf, field_tethered * tet)
916{
917 gchar * str;
918 str = g_strdup_printf ("%4s\t\%d", fkeysw[activef][1][tet -> key], tet -> num);
919 print_info (str, NULL, buf);
920 g_free (str);
921 int i;
922 for (i=0; i<fvalues[activef][0][tmp_ftet -> key]; i++)
923 {
924 str = g_strdup_printf ("\t%15.10f", tmp_ftet -> val[i]);
925 print_info (str, NULL, buf);
926 g_free (str);
927 }
928 print_info ("\n", NULL, buf);
929}
930
939void print_dlp_pmf (GtkTextBuffer * buf, field_pmf * pmf)
940{
941 gchar * str;
942 int i, j;
943 print_info ("PMF", "bold", buf);
944 str = g_strdup_printf ("\t%f\n", pmf -> length);
945 print_info (str, NULL, buf);
946 g_free (str);
947 for (i=0; i<2; i++)
948 {
949 str = g_strdup_printf ("PMF UNIT %d\n", pmf -> num[i]);
950 print_info (str, NULL, buf);
951 g_free (str);
952 for (j=0; j < pmf -> num[i]; j++)
953 {
954 str = g_strdup_printf ("%d\t%f\n", pmf -> list[i][j]+1, pmf -> weight[i][j]);
955 print_info (str, NULL, buf);
956 g_free (str);
957 }
958 }
959}
960
969void print_dlp_cons (GtkTextBuffer * buf, field_constraint * cons)
970{
971 gchar * str;
972 str = g_strdup_printf ("%d\t\%d\t%f\n", cons -> ia[0], cons -> ia[1], cons -> length);
973 print_info (str, NULL, buf);
974 g_free (str);
975}
976
986void print_dlp_shell (GtkTextBuffer * buf, field_molecule * fmol, field_shell * shell)
987{
988 gchar * str;
989 str = g_strdup_printf ("%d\t\%d\t%f\t%f\n", shell -> ia[0], shell -> ia[1], shell -> k2, shell -> k4);
990 print_info (str, NULL, buf);
991 g_free (str);
992}
993
1003void print_dlp_atom (GtkTextBuffer * buf, int at, int numat)
1004{
1005 gchar * str;
1006 if (tmp_fat -> frozen_id[at])
1007 {
1008 str = g_strdup_printf ("%8s %15.10f %15.10f %d %d\n", tmp_fat -> name, tmp_fat -> mass, tmp_fat -> charge, numat, 1);
1009 }
1010 else
1011 {
1012 str = g_strdup_printf ("%8s %15.10f %15.10f %d\n", tmp_fat -> name, tmp_fat -> mass, tmp_fat -> charge, numat);
1013 }
1014 print_info (str, NULL, buf);
1015 g_free (str);
1016}
1017
1027{
1028 int i = 0;
1029 tmp_fstr = fmol -> first_struct[sid];
1030 while (tmp_fstr)
1031 {
1032 if (tmp_fstr -> def -> use)
1033 {
1034 i += tmp_fstr -> num;
1035 }
1036 else if (tmp_fstr -> other)
1037 {
1038 tmp_fprop = tmp_fstr -> other;
1039 while (tmp_fprop)
1040 {
1041 if (tmp_fprop -> use) i ++;
1042 tmp_fprop = tmp_fprop -> next;
1043 }
1044 }
1045 tmp_fstr = tmp_fstr -> next;
1046 }
1047 return i;
1048}
1049
1058void print_dlp_molecule (GtkTextBuffer * buf, field_molecule * fmol)
1059{
1060 gchar * str;
1061 str = g_strdup_printf ("%s", fmol -> name);
1062 print_info (str, "bold_orange", buf);
1063 g_free (str);
1064 print_info ("\nNUMMOLS\t", "bold", buf);
1065 str = g_strdup_printf ("%d", fmol -> multi);
1066 print_info (str, "bold_green", buf);
1067 g_free (str);
1068 int i, j, k, l, m, n, o, p;
1069
1070 j = 0;
1071 tmp_fat = fmol -> first_atom;
1072 for (i=0; i<fmol -> atoms; i++)
1073 {
1074 j += tmp_fat -> num;
1075 if (tmp_fat -> next != NULL) tmp_fat = tmp_fat -> next;
1076 }
1077 j /= fmol -> multi;
1078 if (j != fmol -> mol -> natoms) g_debug ("PRINT:: Error the number of atom(s) is wrong ?!");
1079 print_info ("\nATOMS\t", "bold", buf);
1080 str = g_strdup_printf ("%d\n", fmol -> mol -> natoms);
1081 print_info (str, "bold_blue", buf);
1082 g_free (str);
1083 for (i=0; i < fmol -> mol -> natoms ; i+=(m-i))
1084 {
1085 j = fmol -> atoms_id[i][0].a;
1086 tmp_fat = get_active_atom (fmol -> id, j);
1087 k = fmol -> atoms_id[i][0].b;
1088 l = tmp_fat -> frozen_id[k];
1089 for (m=i+1; m<fmol -> mol -> natoms; m++)
1090 {
1091 n = fmol -> atoms_id[m][0].a;
1092 o = fmol -> atoms_id[m][0].b;
1093 tmp_fbt = get_active_atom (fmol -> id, n);
1094 p = tmp_fbt -> frozen_id[o];
1095 if (j != n || l != p) break;
1096 }
1097 print_dlp_atom (buf, k, m-i);
1098 }
1099 // Shells
1100 int ncs = 0;
1101 if (tmp_field -> afp[10])
1102 {
1103 tmp_fshell = fmol -> first_shell;
1104 while (tmp_fshell)
1105 {
1106 if (tmp_fshell -> use)
1107 {
1108 if (tmp_fshell -> ia[0] && tmp_fshell -> ia[1]) ncs ++;
1109 }
1110 tmp_fshell = tmp_fshell -> next;
1111 }
1112 }
1113 if (ncs)
1114 {
1115 print_info ("SHELLS\t", "bold", buf);
1116 str = g_strdup_printf ("%d\n", ncs);
1117 print_info (str, "bold", buf);
1118 g_free (str);
1119 tmp_fshell = fmol -> first_shell;
1120 while (tmp_fshell)
1121 {
1122 if (tmp_fshell -> use && tmp_fshell -> ia[0] && tmp_fshell -> ia[1])
1123 {
1124 print_dlp_shell (buf, fmol, tmp_fshell);
1125 }
1126 tmp_fshell = tmp_fshell -> next;
1127 }
1128 }
1129
1130 // Constraints
1131 if (tmp_field -> afp[11])
1132 {
1133 j = 0;
1134 tmp_fcons = fmol -> first_constraint;
1135 while (tmp_fcons)
1136 {
1137 if (tmp_fcons -> use) j ++;
1138 tmp_fcons = tmp_fcons -> next;
1139 }
1140 if (j > 0)
1141 {
1142 print_info ("CONSTRAINTS\t", "bold", buf);
1143 str = g_strdup_printf ("%d\n", j);
1144 print_info (str, "bold", buf);
1145 g_free (str);
1146 tmp_fcons = fmol -> first_constraint;
1147 while (tmp_fcons)
1148 {
1149 if (tmp_fcons -> use) print_dlp_cons (buf, tmp_fcons);
1150 tmp_fcons = tmp_fcons -> next;
1151 }
1152 }
1153 }
1154
1155 // PMFs
1156 if (tmp_field -> afp[12])
1157 {
1158 tmp_fpmf = fmol -> first_pmf;
1159 while (tmp_fpmf)
1160 {
1161 if (tmp_fpmf -> use)
1162 {
1163 print_dlp_pmf (buf, tmp_fpmf);
1164 break;
1165 }
1166 tmp_fpmf = tmp_fpmf -> next;
1167 }
1168 }
1169
1170 // Rigid
1171 if (tmp_field -> afp[13])
1172 {
1173 j = 0;
1174 k = 0;
1175 tmp_frig = fmol -> first_rigid;
1176 while (tmp_frig)
1177 {
1178 if (tmp_frig -> use) j ++;
1179 tmp_frig = tmp_frig -> next;
1180 }
1181 if (j > 0)
1182 {
1183 print_info ("RIGID\t", "bold", buf);
1184 str = g_strdup_printf ("%d\n", j);
1185 print_info (str, "bold", buf);
1186 g_free (str);
1187 tmp_frig = fmol -> first_rigid;
1188 while (tmp_frig)
1189 {
1190 if (tmp_frig -> use)
1191 {
1193 }
1194 tmp_frig = tmp_frig -> next;
1195 }
1196 }
1197 }
1198
1199 // Tethering
1200 if (tmp_field -> afp[14])
1201 {
1202 j = 0;
1203 tmp_ftet = fmol -> first_tethered;
1204 while (tmp_ftet)
1205 {
1206 if (tmp_ftet -> use && tmp_ftet -> num) j ++;
1207 tmp_ftet = tmp_ftet -> next;
1208 }
1209 if (j > 0)
1210 {
1211 print_info ("TETH\t", "bold", buf);
1212 str = g_strdup_printf ("%d\n", j);
1213 print_info (str, "bold", buf);
1214 g_free (str);
1215 tmp_ftet = fmol -> first_tethered;
1216 while (tmp_ftet)
1217 {
1218 if (tmp_ftet -> use) print_dlp_tet (buf, tmp_ftet);
1219 tmp_ftet = tmp_ftet -> next;
1220 }
1221 }
1222 }
1223 gchar * str_title[8] = {"BONDS ", "BONDS ", "ANGLES ", "ANGLES ", "DIHEDRALS ", "DIHEDRALS ", "DIHEDRALS ", "INVERSIONS "};
1224 gboolean doprint;
1225 for (i=0; i<8; i++)
1226 {
1227 if (tmp_field -> afp[i+15])
1228 {
1229 j = get_num_struct_to_print (fmol, i);
1230 if ((i == 0 || i == 2 || i == 4) && tmp_field -> afp[i+16])
1231 {
1232 // To add the number of constraints
1233 j += get_num_struct_to_print (fmol, i+1);
1234 }
1235 if ((i == 4 && tmp_field -> afp[i+17]) || (i == 5 && ! tmp_field -> afp[i+14] && tmp_field -> afp[i+16]))
1236 {
1237 // To add the number of impropers
1238 k = (i == 4) ? 2 : 1;
1239 j += get_num_struct_to_print (fmol, i+k);
1240 }
1241 doprint = (i == 0 || i == 2 || i == 4 || i == 7) ? TRUE : FALSE;
1242 if ((i == 1 || i == 3 || i == 5) && ! tmp_field -> afp[i+14])
1243 {
1244 doprint = TRUE;
1245 }
1246 else if (i == 5 && ! tmp_field -> afp[i+14])
1247 {
1248 doprint = TRUE;
1249 }
1250 else if (i == 6 && ! tmp_field -> afp[i+13] && ! tmp_field -> afp[i+14])
1251 {
1252 doprint = TRUE;
1253 }
1254
1255 if (j > 0)
1256 {
1257 if (doprint)
1258 {
1259 print_info (str_title[i], "bold", buf);
1260 str = g_strdup_printf ("%d\n", j);
1261 print_info (str, "bold_blue", buf);
1262 g_free (str);
1263 }
1264 tmp_fstr = fmol -> first_struct[i];
1265 while (tmp_fstr)
1266 {
1267 tmp_fat = get_active_atom (fmol -> id, tmp_fstr -> aid[0]);
1268 tmp_fbt = get_active_atom (fmol -> id, tmp_fstr -> aid[1]);
1269 if (i > 1) tmp_fct = get_active_atom (fmol -> id, tmp_fstr -> aid[2]);
1270 if (i > 3) tmp_fdt = get_active_atom (fmol -> id, tmp_fstr -> aid[3]);
1271 if (i == 0 || i == 1) print_dlp_bond (i, buf, tmp_fstr, fmol -> fragments[0], NULL, NULL);
1272 if (i == 2 || i == 3) print_dlp_angle (i, buf, tmp_fstr, fmol -> fragments[0], NULL, NULL);
1273 if (i == 4 || i == 5) print_dlp_dihedral (i, buf, tmp_fstr, fmol -> fragments[0], NULL, NULL);
1274 if (i == 6 || i == 7) print_dlp_improper_inversion (i, buf, tmp_fstr, fmol -> fragments[0], NULL, NULL);
1275 tmp_fstr = tmp_fstr -> next;
1276 }
1277 }
1278 }
1279 }
1280
1281 print_info ("FINISH\n", "bold_orange", buf);
1282}
1283
1292void print_dlp_body (GtkTextBuffer * buf, field_nth_body * body)
1293{
1294 gchar * str;
1295 int i, j;
1296 j = body_at (body -> bd);
1297 if (! body -> bd)
1298 {
1299 for (i=0; i<j; i++) print_info (g_strdup_printf ("%8s\t", get_body_element_name (body, i, 0)), NULL, buf);
1300 }
1301 else
1302 {
1303 for (i=0; i<j; i++) print_info (g_strdup_printf ("%8s\t", get_active_atom(body -> ma[i][0], body -> a[i][0]) -> name), NULL, buf);
1304 }
1305 str = g_strdup_printf ("%4s",fkeysw[activef][10+ body -> bd][body -> key]);
1306 print_info (str, NULL, buf);
1307 g_free (str);
1308 for (i=0; i<fvalues[activef][9+ body -> bd][body -> key]; i++)
1309 {
1310 str = g_strdup_printf ("\t%15.10f", body -> val[i]);
1311 print_info (str, NULL, buf);
1312 g_free (str);
1313 }
1314 print_info ("\n", NULL, buf);
1315}
1316
1326void print_dlp_tersoff_cross (GtkTextBuffer * buf, field_nth_body * body_a, field_nth_body * body_b)
1327{
1328 gchar * str;
1329 int j;
1330 print_info (g_strdup_printf ("%8s\t", get_active_atom(body_a -> ma[0][0], body_a -> a[0][0]) -> name), NULL, buf);
1331 print_info (g_strdup_printf ("%8s\t", get_active_atom(body_b -> ma[0][0], body_b -> a[0][0]) -> name), NULL, buf);
1332 for (j=0; j<3; j++)
1333 {
1334
1335 str = g_strdup_printf ("%15.10f", tmp_field -> cross[body_a -> id][body_b -> id][j]);
1336 print_info (str, NULL, buf);
1337 g_free (str);
1338 if (j<2) print_info ("\t", NULL, buf);
1339 }
1340 print_info ("\n", NULL, buf);
1341}
1342
1351void print_dlp_tersoff (GtkTextBuffer * buf, field_nth_body * body)
1352{
1353 gchar * str;
1354 int i, j, k;
1355 int num[2]={2, 3};
1356 int nc[2][3]={{5, 6, 3}, {5, 6, 5}};
1357 j = body_at (body -> bd);
1358 k = 0;
1359 for (i=0; i<num[body -> key]; i++)
1360 {
1361 if (i==0)
1362 {
1363 print_info (g_strdup_printf ("%8s\t", get_active_atom(body -> ma[0][0], body -> a[0][0]) -> name), NULL, buf);
1364 str = g_strdup_printf ("%4s\t",fkeysw[activef][10+body -> bd][body -> key]);
1365 print_info (str, NULL, buf);
1366 g_free (str);
1367 }
1368 else
1369 {
1370 print_info (" \t \t", NULL, buf);
1371 }
1372 for (j=0; j<nc[body -> key][i]; j++)
1373 {
1374 if (j > 0) print_info ("\t", NULL, buf);
1375 str = g_strdup_printf ("%15.10f", body -> val[j+k]);
1376 print_info (str, NULL, buf);
1377 g_free (str);
1378 }
1379 print_info ("\n", NULL, buf);
1380 k += nc[body -> key][i];
1381 }
1382 if (! body -> key)
1383 {
1384 field_nth_body * tmp_fbo;
1385 tmp_fbo = tmp_field -> first_body[2];
1386 while (tmp_fbo)
1387 {
1388 print_dlp_tersoff_cross (buf, body, tmp_fbo);
1389 tmp_fbo = tmp_fbo -> next;
1390 }
1391 }
1392}
1393
1401void print_dlp_field (GtkTextBuffer * buf)
1402{
1403 int i, j;
1404 gchar * str;
1405
1406 GtkTextIter bStart;
1407 GtkTextIter bEnd;
1408
1409 gtk_text_buffer_get_start_iter (buf, & bStart);
1410 gtk_text_buffer_get_end_iter (buf, & bEnd);
1411 gtk_text_buffer_delete (buf, & bStart, & bEnd);
1412
1413 str = g_strdup_printf ("# This file was created using %s\n", PACKAGE);
1414 print_info (str, NULL, buf);
1415 g_free (str);
1416 str = g_strdup_printf ("# %s contains:\n", prepare_for_title(tmp_proj -> name));
1417 print_info (str, NULL, buf);
1418 g_free (str);
1419 i = 0;
1420 for (j=0; j<tmp_proj -> modelfc -> mol_by_step[0]; j++)
1421 {
1422 i += tmp_proj -> modelfc -> mols[0][j].multiplicity;
1423 }
1424 str = g_strdup_printf ("# - %d atoms\n"
1425 "# - %d isolated molecular fragments\n"
1426 "# - %d distinct molecules\n",
1427 tmp_proj -> natomes, i, tmp_proj -> modelfc -> mol_by_step[0]);
1428 print_info (str, NULL, buf);
1429 g_free (str);
1430
1431 print_info ("# Energy unit:\n", NULL, buf);
1432 print_info ("UNITS ", "bold", buf);
1433 str = g_strdup_printf ("%s\n", fkeysw[activef][0][tmp_field -> energy_unit]);
1434 print_info (str, "bold_green", buf);
1435 g_free (str);
1436 print_info ("# Number of field molecules:\n", NULL, buf);
1437 print_info ("MOLECULES ", "bold", buf);
1438 str = g_strdup_printf ("%d\n", tmp_field -> molecules);
1439 print_info (str, "bold_red", buf);
1440 g_free (str);
1441 tmp_fmol = tmp_field -> first_molecule;
1442 for (i=0; i<tmp_field -> molecules; i++)
1443 {
1444 str = g_strdup_printf ("# Begin molecule %d\n", i+1);
1445 print_info (str, NULL, buf);
1446 g_free (str);
1448 str = g_strdup_printf ("# End molecule %d\n", i+1);
1449 print_info (str, NULL, buf);
1450 g_free (str);
1451 if (tmp_fmol -> next != NULL) tmp_fmol = tmp_fmol -> next;
1452 }
1453
1454 // Non bonded
1455 gchar * nd_title[5] = {"VDW", "METAL", "TERSOFF", "TBP", "FBP"};
1456 gchar * com_ndb[5] = {"Van der Walls pair", "metal", "Tersoff", "three-body", "four-body"};
1457 for (i=0; i<5; i++)
1458 {
1459 if (tmp_field -> afp[i+23])
1460 {
1461 j=0;
1462 tmp_fbody = tmp_field -> first_body[i];
1463 while (tmp_fbody)
1464 {
1465 if (tmp_fbody -> use) j++;
1466 tmp_fbody = tmp_fbody -> next;
1467 }
1468 if (j > 0)
1469 {
1470 str = g_strdup_printf ("# Non-bonded: %s potential(s)\n", com_ndb[i]);
1471 print_info (str, NULL, buf);
1472 g_free (str);
1473 print_info (nd_title[i], "bold", buf);
1474 str = g_strdup_printf (" %d\n", j);
1475 print_info (str, "bold_red", buf);
1476 tmp_fbody = tmp_field -> first_body[i];
1477 while (tmp_fbody)
1478 {
1479 if (tmp_fbody -> use)
1480 {
1481 if (i == 2)
1482 {
1484 }
1485 else
1486 {
1488 }
1489 }
1490 tmp_fbody = tmp_fbody -> next;
1491 }
1492 }
1493 }
1494 }
1495
1496 // External fields
1497 if (tmp_field -> afp[28])
1498 {
1499 i = 0;
1500 tmp_fext = tmp_field -> first_external;
1501 while (tmp_fext)
1502 {
1503 if (tmp_fext -> use) i ++;
1504 tmp_fext = tmp_fext -> next;
1505 }
1506 if (i == 1)
1507 {
1508 print_info ("EXTERN", "bold", buf);
1509 tmp_fext = tmp_field -> first_external;
1510 while (tmp_fext)
1511 {
1512 if (tmp_fext -> use)
1513 {
1514 str = g_strdup_printf ("\n%4s",fkeysw[activef][15][tmp_fext -> key]);
1515 print_info (str, NULL, buf);
1516 g_free (str);
1517 for (j=0; j<fvalues[activef][SEXTERN-6][tmp_fext -> key]; j++)
1518 {
1519 print_info (g_strdup_printf ("\t%15.10f", tmp_fext -> val[j]), NULL, buf);
1520 }
1521 print_info ("\n", NULL, buf);
1522 break;
1523 }
1524 tmp_fext = tmp_fext -> next;
1525 }
1526 }
1527 }
1528 print_info ("CLOSE", "bold", buf);
1529}
1530
1537{
1538 box_info * box = & tmp_proj -> cell.box[0];
1539 if (box -> param[1][0] == box -> param[1][1] && box -> param[1][0] == box -> param[1][2] && box -> param[1][0] == 90.0)
1540 {
1541 if (box -> param[0][0] == box -> param[0][1] && box -> param[0][0] == box -> param[0][2])
1542 {
1543 return 1;
1544 }
1545 else
1546 {
1547 return 2;
1548 }
1549 }
1550 else if (box -> vect[0][1] == 0.0 && box -> vect[0][2] == 0.0 && box -> vect[1][0] == 0.0
1551 && box -> vect[1][2] == 0.0 && box -> vect[2][0] == 0.0 && box -> vect[2][1] == 0.0)
1552 {
1553 if (box -> vect[0][0] == box -> vect[1][1] && box -> vect[0][0] == box -> vect[2][2])
1554 {
1555 return 1;
1556 }
1557 else
1558 {
1559 return 2;
1560 }
1561 }
1562 else
1563 {
1564 return 3;
1565 }
1566}
1567
1575void print_dlp_config (GtkTextBuffer * buf)
1576{
1577 int h, i, j, k, l, m, n;
1578 int pbc;
1579 gchar * str;
1580
1581 GtkTextIter bStart;
1582 GtkTextIter bEnd;
1583
1584 gtk_text_buffer_get_start_iter (buf, & bStart);
1585 gtk_text_buffer_get_end_iter (buf, & bEnd);
1586 gtk_text_buffer_delete (buf, & bStart, & bEnd);
1587
1588 str = g_strdup_printf ("# DL-POLY CONFIG file created by %s, %s - %d atoms\n",
1589 PACKAGE,
1590 prepare_for_title(tmp_proj -> name),
1591 tmp_proj -> natomes);
1592 print_info (str, "bold", buf);
1593 g_free (str);
1594 if (tmp_proj -> cell.pbc)
1595 {
1596 pbc = get_pbc ();
1597 }
1598 else
1599 {
1600 pbc = 0;
1601 }
1602 str = g_strdup_printf ("%d", 0);
1603 print_info (str, "bold_red", buf);
1604 g_free (str);
1605 str = g_strdup_printf ("\t%d", pbc);
1606 print_info (str, "bold_green", buf);
1607 g_free (str);
1608 str = g_strdup_printf ("\t%d\n", tmp_proj -> natomes);
1609 print_info (str, "bold_blue", buf);
1610 g_free (str);
1611 if (pbc > 0)
1612 {
1613 for (i=0; i<3; i++)
1614 {
1615 str = g_strdup_printf ("%f\t%f\t%f\n",
1616 tmp_proj -> cell.box[0].vect[i][0],
1617 tmp_proj -> cell.box[0].vect[i][1],
1618 tmp_proj -> cell.box[0].vect[i][2]);
1619 print_info (str, NULL, buf);
1620 g_free (str);
1621
1622 }
1623 }
1624 tmp_fmol = tmp_field -> first_molecule;
1625 h = 0;
1626 for (i=0; i<tmp_field -> molecules; i++)
1627 {
1628 for (j=0; j<tmp_fmol -> multi; j++)
1629 {
1630 for (k=0; k<tmp_fmol -> mol -> natoms; k++)
1631 {
1632 l = tmp_fmol -> atoms_id[k][j].a;
1633 m = tmp_fmol -> atoms_id[k][j].b;
1634 tmp_fat = get_active_atom (tmp_fmol -> id, l);
1635 str = g_strdup_printf ("%8s", tmp_fat -> name);
1636 print_info (str, "bold", buf);
1637 g_free (str);
1638 if (tmp_field -> sys_opts[2])
1639 {
1640 print_info ("\n", NULL, buf);
1641 }
1642 else
1643 {
1644 str = g_strdup_printf (" %d\n", h+1);
1645 print_info (str, "bold_red", buf);
1646 g_free (str);
1647 h ++;
1648 }
1649 n = tmp_fat -> list[m];
1650 str = g_strdup_printf ("%f\t%f\t%f\n", tmp_proj -> atoms[0][n].x, tmp_proj -> atoms[0][n].y, tmp_proj -> atoms[0][n].z);
1651 print_info (str, NULL, buf);
1652 g_free (str);
1653 }
1654 }
1655 if (tmp_fmol -> next != NULL) tmp_fmol = tmp_fmol -> next;
1656 }
1657}
1658
1659
1660gchar * ens_keyw[4] = {"nve", "nvt", "npt", "nst"};
1661gchar * thermo_keyw[10] = {"evans", "lang", "ander", "ber", "hoover", "gst", "dpd", "mtk", "ttm", "inhomo"};
1662gchar * pseudo_thermo[3] = {"langevin", "gauss", "direct"};
1663gchar * area_keyw[5]={"area", "tens", "tens", "orth", "orth"};
1664gchar * md_keyw[4]={"temp ", "steps ", "integrat ", "pres "};
1665gchar * md_text[4]={"# Target temperature in K\n", "# Number of MD steps\n", "# Integration time step in ps\n", "# Target presssure in k atms\n"};
1666gchar * min_key[3]={"force", "energy", "distance"};
1667//gchar * md_legend[3]={"# Target temperature", "# Number of MD steps", "# MD time step d(t)"};
1668
1674gboolean print_ana ()
1675{
1676 if ((int)tmp_field -> ana_opts[0] || (int)tmp_field -> ana_opts[4] || (int)tmp_field -> ana_opts[8] || (int)tmp_field -> ana_opts[11] || (int)tmp_field -> ana_opts[14])
1677 {
1678 return TRUE;
1679 }
1680 else
1681 {
1682 return FALSE;
1683 }
1684}
1685
1686extern gchar * eval_m[10];
1687extern gchar * eval_vdw[6];
1688extern gchar * io_rw_m[4];
1689extern gchar * io_pres[2];
1690
1691gchar * elec_key[10]={"coul ",
1692 "distan ",
1693 "ewald precision ",
1694 "ewald ",
1695 "reaction ",
1696 "reaction damp ",
1697 "reaction precision ",
1698 "shift ",
1699 "shift damp ",
1700 "shift precision "};
1701gchar * vdw_key[6]={"lor ", "fend", "hoge",
1702 "halg", "tang", "func"};
1703gchar * sys_info[8]={"\n# Ignore particle indices from CONFIG file and set indices by order of reading",
1704 "\n# Ignore strict checks when reading CONFIG file, warning messages and assume safe simulation parameters",
1705 "\n# Skip detailed topology reporting when reading FIELD file",
1706 "\n# Ignore center of mass momentum removal during the simulation",
1707 "\n# Tolerance for the relaxed shell model\n",
1708 "\n# Subcelling threshold density of particle per link cell\n",
1709 "\n# Create an expanded version of the current model\n",
1710 "\n# Restart job from previous run: "};
1711gchar * sys_key[8]={NULL, NULL, NULL, NULL, "rlxtol ", "subcell ", "nfold ", "restart "};
1712gchar * sys_string[8]={"ind", "str", "top", "vom", NULL, NULL, NULL, NULL};
1713gchar * rest_inf[3]={"\n# Continue current simulation - require REVOLD file",
1714 "\n# Start new simulation from older run without temperature reset",
1715 "\n# Start new simulation from older run with temperature reset"};
1716gchar * rest_key[3]={NULL, "noscale", "scale"};
1717
1718gchar * time_inf[2]={"\n\n# Set job time to ", "\n\n# Set job closure time to "};
1719gchar * time_key[2]={"job time ", "close time "};
1720gchar * io_inf[2]={"\n# I/O read interface, with:\n", "\n# I/O write interface, with:\n"};
1721gchar * io_key[2]={"\nio read ", "\nio writ "};
1722gchar * io_meth[4]={"mpiio", "direct", "master", "netcdf"};
1723gchar * io_pec[2]={"off", "on"};
1724gchar * io_typ[2]={"sorted", "unsorted"};
1725
1734void print_int (GtkTextBuffer * buf, int data)
1735{
1736 gchar * str = g_strdup_printf (" %d", data);
1737 print_info (str, "bold_blue", buf);
1738 g_free (str);
1739}
1740
1752void print_control_int (GtkTextBuffer * buf, int data, gchar * info_a, gchar * info_b, gchar * key)
1753{
1754 gchar * str = g_strdup_printf ("%d", data);
1755 print_info (info_a, NULL, buf);
1756 print_info (str, NULL, buf);
1757 g_free (str);
1758 if (info_b != NULL) print_info (info_b, NULL, buf);
1759 print_info ("\n", NULL, buf);
1760 print_info (key, "bold", buf);
1761 print_int (buf, data);
1762}
1763
1772void print_float (GtkTextBuffer * buf, double data)
1773{
1774 gchar * str = g_strdup_printf (" %f", data);
1775 print_info (str, "bold_red", buf);
1776 g_free (str);
1777}
1778
1790void print_control_float (GtkTextBuffer * buf, double data, gchar * info_a, gchar * info_b, gchar * key)
1791{
1792 gchar * str = g_strdup_printf ("%f", data);
1793 print_info (info_a, NULL, buf);
1794 print_info (str, NULL, buf);
1795 g_free (str);
1796 if (info_b != NULL) print_info (info_b, NULL, buf);
1797 print_info ("\n", NULL, buf);
1798 print_info (key, "bold", buf);
1799 print_float (buf, data);
1800}
1801
1810void print_sci (GtkTextBuffer * buf, double data)
1811{
1812 gchar * str = g_strdup_printf (" %e", data);
1813 print_info (str, "bold_orange", buf);
1814 g_free (str);
1815}
1816
1828void print_control_sci (GtkTextBuffer * buf, double data, gchar * info_a, gchar * info_b, gchar * key)
1829{
1830 gchar * str = g_strdup_printf ("%e", data);
1831 print_info (info_a, NULL, buf);
1832 print_info (str, NULL, buf);
1833 if (info_b != NULL) print_info (info_b, NULL, buf);
1834 print_info ("\n", NULL, buf);
1835 print_info (key, "bold", buf);
1836 print_info (str, "bold_orange", buf);
1837 g_free (str);
1838}
1839
1848void print_string (GtkTextBuffer * buf, gchar * string)
1849{
1850 print_info (" ", NULL, buf);
1851 print_info (string, "bold_green", buf);
1852}
1853
1865void print_control_string (GtkTextBuffer * buf, gchar * string, gchar * info_a, gchar * info_b, gchar * key)
1866{
1867 if (info_a != NULL) print_info (info_a, NULL, buf);
1868 if (info_b != NULL) print_info (info_b, NULL, buf);
1869 if (info_a != NULL) print_info ("\n", NULL, buf);
1870 print_info (key, "bold", buf);
1871 if (string) print_string (buf, string);
1872}
1873
1883void print_control_key (GtkTextBuffer * buf, gchar * info, gchar * key)
1884{
1885 if (info != NULL) print_info (info, NULL, buf);
1886 print_info (key, "bold", buf);
1887}
1888
1889
1897void print_dlp_control (GtkTextBuffer * buf)
1898{
1899 int i, j, k;
1900 gchar * str;
1901 gchar * str_a, * str_b, * str_c;
1902
1903 GtkTextIter bStart;
1904 GtkTextIter bEnd;
1905
1906 gtk_text_buffer_get_start_iter (buf, & bStart);
1907 gtk_text_buffer_get_end_iter (buf, & bEnd);
1908 gtk_text_buffer_delete (buf, & bStart, & bEnd);
1909
1910 str = g_strdup_printf ("# DL-POLY CONTROL file created by %s, %s - %d atoms\n\n",
1911 PACKAGE,
1912 prepare_for_title(tmp_proj -> name),
1913 tmp_proj -> natomes);
1914 print_info (str, "bold", buf);
1915 g_free (str);
1916
1917
1918 if (tmp_field -> sys_opts[0] != 1.0)
1919 {
1920 print_control_float (buf, tmp_field -> sys_opts[0], "# Relative dielectric constant = ", NULL, "eps ");
1921 }
1922 if (tmp_field -> sys_opts[1] != 0.0)
1923 {
1924 print_control_float (buf, tmp_field -> sys_opts[1], "\n# Allowing local variation of system density of : ", " \%", "densvar ");
1925 }
1926 for (i=2; i<10; i++)
1927 {
1928 j = (i < 7) ? i : (i == 7) ? 8 : (i == 8) ? 10 : 14;
1929 if (tmp_field -> sys_opts[j] == 1.0)
1930 {
1931 if (i == 9)
1932 {
1933 k = (int)tmp_field -> sys_opts[15];
1934 print_control_string (buf, rest_key[k], sys_info[i-2], rest_inf[k], sys_key[i-2]);
1935 }
1936 else if (i < 6)
1937 {
1938 print_control_string (buf, sys_string[i-2], sys_info[i-2], NULL, "no ");
1939 }
1940 else
1941 {
1942 print_control_key (buf, sys_info[i-2], sys_key[i-2]);
1943 }
1944 if (i == 6 || i == 7)
1945 {
1946 print_float (buf, tmp_field -> sys_opts[j+1]);
1947 }
1948 else if (i == 8)
1949 {
1950 for (k=1; k<4; k++) print_int (buf, (int)tmp_field -> sys_opts[j+k]);
1951 }
1952 print_info ("\n", NULL, buf);
1953 }
1954 }
1955
1956 if (tmp_field -> vdw_opts[0] == 1.0)
1957 {
1958 print_info ("\n# Non bonded short range interactions - type vdW", NULL, buf);
1959 print_control_float (buf, tmp_field -> vdw_opts[1], "\n# van Der Waals short range cutoff = ", " Ang.", "rvdw ");
1960 if (tmp_field -> vdw_opts[2] == 1.0)
1961 {
1962 print_control_string (buf, "direct", "\n# Enforce direct calculation of vdW interactions",
1963 "\n# Do not work with system using tabulated potentials", "vdw ");
1964 }
1965 if (tmp_field -> vdw_opts[3] == 1.0)
1966 {
1967 print_control_string (buf, "shift", "\n# Apply force-shifting for vdW interactions", NULL, "vdw ");
1968 }
1969 if (tmp_field -> vdw_opts[4] == 1.0)
1970 {
1971 print_control_string (buf, vdw_key[(int)tmp_field -> vdw_opts[5]], "\n# Apply mixing rule of type: ", eval_vdw[(int)tmp_field -> vdw_opts[5]], "vdw mix ");
1972 }
1973 }
1974 else
1975 {
1976 print_control_string (buf, "vdw", "\n# No van der Waals interactions (short range)", NULL, "no ");
1977 }
1978 print_info ("\n\n", NULL, buf);
1979
1980 if (tmp_field -> elec_opts[0] == 1.0)
1981 {
1982 print_info ("\n# Non bonded long range interactions", NULL, buf);
1983 print_control_float (buf, tmp_field -> elec_opts[1], "\n# Electrostatics long range cutoff = ", " Ang.", "cut ");
1984 if (tmp_field -> elec_opts[2] == 1.0)
1985 {
1986 print_control_float (buf, tmp_field -> elec_opts[3], "\n# Use optional padding to the cutoff = ", " Ang.", "pad ");
1987 }
1988 if (tmp_field -> elec_opts[4] == 1.0)
1989 {
1990 print_control_key (buf, "\n# Use extended coulombic exclusion\n", "exclu");
1991 }
1992 print_info ("\n# Electrostatics calculated using ", NULL, buf);
1993 print_info (eval_m[(int)tmp_field -> elec_opts[5]], NULL, buf);
1994 print_info ("\n", NULL, buf);
1995 print_info (elec_key[(int)tmp_field -> elec_opts[5]], "bold", buf);
1996 if (tmp_field -> elec_opts[5] == 2.0 || tmp_field -> elec_opts[5] == 6.0 || tmp_field -> elec_opts[5] == 9.0)
1997 {
1998 print_sci (buf, tmp_field -> elec_opts[6]);
1999 }
2000 else if (tmp_field -> elec_opts[5] == 3.0)
2001 {
2002 print_float (buf, tmp_field -> elec_opts[6]);
2003 for (k=7; k<10; k++) print_int (buf, (int)tmp_field -> elec_opts[k]);
2004 }
2005 else if (tmp_field -> elec_opts[5] == 5.0 || tmp_field -> elec_opts[5] == 8.0)
2006 {
2007 print_float (buf, tmp_field -> elec_opts[6]);
2008 }
2009
2010 if (tmp_field -> elec_opts[5] == 2.0 || tmp_field -> elec_opts[5] == 3.0)
2011 {
2012 print_control_int (buf, (int)tmp_field -> elec_opts[10], "\n# Evaluate k space contribution to the Ewald sum every: ", " MD step(s)", "ewald evalu ");
2013 }
2014 }
2015 else
2016 {
2017 print_control_string (buf, "elec", "# No electrostatics interactions (long range)", NULL, "no ");
2018 }
2019 print_info ("\n", NULL, buf);
2020
2021 if (tmp_field -> met_opts[0] == 1.0 || tmp_field -> met_opts[1] == 1.0)
2022 {
2023 print_info ("\n# Metallic interactions", NULL, buf);
2024 }
2025 if (tmp_field -> met_opts[0] == 1.0)
2026 {
2027 print_control_string (buf, "direct", "\n# Enforce direct calculation of metal interactions", "\n# This does not work with metal alloy systems using the *EAM* potentials", "metal ");
2028 }
2029 if (tmp_field -> met_opts[1] == 1.0)
2030 {
2031 print_control_string (buf, "sqrtrho", "\n# Switch the TABEAM default embedding functions, F, from F(ρ) to F(√ρ)", NULL, "metal ");
2032 }
2033 if (tmp_field -> met_opts[0] == 1.0 || tmp_field -> met_opts[1] == 1.0) print_info ("\n", NULL, buf);
2034
2035 print_control_string (buf, ens_keyw[tmp_field -> ensemble], "\n# Thermostat information", NULL, "ensemble ");
2036 if (tmp_field -> ensemble)
2037 {
2038 switch (tmp_field -> ensemble)
2039 {
2040 case 1:
2041 i = (tmp_field -> thermostat > 6) ? tmp_field -> thermostat + 1 : tmp_field -> thermostat;
2042 break;
2043 default:
2044 i = (! tmp_field -> thermostat) ? 1 : (tmp_field -> thermostat == 3) ? 7 : tmp_field -> thermostat + 2;
2045 break;
2046 }
2047 print_string (buf, thermo_keyw[i]);
2048 if (tmp_field -> ensemble > 1 || tmp_field -> thermostat)
2049 {
2050 if (tmp_field -> ensemble == 1 && tmp_field -> thermostat == 6)
2051 {
2052 str = g_strdup_printf ("s%1d", (int)tmp_field -> thermo_opts[0]+1);
2053 print_string (buf, str);
2054 g_free (str);
2055 print_float (buf, tmp_field -> thermo_opts[1]);
2056 }
2057 else
2058 {
2059 if (tmp_field -> ensemble == 3 && tmp_field -> thermo_opts[3] > 0.0) print_string (buf, "Q");
2060 print_float (buf, tmp_field -> thermo_opts[0]);
2061 if (tmp_field -> ensemble != 1 || (tmp_field -> thermostat == 2 || tmp_field -> thermostat == 5 || tmp_field -> thermostat > 6))
2062 {
2063 print_float (buf, tmp_field -> thermo_opts[1]);
2064 if (tmp_field -> ensemble == 1 && tmp_field -> thermostat > 6) print_float (buf, tmp_field -> thermo_opts[2]);
2065 }
2066 }
2067 if (tmp_field -> ensemble == 3 && tmp_field -> thermo_opts[3] > 0.0)
2068 {
2069 i = (int)tmp_field -> thermo_opts[3] - 1;
2070 print_string (buf, area_keyw[i]);
2071 if (tmp_field -> thermo_opts[3] == 2.0) print_float (buf, tmp_field -> thermo_opts[4]);
2072 }
2073 if (tmp_field -> thermo_opts[5] == 1.0) print_string (buf, "semi");
2074 }
2075 }
2076
2077 print_info ("\n\n", NULL, buf);
2078 if (tmp_field -> thermo_opts[6] == 1.0)
2079 {
2080 print_info ("# Attach a pseudo thermal bath with:\n", NULL, buf);
2081 if (tmp_field -> thermo_opts[7] > 0.0)
2082 {
2083 str = g_strdup_printf ("# - thermostat of type: %s\n", pseudo_thermo[(int)tmp_field -> thermo_opts[7] - 1]);
2084 }
2085 else
2086 {
2087 str = g_strdup_printf ("# - thermostats of type Langevin and Direct applied successively\n");
2088 }
2089 print_info (str, NULL, buf);
2090 g_free (str);
2091 str = g_strdup_printf ("# - thickness of thermostat layer to MD cell boundaries: %f Ang.\n", tmp_field -> thermo_opts[8]);
2092 print_info (str, NULL, buf);
2093 g_free (str);
2094 if (tmp_field -> thermo_opts[9] > 0.0)
2095 {
2096 str = g_strdup_printf ("# - Target temperature: %f K\n", tmp_field -> thermo_opts[9]);
2097 print_info (str, NULL, buf);
2098 g_free (str);
2099 }
2100 else
2101 {
2102 print_info ("# - Target temperature: system target temperature\n", NULL, buf);
2103 }
2104 print_info ("pseudo ", "bold", buf);
2105 if (tmp_field -> thermo_opts[7] > 0.0)
2106 {
2107 print_info (pseudo_thermo[(int)tmp_field -> thermo_opts[7] - 1], "bold_green", buf);
2108 }
2109 print_float (buf, tmp_field -> thermo_opts[8]);
2110 if (tmp_field -> thermo_opts[9] > 0.0) print_float (buf, tmp_field -> thermo_opts[9]);
2111 print_info ("\n\n", NULL, buf);
2112 }
2113
2114 // MD information
2115 print_info ("# Molecular dynamics information\n", NULL, buf);
2116 for (i=0; i<2+(int)tmp_field -> md_opts[1]; i++)
2117 {
2118 print_control_key (buf, md_text[i], md_keyw[i]);
2119 switch (i)
2120 {
2121 case 0:
2122 print_float (buf, tmp_field -> md_opts[0]);
2123 if (tmp_field -> ensemble > 1)
2124 {
2125 print_info ("\n", NULL, buf);
2126 print_info (md_keyw[3], "bold", buf);
2127 print_float (buf, tmp_field -> md_opts[5]);
2128 }
2129 break;
2130 case 1:
2131 print_int (buf, (int)tmp_field -> md_opts[2]);
2132 break;
2133 case 2:
2134 print_info ("leapfrog", "bold_green", buf);
2135 break;
2136 }
2137 print_info ("\n", NULL, buf);
2138 }
2139
2140 if (tmp_field -> md_opts[3] == 1.0)
2141 {
2142 print_control_float (buf, tmp_field -> md_opts[4], "# Variable time step, initial time step = ", " ps", "variable timestep ");
2143 print_control_float (buf, tmp_field -> md_opts[6], "\n# Maximum time step allowed = ", " ps", "mxstep ");
2144 print_control_float (buf, tmp_field -> md_opts[7], "\n# Maximum move allowed = ", " Ang.", "maxdis ");
2145 print_control_float (buf, tmp_field -> md_opts[8], "\n# Minimum move allowed = ", " Ang.", "mindis ");
2146 }
2147 else
2148 {
2149 print_control_float (buf, tmp_field -> md_opts[4], "# MD time step = ", " fs", "timestep ");
2150 }
2151
2152 print_control_int (buf, (int)tmp_field -> md_opts[9], "\n# Shake / Rattle iterations limit: ", " cycle(s)", "mxshak ");
2153 print_control_sci (buf, tmp_field -> md_opts[10], "\n# Shake / Rattle tolerance: ", NULL, "shake ");
2154
2155 if (tmp_field -> md_opts[1] == 1.0)
2156 {
2157 print_control_int (buf, (int)tmp_field -> md_opts[11], "\n# FIQA iterations limit: ", " cycle(s)", "mxquat ");
2158 print_control_sci (buf, tmp_field -> md_opts[12], "\n# Quaternion tolerance: ", NULL, "quater ");
2159 }
2160
2161 if (tmp_field -> md_opts[13] == 1.0)
2162 {
2163 print_info ("\n\n# Initiate impact on particle\n# - with particle index: ", NULL, buf);
2164 str = g_strdup_printf ("%d", (int)tmp_field -> md_opts[14]);
2165 print_info (str, NULL, buf);
2166 print_info ("\n# - at MD step: ", NULL, buf);
2167 str = g_strdup_printf ("%d", (int)tmp_field -> md_opts[15]);
2168 print_info (str, NULL, buf);
2169 g_free (str);
2170 print_info ("\n# - with energy (k eV): ", NULL, buf);
2171 str = g_strdup_printf ("%f", tmp_field -> md_opts[16]);
2172 print_info (str, NULL, buf);
2173 g_free (str);
2174 print_info ("\n# - direction (x, y, z): ", NULL, buf);
2175 str = g_strdup_printf ("%f %f %f", tmp_field -> md_opts[17], tmp_field -> md_opts[18], tmp_field -> md_opts[19]);
2176 print_info (str, NULL, buf);
2177 g_free (str);
2178 print_info ("\n", NULL, buf);
2179 print_info ("impact ", "bold", buf);
2180 for (k=14; k<16; k++) print_int (buf, (int)tmp_field -> md_opts[k]);
2181 for (k=16; k<20; k++) print_float (buf, tmp_field -> md_opts[k]);
2182 }
2183
2184
2185 if (tmp_field -> equi_opts[0] == 1.0)
2186 {
2187 // Equilibration information
2188 print_info ("\n\n# Equilibration information", NULL, buf);
2189 print_control_int (buf, (int)tmp_field -> equi_opts[1], "\n# Equilibrate during: ", " MD step(s)", "equil ");
2190 if (tmp_field -> equi_opts[2] == 1.0)
2191 {
2192 print_control_int (buf, (int)tmp_field -> equi_opts[3], "\n# During equilibration: rescale system temperature every: ", " MD step(s)", "scale ");
2193 }
2194 if (tmp_field -> equi_opts[4] == 1.0)
2195 {
2196 print_control_float (buf, tmp_field -> equi_opts[5], "\n# During equilibration: cap forces, with fmax= ", " Kb T Ang-1", "cap ");
2197 }
2198 if (tmp_field -> equi_opts[6] == 1.0)
2199 {
2200 print_control_int (buf, (int)tmp_field -> equi_opts[7], "\n# During equilibration: resample the instantaneous momenta distribution every: ", " MD step(s)", "regaus ");
2201 }
2202 if (tmp_field -> equi_opts[8] == 1.0)
2203 {
2204 str = g_strdup_printf ("\n# Every %d step(s) during equilibration: minimize %s with target %s= %f\n",
2205 (int)tmp_field -> equi_opts[11], min_key[(int)tmp_field -> equi_opts[9]], min_key[(int)tmp_field -> equi_opts[9]], tmp_field -> equi_opts[10]);
2206 print_control_key (buf, str, "minim ");
2207 g_free (str);
2208 print_string (buf, min_key[(int)tmp_field -> equi_opts[9]]);
2209 print_int (buf, (int)tmp_field -> equi_opts[11]);
2210 print_float (buf, tmp_field -> equi_opts[10]);
2211 print_info ("\n", NULL, buf);
2212 }
2213 if (tmp_field -> equi_opts[12] == 1.0)
2214 {
2215 str = g_strdup_printf ("# At the start of the equilibration: minimize %s with target %s= %f\n",
2216 min_key[(int)tmp_field -> equi_opts[13]], min_key[(int)tmp_field -> equi_opts[13]], tmp_field -> equi_opts[14]);
2217 print_control_key (buf, str, "optim ");
2218 g_free (str);
2219 print_string (buf, min_key[(int)tmp_field -> equi_opts[13]]);
2220 print_float (buf, tmp_field -> equi_opts[14]);
2221 print_info ("\n", NULL, buf);
2222 }
2223 if (tmp_field -> equi_opts[15] == 1.0)
2224 {
2225 print_control_key (buf, "# During equilibration: perform a zero temperature MD minimization\n", "zero");
2226 print_info ("\n", NULL, buf);
2227 }
2228 if (tmp_field -> equi_opts[16] == 1.0)
2229 {
2230 print_control_key (buf, "# Include equilibration data in overall statistics\n", "collect");
2231 print_info ("\n", NULL, buf);
2232 }
2233 }
2234
2235 if (print_ana())
2236 {
2237 print_info ("\n# Analysis information", NULL, buf);
2238 if (tmp_field -> ana_opts[0] == 1.0)
2239 {
2240 print_control_string (buf, "all", "\n# Calculate and collect all intra-molecular PDFs", NULL, "ana ");
2241 for (k=1; k<3; k++) print_int (buf, (int)tmp_field -> ana_opts[k]);
2242 print_float (buf, tmp_field -> ana_opts[3]);
2243 }
2244 if (tmp_field -> ana_opts[4] == 1.0)
2245 {
2246 print_control_string (buf, "bon", "\n# Calculate and collect bonds PDFs", NULL, "ana ");
2247 for (k=5; k<6; k++) print_int (buf, (int)tmp_field -> ana_opts[k]);
2248 print_float (buf, tmp_field -> ana_opts[7]);
2249 }
2250 if (tmp_field -> ana_opts[8] == 1.0)
2251 {
2252 print_control_string (buf, "ang", "\n# Calculate and collect angles PDFs", NULL, "ana ");
2253 for (k=9; k<11; k++) print_int (buf, (int)tmp_field -> ana_opts[k]);
2254 }
2255 if (tmp_field -> ana_opts[11] == 1.0)
2256 {
2257 print_control_string (buf, "dih", "\n# Calculate and collect dihedrals PDFs", NULL, "ana ");
2258 for (k=12; k<14; k++) print_int (buf, (int)tmp_field -> ana_opts[k]);
2259 }
2260 if (tmp_field -> ana_opts[14] == 1.0)
2261 {
2262 print_control_string (buf, "inv", "\n# Calculate and collect inversions PDFs", NULL, "ana ");
2263 for (k=15; k<17; k++) print_int (buf, (int)tmp_field -> ana_opts[k]);
2264 }
2265 print_control_string (buf, "ana", "\n# Print any opted for analysis inter and intra-molecular PDFs", NULL, "print ");
2266 }
2267
2268 print_info ("\n", NULL, buf);
2269
2270 if (tmp_field -> out_opts[21] == 1.0 || tmp_field -> out_opts[27] == 1.0)
2271 {
2272 print_control_float (buf, tmp_field -> out_opts[23], "\n# Bin size for RDfs and Z-density distribution: ", " Ang.", "binsize ");
2273 }
2274 if (tmp_field -> out_opts[21] == 1.0)
2275 {
2276 print_control_int (buf, (int)tmp_field -> out_opts[22], "\n# Calculate and collect radial distribution functions every: ", " MD step(s)", "rdf ");
2277 print_info ("\n", NULL, buf);
2278 print_control_string (buf, "rdf", NULL, NULL, "print ");
2279 }
2280 if (tmp_field -> out_opts[27] == 1.0)
2281 {
2282 print_control_int (buf, (int)tmp_field -> out_opts[28], "\n# Calculate and collect Z-density profile every: ", " MD step(s)", "zden ");
2283 print_info ("\n", NULL, buf);
2284 print_control_string (buf, "zden", NULL, NULL, "print ");
2285 }
2286 if (tmp_field -> out_opts[24] == 1.0)
2287 {
2288 print_control_key (buf, "\n# Velocity autocorrelation functions, VAFs\n", "vaf ");
2289 for (k=25; k<27; k++) print_int (buf, (int)tmp_field -> out_opts[k]);
2290 print_info ("\n", NULL, buf);
2291 print_control_string (buf, "vaf", NULL, NULL, "print ");
2292 if (tmp_field -> out_opts[29] == 1.0)
2293 {
2294 print_control_string (buf, "vafav", "\n# Ignore time averaging for the VAFs", NULL, "no ");
2295 }
2296 }
2297
2298 if ((int)tmp_field -> out_opts[0] || (int)tmp_field -> out_opts[4] || (int)tmp_field -> out_opts[8]
2299 || (int)tmp_field -> out_opts[12] || (int)tmp_field -> out_opts[15] || (int)tmp_field -> out_opts[17] || (int)tmp_field -> out_opts[19])
2300 {
2301 print_info ("\n\n# Output information", NULL, buf);
2302 if ((int)tmp_field -> out_opts[0])
2303 {
2304 print_control_key (buf, "\n# Write defects trajectory file, DEFECTS\n", "defe ");
2305 for (k=1; k<3; k++) print_int (buf, (int)tmp_field -> out_opts[k]);
2306 print_float (buf, tmp_field -> out_opts[3]);
2307 }
2308 if ((int)tmp_field -> out_opts[4])
2309 {
2310 print_control_key (buf, "\n# Write displacement trajectory file, RSDDAT\n", "disp ");
2311 for (k=5; k<7; k++) print_int (buf, (int)tmp_field -> out_opts[k]);
2312 print_float (buf, tmp_field -> out_opts[7]);
2313 }
2314 if ((int)tmp_field -> out_opts[8])
2315 {
2316 print_control_key (buf, "\n# Write HISTORY file\n", "traj ");
2317 for (k=9; k<11; k++) print_int (buf, (int)tmp_field -> out_opts[k]);
2318 print_float (buf, tmp_field -> out_opts[11]);
2319 }
2320 if ((int)tmp_field -> out_opts[12])
2321 {
2322 print_control_key (buf, "\n# Write MSDTMP file\n", "msdtmp ");
2323 for (k=13; k<15; k++) print_int (buf, (int)tmp_field -> out_opts[k]);
2324 }
2325 if ((int)tmp_field -> out_opts[15])
2326 {
2327 print_control_int (buf, (int)tmp_field -> out_opts[16], "\n# Print system data every: ", " MD step(s)", "print ");
2328 }
2329 if ((int)tmp_field -> out_opts[17])
2330 {
2331 print_control_int (buf, (int)tmp_field -> out_opts[18], "\n# Accumulate statics data every: ", " MD step(s)", "stats ");
2332 }
2333 if ((int)tmp_field -> out_opts[19])
2334 {
2335 print_control_int (buf, (int)tmp_field -> out_opts[20], "\n# Set rolling average stack to: ", " MD step(s)", "stack ");
2336 }
2337 }
2338
2339 print_control_int (buf, (int)tmp_field -> out_opts[30], "\n# Dump restart information every: ", " MD step(s)", "dump ");
2340
2341 for (i=0; i<2; i++)
2342 {
2343 if (tmp_field -> io_opts[2*i] == 1.0)
2344 {
2345 print_control_float (buf, tmp_field -> io_opts[2*i+1], time_inf[i], " s", time_key[i]);
2346 }
2347 }
2348 print_info ("\n", NULL, buf);
2349 for (i=0; i<2; i++)
2350 {
2351 j=4 + i*6;
2352 if (tmp_field -> io_opts[j] == 1.0)
2353 {
2354 j ++;
2355 print_info (io_inf[i], NULL, buf);
2356 print_info ("# - method = ", NULL, buf);
2357 print_info (io_rw_m[(int)tmp_field -> io_opts[j]], NULL, buf);
2358 j++;
2359 if (i)
2360 {
2361 if (tmp_field -> io_opts[j-1] == 3.0)
2362 {
2363 print_info ("\n# - precision = ", NULL, buf);
2364 print_info (io_pres[(int)tmp_field -> io_opts[j]], NULL, buf);
2365 }
2366 j ++;
2367 print_info ("\n# - type = ", NULL, buf);
2368 print_info (io_typ[(int)tmp_field -> io_opts[j]], NULL, buf);
2369 j++;
2370 }
2371 if (tmp_field -> io_opts[4+7*i] != 2.0)
2372 {
2373 print_info ("\n# - j, reader count = ", NULL, buf);
2374 str_a = g_strdup_printf ("%d", (int)tmp_field -> io_opts[j]);
2375 print_info (str_a, NULL, buf);
2376 }
2377 j++;
2378 if (tmp_field -> io_opts[4+7*i] != 2.0)
2379 {
2380 print_info ("\n# - k, batch size = ", NULL, buf);
2381 str_b = g_strdup_printf ("%d", (int)tmp_field -> io_opts[j]);
2382 print_info (str_b, NULL, buf);
2383 }
2384 j++;
2385 print_info ("\n# - l, buffer size = ", NULL, buf);
2386 str_c = g_strdup_printf ("%d", (int)tmp_field -> io_opts[j]);
2387 print_info (str_c, NULL, buf);
2388 j++;
2389 if (tmp_field -> io_opts[4+7*i] != 2.0)
2390 {
2391 print_info ("\n# - e, parallel error check is ", NULL, buf);
2392 print_info (io_pec[(int)tmp_field -> io_opts[j]], NULL, buf);
2393 }
2394 print_info (io_key[i], "bold", buf);
2395 print_info (io_meth[(int)tmp_field -> io_opts[5+6*i]], "bold_green", buf);
2396 if (i)
2397 {
2398 if (tmp_field -> io_opts[11] == 3.0)
2399 {
2400 print_info (io_pres[(int)tmp_field -> io_opts[12]], "bold_green", buf);
2401 }
2402 print_info (" ", NULL, buf);
2403 print_info (io_typ[(int)tmp_field -> io_opts[13]], "bold_green", buf);
2404 }
2405 if (tmp_field -> io_opts[4+7*i] != 2.0)
2406 {
2407 print_info (" ", NULL, buf);
2408 print_info (str_a, "bold_blue", buf);
2409 g_free (str_a);
2410 print_info (" ", NULL, buf);
2411 print_info (str_b, "bold_blue", buf);
2412 g_free (str_b);
2413 }
2414 print_info (" ", NULL, buf);
2415 print_info (str_c, "bold_blue", buf);
2416 g_free (str_c);
2417 if (tmp_field -> io_opts[4+7*i] != 2.0)
2418 {
2419 (tmp_field -> io_opts[j] == 0.0) ? print_string (buf, "N") : print_string (buf, "Y");
2420 }
2421 print_info ("\n", NULL, buf);
2422 j++;
2423 }
2424 }
2425 if (tmp_field -> io_opts[18] == 1.0)
2426 {
2427 print_control_key (buf, "\n# Seeds for the random number generators\n", "seed ");
2428 for (i=19; i<22; i++) print_int (buf, (int)tmp_field -> io_opts[i]);
2429 print_info ("\n", NULL, buf);
2430 }
2431 if (tmp_field -> io_opts[22] == 1.0)
2432 {
2433 print_control_key (buf, "\n# Limits to 2 the number of processors in z direction for slab simulations\n", "slab");
2434 }
2435
2436 print_info ("\n\n", NULL, buf);
2437 print_info ("finish", "bold", buf); // Close the CONTROL file
2438}
insertion_menu mol[]
Definition w_library.c:193
gchar * param[2]
int atoms[NUM_STYLES][2]
field_prop * get_active_prop_using_atoms(struct field_prop *pr, int ti, int *ids)
retrieve field molecule structural property using atoms
Definition dlp_active.c:286
field_atom * get_active_atom(int a, int b)
retrieve field atom
Definition dlp_active.c:145
float val
Definition dlp_init.c:117
gchar * sys_opts[10]
int * atoms_id
double *** cross
Definition dlp_edit.c:418
char * fvars_bond[2][FBONDS][FBONDS_P]
Definition dlp_field.c:306
char * fvars_dihedral[2][FDIHEDRAL][FDIHEDRAL_P]
Definition dlp_field.c:450
char * fvars_angle[2][FANGLES][FANGLES_P]
Definition dlp_field.c:372
field_constraint * tmp_fcons
Definition dlp_field.c:959
int fvalues[2][15][21]
Definition dlp_field.c:255
field_nth_body * tmp_fbody
Definition dlp_field.c:965
project * tmp_proj
Definition dlp_field.c:953
field_shell * tmp_fshell
Definition dlp_field.c:958
field_rigid * tmp_frig
Definition dlp_field.c:961
int body_at(int b)
find the number of atom(s) in a non bonded interaction
Definition dlp_field.c:1022
field_pmf * tmp_fpmf
Definition dlp_field.c:960
gchar * fnames[2][16][21]
Definition dlp_field.c:218
field_external * tmp_fext
Definition dlp_field.c:967
field_atom * tmp_fat
Definition dlp_field.c:957
int struct_id(int f)
number of atoms in a structural element
Definition dlp_field.c:999
field_struct * tmp_fstr
Definition dlp_field.c:964
gchar * parameters_info(int obj, int key, gchar **words, float *data)
prepare classical force field parameter description string
Definition dlp_field.c:1097
field_molecule * tmp_fmol
Definition dlp_field.c:955
field_atom * tmp_fdt
Definition dlp_field.c:957
field_prop * tmp_fprop
Definition dlp_field.c:963
field_tethered * tmp_ftet
Definition dlp_field.c:962
field_atom * tmp_fct
Definition dlp_field.c:957
gchar * fkeysw[2][16][21]
Definition dlp_field.c:185
field_atom * tmp_fbt
Definition dlp_field.c:957
classical_field * tmp_field
Definition dlp_field.c:951
Variable declarations for the creation of the DL_POLY input file(s)
angle dihedral_3d(cell_info *cell, int mdstep, atom *at, atom *bt, atom *ct, atom *dt)
dihedral between atom a, b, c and d in 3D
Definition ogl_utils.c:204
angle inversion_3d(cell_info *cell, int mdstep, atom *at, atom *bt, atom *ct, atom *dt)
inversion angle between atom a, b, c and d in 3D
Definition ogl_utils.c:240
#define SEXTERN
Definition dlp_field.h:118
angle angle_3d(cell_info *cell, int mdstep, atom *at, atom *bt, atom *ct)
angle between atom a, b and c in 3D
Definition ogl_utils.c:179
gboolean afp[MAXDATA]
distance distance_3d(cell_info *cell, int mdstep, atom *at, atom *bt)
distance between atom a and b in 3D
Definition ogl_utils.c:81
int multi
Definition dlp_init.c:121
gchar * eval_vdw[6]
gchar * vdw_key[6]
Definition dlp_print.c:1701
gboolean print_this_imp_inv(imp_inv *inv, int di, int a, int b, int c, int d)
print this improper / inversion structure or not (already printed) ?
Definition dlp_print.c:239
gchar * md_keyw[4]
Definition dlp_print.c:1664
gchar * pseudo_thermo[3]
Definition dlp_print.c:1662
gchar * io_typ[2]
Definition dlp_print.c:1724
void print_dlp_tet(GtkTextBuffer *buf, field_tethered *tet)
print force field tethered potential
Definition dlp_print.c:915
void print_dlp_bond(int bi, GtkTextBuffer *buf, field_struct *bd, int fi, GtkTreeStore *store, GtkTreeIter *iter)
print / fill tree store with force field bond(s) information
Definition dlp_print.c:774
void print_sci(GtkTextBuffer *buf, double data)
print float in scientific format
Definition dlp_print.c:1810
void print_dlp_dihedral(int dih, GtkTextBuffer *buf, field_struct *dh, int fi, GtkTreeStore *store, GtkTreeIter *iter)
print / fill tree store with force field dihedral(s) information
Definition dlp_print.c:499
void print_all_field_struct(field_molecule *mol, int str)
print all field structural element(s)
Definition dlp_print.c:195
gchar * io_inf[2]
Definition dlp_print.c:1720
gchar * ens_keyw[4]
Definition dlp_print.c:1660
void print_dlp_atom(GtkTextBuffer *buf, int at, int numat)
print force field atom
Definition dlp_print.c:1003
void print_dlp_tersoff(GtkTextBuffer *buf, field_nth_body *body)
print force field Tersoff potential
Definition dlp_print.c:1351
void print_dlp_pmf(GtkTextBuffer *buf, field_pmf *pmf)
print force field mean force potential
Definition dlp_print.c:939
void print_control_int(GtkTextBuffer *buf, int data, gchar *info_a, gchar *info_b, gchar *key)
print CONTROL file print integer value
Definition dlp_print.c:1752
void print_control_string(GtkTextBuffer *buf, gchar *string, gchar *info_a, gchar *info_b, gchar *key)
print CONTROL file print string
Definition dlp_print.c:1865
void print_dlp_improper_inversion(int di, GtkTextBuffer *buf, field_struct *dhii, int fi, GtkTreeStore *store, GtkTreeIter *iter)
print / fill tree store with force field improper(s)/inversion(s) information
Definition dlp_print.c:299
void print_dlp_angle(int ai, GtkTextBuffer *buf, field_struct *an, int fi, GtkTreeStore *store, GtkTreeIter *iter)
print / fill tree store with force field angle(s) information
Definition dlp_print.c:649
void print_dlp_tersoff_cross(GtkTextBuffer *buf, field_nth_body *body_a, field_nth_body *body_b)
print Tersoff potential cross term
Definition dlp_print.c:1326
void print_dlp_molecule(GtkTextBuffer *buf, field_molecule *fmol)
print force field molecule
Definition dlp_print.c:1058
gchar * io_pres[2]
gchar * sys_key[8]
Definition dlp_print.c:1711
void print_dlp_control(GtkTextBuffer *buf)
print DL-POLY CONTROL file
Definition dlp_print.c:1897
gchar * io_rw_m[4]
gchar * time_key[2]
Definition dlp_print.c:1719
gchar * sys_string[8]
Definition dlp_print.c:1712
gboolean member_of_atom(field_atom *fat, int id)
is the id atom from the model in target field atom
Definition dlp_print.c:277
void print_control_float(GtkTextBuffer *buf, double data, gchar *info_a, gchar *info_b, gchar *key)
print CONTROL file print float value
Definition dlp_print.c:1790
int get_num_vdw_max()
Get the number of field shell interactions.
Definition dlp_edit.c:2002
int get_pbc()
get the PBC DL-POLY lattice type
Definition dlp_print.c:1536
void print_float(GtkTextBuffer *buf, double data)
print float value
Definition dlp_print.c:1772
gchar * sys_info[8]
Definition dlp_print.c:1703
gchar * io_pec[2]
Definition dlp_print.c:1723
gchar * eval_m[10]
gchar * thermo_keyw[10]
Definition dlp_print.c:1661
gchar * time_inf[2]
Definition dlp_print.c:1718
gchar * io_meth[4]
Definition dlp_print.c:1722
gchar * rest_key[3]
Definition dlp_print.c:1716
void print_field_prop(field_prop *pro, int st, field_molecule *mol)
print force field property
Definition dlp_print.c:94
void print_dlp_field(GtkTextBuffer *buf)
print DL-POLY classical force field
Definition dlp_print.c:1401
void print_dlp_config(GtkTextBuffer *buf)
print DL-POLY CONFIG file
Definition dlp_print.c:1575
gchar * area_keyw[5]
Definition dlp_print.c:1663
int get_num_struct_to_print(field_molecule *fmol, int sid)
find the number of structural element(s) to print
Definition dlp_print.c:1026
void print_dlp_body(GtkTextBuffer *buf, field_nth_body *body)
print force field non bonded potential
Definition dlp_print.c:1292
gboolean in_bond(int at, int bd[2])
is atom at in bond bd
Definition dlp_init.c:892
void print_string(GtkTextBuffer *buf, gchar *string)
print string
Definition dlp_print.c:1848
gboolean print_ana()
determine if the analysis information section is required
Definition dlp_print.c:1674
void print_field_struct(field_struct *stru, field_molecule *mol)
print force field structural element
Definition dlp_print.c:154
gchar * rest_inf[3]
Definition dlp_print.c:1713
void print_control_sci(GtkTextBuffer *buf, double data, gchar *info_a, gchar *info_b, gchar *key)
print CONTROL file print float value in scientific format
Definition dlp_print.c:1828
void print_int(GtkTextBuffer *buf, int data)
print integer value
Definition dlp_print.c:1734
void print_dlp_shell(GtkTextBuffer *buf, field_molecule *fmol, field_shell *shell)
print force field core shell
Definition dlp_print.c:986
void print_dlp_rigid(GtkTextBuffer *buf, field_rigid *rig)
print force field rigid
Definition dlp_print.c:877
gchar * md_text[4]
Definition dlp_print.c:1665
void print_dlp_cons(GtkTextBuffer *buf, field_constraint *cons)
print force field constraint
Definition dlp_print.c:969
gchar * elec_key[10]
Definition dlp_print.c:1691
gchar * get_body_element_name(field_nth_body *body, int aid, int nbd)
get field body potential element name
Definition dlp_edit.c:1934
gchar * io_key[2]
Definition dlp_print.c:1721
void print_control_key(GtkTextBuffer *buf, gchar *info, gchar *key)
print CONTROL file print key
Definition dlp_print.c:1883
gchar * min_key[3]
Definition dlp_print.c:1666
gboolean * allocbool(int val)
allocate a gboolean * pointer
Definition global.c:254
int activef
Definition global.c:161
int * allocint(int val)
allocate an int * pointer
Definition global.c:314
struct thermostat thermostat
Definition global.h:733
gchar * prepare_for_title(gchar *init)
prepare a string for a window title, getting rid of all markup
Definition tools.c:71
void print_info(gchar *str, gchar *stag, GtkTextBuffer *buffer)
print information in GtkTextBuffer
Definition interface.c:738
Messaging function declarations.
integer(kind=c_int) function molecules(frag_and_mol, allbonds)
double z
Definition ogl_draw.c:57
double y
Definition ogl_draw.c:57
double x
Definition ogl_draw.c:57
double angle
Definition glwin.h:115
double length
Definition glwin.h:123
imp_inv * prev
Definition dlp_print.c:224
imp_inv * next
Definition dlp_print.c:223
int b
Definition tab-1.c:95
int c
Definition tab-1.c:95
int d
Definition tab-1.c:95
int a
Definition tab-1.c:95