2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * Copyright (c) 2013, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gmx_fatal.h"
51 /* The source code in this file should be thread-safe.
52 Please keep it that way. */
54 /* This number should be increased whenever the file format changes! */
55 static const int enx_version = 5;
57 const char *enx_block_id_name[] = {
58 "Averaged orientation restraints",
59 "Instantaneous orientation restraints",
60 "Orientation restraint order tensor(s)",
61 "Distance restraints",
68 /* Stuff for reading pre 4.1 energy files */
70 gmx_bool bOldFileOpen; /* Is this an open old file? */
71 gmx_bool bReadFirstStep; /* Did we read the first step? */
72 int first_step; /* First step in the energy file */
73 int step_prev; /* Previous step */
74 int nsum_prev; /* Previous step sum length */
75 t_energy *ener_prev; /* Previous energy sums */
86 static void enxsubblock_init(t_enxsubblock *sb)
90 sb->type = xdr_datatype_double;
92 sb->type = xdr_datatype_float;
108 static void enxsubblock_free(t_enxsubblock *sb)
144 for (i = 0; i < sb->sval_alloc; i++)
157 /* allocate the appropriate amount of memory for the given type and nr */
158 static void enxsubblock_alloc(t_enxsubblock *sb)
160 /* allocate the appropriate amount of memory */
163 case xdr_datatype_float:
164 if (sb->nr > sb->fval_alloc)
166 srenew(sb->fval, sb->nr);
167 sb->fval_alloc = sb->nr;
170 case xdr_datatype_double:
171 if (sb->nr > sb->dval_alloc)
173 srenew(sb->dval, sb->nr);
174 sb->dval_alloc = sb->nr;
177 case xdr_datatype_int:
178 if (sb->nr > sb->ival_alloc)
180 srenew(sb->ival, sb->nr);
181 sb->ival_alloc = sb->nr;
184 case xdr_datatype_int64:
185 if (sb->nr > sb->lval_alloc)
187 srenew(sb->lval, sb->nr);
188 sb->lval_alloc = sb->nr;
191 case xdr_datatype_char:
192 if (sb->nr > sb->cval_alloc)
194 srenew(sb->cval, sb->nr);
195 sb->cval_alloc = sb->nr;
198 case xdr_datatype_string:
199 if (sb->nr > sb->sval_alloc)
203 srenew(sb->sval, sb->nr);
204 for (i = sb->sval_alloc; i < sb->nr; i++)
208 sb->sval_alloc = sb->nr;
212 gmx_incons("Unknown block type: this file is corrupted or from the future");
216 static void enxblock_init(t_enxblock *eb)
224 static void enxblock_free(t_enxblock *eb)
226 if (eb->nsub_alloc > 0)
229 for (i = 0; i < eb->nsub_alloc; i++)
231 enxsubblock_free(&(eb->sub[i]));
239 void init_enxframe(t_enxframe *fr)
250 fr->nblock_alloc = 0;
255 void free_enxframe(t_enxframe *fr)
263 for (b = 0; b < fr->nblock_alloc; b++)
265 enxblock_free(&(fr->block[b]));
270 void add_blocks_enxframe(t_enxframe *fr, int n)
273 if (n > fr->nblock_alloc)
277 srenew(fr->block, n);
278 for (b = fr->nblock_alloc; b < fr->nblock; b++)
280 enxblock_init(&(fr->block[b]));
282 fr->nblock_alloc = n;
286 t_enxblock *find_block_id_enxframe(t_enxframe *ef, int id, t_enxblock *prev)
288 gmx_off_t starti = 0;
293 starti = (prev - ef->block) + 1;
295 for (i = starti; i < ef->nblock; i++)
297 if (ef->block[i].id == id)
299 return &(ef->block[i]);
305 void add_subblocks_enxblock(t_enxblock *eb, int n)
308 if (eb->nsub > eb->nsub_alloc)
313 for (b = eb->nsub_alloc; b < n; b++)
315 enxsubblock_init(&(eb->sub[b]));
321 static void enx_warning(const char *msg)
323 if (getenv("GMX_ENX_NO_FATAL") != NULL)
329 gmx_fatal(FARGS, "%s\n%s",
331 "If you want to use the correct frames before the corrupted frame and avoid this fatal error set the env.var. GMX_ENX_NO_FATAL");
335 static void edr_strings(XDR *xdr, gmx_bool bRead, int file_version,
336 int n, gmx_enxnm_t **nms)
345 for (i = 0; i < n; i++)
361 if (!xdr_string(xdr, &(nm->name), STRLEN))
363 gmx_file("Cannot write energy names to file; maybe you are out of disk space?");
365 if (file_version >= 2)
367 if (!xdr_string(xdr, &(nm->unit), STRLEN))
369 gmx_file("Cannot write energy names to file; maybe you are out of disk space?");
374 nm->unit = strdup("kJ/mol");
379 void do_enxnms(ener_file_t ef, int *nre, gmx_enxnm_t **nms)
383 gmx_bool bRead = gmx_fio_getread(ef->fio);
387 gmx_fio_checktype(ef->fio);
389 xdr = gmx_fio_getxdr(ef->fio);
391 if (!xdr_int(xdr, &magic))
395 gmx_file("Cannot write energy names to file; maybe you are out of disk space?");
402 /* Assume this is an old edr format */
405 ef->eo.bOldFileOpen = TRUE;
406 ef->eo.bReadFirstStep = FALSE;
407 srenew(ef->eo.ener_prev, *nre);
411 ef->eo.bOldFileOpen = FALSE;
415 gmx_fatal(FARGS, "Energy names magic number mismatch, this is not a GROMACS edr file");
417 file_version = enx_version;
418 xdr_int(xdr, &file_version);
419 if (file_version > enx_version)
421 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program", gmx_fio_getname(ef->fio), file_version, enx_version);
425 if (file_version != enx_version)
427 fprintf(stderr, "Note: enx file_version %d, software version %d\n",
428 file_version, enx_version);
431 edr_strings(xdr, bRead, file_version, *nre, nms);
434 static gmx_bool do_eheader(ener_file_t ef, int *file_version, t_enxframe *fr,
435 int nre_test, gmx_bool *bWrongPrecision, gmx_bool *bOK)
437 int magic = -7777777;
438 real first_real_to_check;
439 int b, i, zero = 0, dum = 0;
440 gmx_bool bRead = gmx_fio_getread(ef->fio);
445 xdr_datatype dtreal = xdr_datatype_float;
447 xdr_datatype dtreal = xdr_datatype_double;
452 *bWrongPrecision = FALSE;
456 /* The original energy frame started with a real,
457 * so we have to use a real for compatibility.
458 * This is VERY DIRTY code, since do_eheader can be called
459 * with the wrong precision set and then we could read r > -1e10,
460 * while actually the intention was r < -1e10.
461 * When nre_test >= 0, do_eheader should therefore terminate
462 * before the number of i/o calls starts depending on what has been read
463 * (which is the case for for instance the block sizes for variable
464 * number of blocks, where this number is read before).
466 first_real_to_check = -2e10;
467 if (!gmx_fio_do_real(ef->fio, first_real_to_check))
471 if (first_real_to_check > -1e10)
473 /* Assume we are reading an old format */
475 fr->t = first_real_to_check;
476 if (!gmx_fio_do_int(ef->fio, dum))
484 if (!gmx_fio_do_int(ef->fio, magic))
488 if (magic != -7777777)
490 enx_warning("Energy header magic number mismatch, this is not a GROMACS edr file");
494 *file_version = enx_version;
495 if (!gmx_fio_do_int(ef->fio, *file_version))
499 if (*bOK && *file_version > enx_version)
501 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program", gmx_fio_getname(ef->fio), file_version, enx_version);
503 if (!gmx_fio_do_double(ef->fio, fr->t))
507 if (!gmx_fio_do_int64(ef->fio, fr->step))
511 if (!bRead && fr->nsum == 1)
513 /* Do not store sums of length 1,
514 * since this does not add information.
516 if (!gmx_fio_do_int(ef->fio, zero))
523 if (!gmx_fio_do_int(ef->fio, fr->nsum))
528 if (*file_version >= 3)
530 if (!gmx_fio_do_int64(ef->fio, fr->nsteps))
537 fr->nsteps = max(1, fr->nsum);
539 if (*file_version >= 5)
541 if (!gmx_fio_do_double(ef->fio, fr->dt))
551 if (!gmx_fio_do_int(ef->fio, fr->nre))
555 if (*file_version < 4)
557 if (!gmx_fio_do_int(ef->fio, ndisre))
564 /* now reserved for possible future use */
565 if (!gmx_fio_do_int(ef->fio, dum))
571 if (!gmx_fio_do_int(ef->fio, fr->nblock))
582 if (*file_version >= 4)
584 enx_warning("Distance restraint blocks in old style in new style file");
592 /* Frames could have nre=0, so we can not rely only on the fr->nre check */
593 if (bRead && nre_test >= 0 &&
594 ((fr->nre > 0 && fr->nre != nre_test) ||
595 fr->nre < 0 || ndisre < 0 || fr->nblock < 0))
597 *bWrongPrecision = TRUE;
601 /* we now know what these should be, or we've already bailed out because
602 of wrong precision */
603 if (*file_version == 1 && (fr->t < 0 || fr->t > 1e20 || fr->step < 0 ) )
605 enx_warning("edr file with negative step number or unreasonable time (and without version number).");
613 add_blocks_enxframe(fr, fr->nblock);
619 /* sub[0] is the instantaneous data, sub[1] is time averaged */
620 add_subblocks_enxblock(&(fr->block[0]), 2);
621 fr->block[0].id = enxDISRE;
622 fr->block[0].sub[0].nr = ndisre;
623 fr->block[0].sub[1].nr = ndisre;
624 fr->block[0].sub[0].type = dtreal;
625 fr->block[0].sub[1].type = dtreal;
629 /* read block header info */
630 for (b = startb; b < fr->nblock; b++)
632 if (*file_version < 4)
634 /* blocks in old version files always have 1 subblock that
635 consists of reals. */
640 add_subblocks_enxblock(&(fr->block[b]), 1);
644 if (fr->block[b].nsub != 1)
646 gmx_incons("Writing an old version .edr file with too many subblocks");
648 if (fr->block[b].sub[0].type != dtreal)
650 gmx_incons("Writing an old version .edr file the wrong subblock type");
653 nrint = fr->block[b].sub[0].nr;
655 if (!gmx_fio_do_int(ef->fio, nrint))
659 fr->block[b].id = b - startb;
660 fr->block[b].sub[0].nr = nrint;
661 fr->block[b].sub[0].type = dtreal;
666 /* in the new version files, the block header only contains
667 the ID and the number of subblocks */
668 int nsub = fr->block[b].nsub;
669 *bOK = *bOK && gmx_fio_do_int(ef->fio, fr->block[b].id);
670 *bOK = *bOK && gmx_fio_do_int(ef->fio, nsub);
672 fr->block[b].nsub = nsub;
675 add_subblocks_enxblock(&(fr->block[b]), nsub);
678 /* read/write type & size for each subblock */
679 for (i = 0; i < nsub; i++)
681 t_enxsubblock *sub = &(fr->block[b].sub[i]); /* shortcut */
682 int typenr = sub->type;
684 *bOK = *bOK && gmx_fio_do_int(ef->fio, typenr);
685 *bOK = *bOK && gmx_fio_do_int(ef->fio, sub->nr);
687 sub->type = (xdr_datatype)typenr;
691 if (!gmx_fio_do_int(ef->fio, fr->e_size))
696 /* now reserved for possible future use */
697 if (!gmx_fio_do_int(ef->fio, dum))
702 /* Do a dummy int to keep the format compatible with the old code */
703 if (!gmx_fio_do_int(ef->fio, dum))
708 if (*bOK && *file_version == 1 && nre_test < 0)
710 if (!ef->eo.bReadFirstStep)
712 ef->eo.bReadFirstStep = TRUE;
713 ef->eo.first_step = fr->step;
714 ef->eo.step_prev = fr->step;
715 ef->eo.nsum_prev = 0;
718 fr->nsum = fr->step - ef->eo.first_step + 1;
719 fr->nsteps = fr->step - ef->eo.step_prev;
726 void free_enxnms(int n, gmx_enxnm_t *nms)
730 for (i = 0; i < n; i++)
739 void close_enx(ener_file_t ef)
741 if (gmx_fio_close(ef->fio) != 0)
743 gmx_file("Cannot close energy file; it might be corrupt, or maybe you are out of disk space?");
747 static gmx_bool empty_file(const char *fn)
754 fp = gmx_fio_fopen(fn, "r");
755 ret = fread(&dum, sizeof(dum), 1, fp);
763 ener_file_t open_enx(const char *fn, const char *mode)
766 gmx_enxnm_t *nms = NULL;
767 int file_version = -1;
769 gmx_bool bWrongPrecision, bOK = TRUE;
770 struct ener_file *ef;
776 ef->fio = gmx_fio_open(fn, mode);
777 gmx_fio_checktype(ef->fio);
778 gmx_fio_setprecision(ef->fio, FALSE);
779 do_enxnms(ef, &nre, &nms);
781 do_eheader(ef, &file_version, fr, nre, &bWrongPrecision, &bOK);
784 gmx_file("Cannot read energy file header. Corrupt file?");
787 /* Now check whether this file is in single precision */
788 if (!bWrongPrecision &&
789 ((fr->e_size && (fr->nre == nre) &&
790 (nre*4*(long int)sizeof(float) == fr->e_size)) ) )
792 fprintf(stderr, "Opened %s as single precision energy file\n", fn);
793 free_enxnms(nre, nms);
797 gmx_fio_rewind(ef->fio);
798 gmx_fio_checktype(ef->fio);
799 gmx_fio_setprecision(ef->fio, TRUE);
800 do_enxnms(ef, &nre, &nms);
801 do_eheader(ef, &file_version, fr, nre, &bWrongPrecision, &bOK);
804 gmx_file("Cannot write energy file header; maybe you are out of disk space?");
807 if (((fr->e_size && (fr->nre == nre) &&
808 (nre*4*(long int)sizeof(double) == fr->e_size)) ))
810 fprintf(stderr, "Opened %s as double precision energy file\n",
817 gmx_fatal(FARGS, "File %s is empty", fn);
821 gmx_fatal(FARGS, "Energy file %s not recognized, maybe different CPU?",
825 free_enxnms(nre, nms);
829 gmx_fio_rewind(ef->fio);
833 ef->fio = gmx_fio_open(fn, mode);
841 t_fileio *enx_file_pointer(const ener_file_t ef)
846 static void convert_full_sums(ener_old_t *ener_old, t_enxframe *fr)
850 double esum_all, eav_all;
856 for (i = 0; i < fr->nre; i++)
858 if (fr->ener[i].e != 0)
862 if (fr->ener[i].esum != 0)
867 if (ne > 0 && ns == 0)
869 /* We do not have all energy sums */
874 /* Convert old full simulation sums to sums between energy frames */
875 nstep_all = fr->step - ener_old->first_step + 1;
876 if (fr->nsum > 1 && fr->nsum == nstep_all && ener_old->nsum_prev > 0)
878 /* Set the new sum length: the frame step difference */
879 fr->nsum = fr->step - ener_old->step_prev;
880 for (i = 0; i < fr->nre; i++)
882 esum_all = fr->ener[i].esum;
883 eav_all = fr->ener[i].eav;
884 fr->ener[i].esum = esum_all - ener_old->ener_prev[i].esum;
885 fr->ener[i].eav = eav_all - ener_old->ener_prev[i].eav
886 - dsqr(ener_old->ener_prev[i].esum/(nstep_all - fr->nsum)
887 - esum_all/nstep_all)*
888 (nstep_all - fr->nsum)*nstep_all/(double)fr->nsum;
889 ener_old->ener_prev[i].esum = esum_all;
890 ener_old->ener_prev[i].eav = eav_all;
892 ener_old->nsum_prev = nstep_all;
894 else if (fr->nsum > 0)
896 if (fr->nsum != nstep_all)
898 fprintf(stderr, "\nWARNING: something is wrong with the energy sums, will not use exact averages\n");
899 ener_old->nsum_prev = 0;
903 ener_old->nsum_prev = nstep_all;
905 /* Copy all sums to ener_prev */
906 for (i = 0; i < fr->nre; i++)
908 ener_old->ener_prev[i].esum = fr->ener[i].esum;
909 ener_old->ener_prev[i].eav = fr->ener[i].eav;
913 ener_old->step_prev = fr->step;
916 gmx_bool do_enx(ener_file_t ef, t_enxframe *fr)
918 int file_version = -1;
920 gmx_bool bRead, bOK, bOK1, bSane;
921 real tmp1, tmp2, rdum;
926 bRead = gmx_fio_getread(ef->fio);
929 fr->e_size = fr->nre*sizeof(fr->ener[0].e)*4;
930 /*d_size = fr->ndisre*(sizeof(real)*2);*/
932 gmx_fio_checktype(ef->fio);
934 if (!do_eheader(ef, &file_version, fr, -1, NULL, &bOK))
938 fprintf(stderr, "\rLast energy frame read %d time %8.3f ",
939 ef->framenr-1, ef->frametime);
943 "\nWARNING: Incomplete energy frame: nr %d time %8.3f\n",
949 gmx_file("Cannot write energy file header; maybe you are out of disk space?");
955 if ((ef->framenr < 20 || ef->framenr % 10 == 0) &&
956 (ef->framenr < 200 || ef->framenr % 100 == 0) &&
957 (ef->framenr < 2000 || ef->framenr % 1000 == 0))
959 fprintf(stderr, "\rReading energy frame %6d time %8.3f ",
963 ef->frametime = fr->t;
965 /* Check sanity of this header */
967 for (b = 0; b < fr->nblock; b++)
969 bSane = bSane || (fr->block[b].nsub > 0);
971 if (!((fr->step >= 0) && bSane))
973 fprintf(stderr, "\nWARNING: there may be something wrong with energy file %s\n",
974 gmx_fio_getname(ef->fio));
975 fprintf(stderr, "Found: step=%s, nre=%d, nblock=%d, time=%g.\n"
976 "Trying to skip frame expect a crash though\n",
977 gmx_step_str(fr->step, buf), fr->nre, fr->nblock, fr->t);
979 if (bRead && fr->nre > fr->e_alloc)
981 srenew(fr->ener, fr->nre);
982 for (i = fr->e_alloc; (i < fr->nre); i++)
986 fr->ener[i].esum = 0;
988 fr->e_alloc = fr->nre;
991 for (i = 0; i < fr->nre; i++)
993 bOK = bOK && gmx_fio_do_real(ef->fio, fr->ener[i].e);
995 /* Do not store sums of length 1,
996 * since this does not add information.
998 if (file_version == 1 ||
999 (bRead && fr->nsum > 0) || fr->nsum > 1)
1001 tmp1 = fr->ener[i].eav;
1002 bOK = bOK && gmx_fio_do_real(ef->fio, tmp1);
1005 fr->ener[i].eav = tmp1;
1008 /* This is to save only in single precision (unless compiled in DP) */
1009 tmp2 = fr->ener[i].esum;
1010 bOK = bOK && gmx_fio_do_real(ef->fio, tmp2);
1013 fr->ener[i].esum = tmp2;
1016 if (file_version == 1)
1018 /* Old, unused real */
1020 bOK = bOK && gmx_fio_do_real(ef->fio, rdum);
1025 /* Here we can not check for file_version==1, since one could have
1026 * continued an old format simulation with a new one with mdrun -append.
1028 if (bRead && ef->eo.bOldFileOpen)
1030 /* Convert old full simulation sums to sums between energy frames */
1031 convert_full_sums(&(ef->eo), fr);
1033 /* read the blocks */
1034 for (b = 0; b < fr->nblock; b++)
1036 /* now read the subblocks. */
1037 int nsub = fr->block[b].nsub; /* shortcut */
1040 for (i = 0; i < nsub; i++)
1042 t_enxsubblock *sub = &(fr->block[b].sub[i]); /* shortcut */
1046 enxsubblock_alloc(sub);
1049 /* read/write data */
1053 case xdr_datatype_float:
1054 bOK1 = gmx_fio_ndo_float(ef->fio, sub->fval, sub->nr);
1056 case xdr_datatype_double:
1057 bOK1 = gmx_fio_ndo_double(ef->fio, sub->dval, sub->nr);
1059 case xdr_datatype_int:
1060 bOK1 = gmx_fio_ndo_int(ef->fio, sub->ival, sub->nr);
1062 case xdr_datatype_int64:
1063 bOK1 = gmx_fio_ndo_int64(ef->fio, sub->lval, sub->nr);
1065 case xdr_datatype_char:
1066 bOK1 = gmx_fio_ndo_uchar(ef->fio, sub->cval, sub->nr);
1068 case xdr_datatype_string:
1069 bOK1 = gmx_fio_ndo_string(ef->fio, sub->sval, sub->nr);
1072 gmx_incons("Reading unknown block data type: this file is corrupted or from the future");
1080 if (gmx_fio_flush(ef->fio) != 0)
1082 gmx_file("Cannot write energy file; maybe you are out of disk space?");
1090 fprintf(stderr, "\nLast energy frame read %d",
1092 fprintf(stderr, "\nWARNING: Incomplete energy frame: nr %d time %8.3f\n",
1093 ef->framenr, fr->t);
1097 gmx_fatal(FARGS, "could not write energies");
1105 static real find_energy(const char *name, int nre, gmx_enxnm_t *enm,
1110 for (i = 0; i < nre; i++)
1112 if (strcmp(enm[i].name, name) == 0)
1114 return fr->ener[i].e;
1118 gmx_fatal(FARGS, "Could not find energy term named '%s'", name);
1124 void get_enx_state(const char *fn, real t, gmx_groups_t *groups, t_inputrec *ir,
1127 /* Should match the names in mdebin.c */
1128 static const char *boxvel_nm[] = {
1129 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
1130 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
1133 static const char *pcouplmu_nm[] = {
1134 "Pcoupl-Mu-XX", "Pcoupl-Mu-YY", "Pcoupl-Mu-ZZ",
1135 "Pcoupl-Mu-YX", "Pcoupl-Mu-ZX", "Pcoupl-Mu-ZY"
1137 static const char *baro_nm[] = {
1142 int ind0[] = { XX, YY, ZZ, YY, ZZ, ZZ };
1143 int ind1[] = { XX, YY, ZZ, XX, XX, YY };
1144 int nre, nfr, i, j, ni, npcoupl;
1147 gmx_enxnm_t *enm = NULL;
1151 in = open_enx(fn, "r");
1152 do_enxnms(in, &nre, &enm);
1155 while ((nfr == 0 || fr->t != t) && do_enx(in, fr))
1160 fprintf(stderr, "\n");
1162 if (nfr == 0 || fr->t != t)
1164 gmx_fatal(FARGS, "Could not find frame with time %f in '%s'", t, fn);
1167 npcoupl = TRICLINIC(ir->compress) ? 6 : 3;
1168 if (ir->epc == epcPARRINELLORAHMAN)
1170 clear_mat(state->boxv);
1171 for (i = 0; i < npcoupl; i++)
1173 state->boxv[ind0[i]][ind1[i]] =
1174 find_energy(boxvel_nm[i], nre, enm, fr);
1176 fprintf(stderr, "\nREAD %d BOX VELOCITIES FROM %s\n\n", npcoupl, fn);
1179 if (ir->etc == etcNOSEHOOVER)
1185 for (i = 0; i < state->ngtc; i++)
1187 ni = groups->grps[egcTC].nm_ind[i];
1188 bufi = *(groups->grpname[ni]);
1189 for (j = 0; (j < state->nhchainlength); j++)
1191 if (IR_NVT_TROTTER(ir))
1193 sprintf(cns, "-%d", j);
1195 sprintf(buf, "Xi%s-%s", cns, bufi);
1196 state->nosehoover_xi[i] = find_energy(buf, nre, enm, fr);
1197 sprintf(buf, "vXi%s-%s", cns, bufi);
1198 state->nosehoover_vxi[i] = find_energy(buf, nre, enm, fr);
1202 fprintf(stderr, "\nREAD %d NOSE-HOOVER Xi chains FROM %s\n\n", state->ngtc, fn);
1204 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir))
1206 for (i = 0; i < state->nnhpres; i++)
1208 bufi = baro_nm[0]; /* All barostat DOF's together for now */
1209 for (j = 0; (j < state->nhchainlength); j++)
1211 sprintf(buf, "Xi-%d-%s", j, bufi);
1212 state->nhpres_xi[i] = find_energy(buf, nre, enm, fr);
1213 sprintf(buf, "vXi-%d-%s", j, bufi);
1214 state->nhpres_vxi[i] = find_energy(buf, nre, enm, fr);
1217 fprintf(stderr, "\nREAD %d NOSE-HOOVER BAROSTAT Xi chains FROM %s\n\n", state->nnhpres, fn);
1221 free_enxnms(nre, enm);