1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
42 #include "gmx_fatal.h"
50 /* The source code in this file should be thread-safe.
51 Please keep it that way. */
53 /* This number should be increased whenever the file format changes! */
54 static const int enx_version = 5;
56 const char *enx_block_id_name[] = {
57 "Averaged orientation restraints",
58 "Instantaneous orientation restraints",
59 "Orientation restraint order tensor(s)",
60 "Distance restraints",
67 /* Stuff for reading pre 4.1 energy files */
69 gmx_bool bOldFileOpen; /* Is this an open old file? */
70 gmx_bool bReadFirstStep; /* Did we read the first step? */
71 int first_step; /* First step in the energy file */
72 int step_prev; /* Previous step */
73 int nsum_prev; /* Previous step sum length */
74 t_energy *ener_prev; /* Previous energy sums */
85 static void enxsubblock_init(t_enxsubblock *sb)
89 sb->type = xdr_datatype_double;
91 sb->type = xdr_datatype_float;
107 static void enxsubblock_free(t_enxsubblock *sb)
143 for (i = 0; i < sb->sval_alloc; i++)
156 /* allocate the appropriate amount of memory for the given type and nr */
157 static void enxsubblock_alloc(t_enxsubblock *sb)
159 /* allocate the appropriate amount of memory */
162 case xdr_datatype_float:
163 if (sb->nr > sb->fval_alloc)
165 srenew(sb->fval, sb->nr);
166 sb->fval_alloc = sb->nr;
169 case xdr_datatype_double:
170 if (sb->nr > sb->dval_alloc)
172 srenew(sb->dval, sb->nr);
173 sb->dval_alloc = sb->nr;
176 case xdr_datatype_int:
177 if (sb->nr > sb->ival_alloc)
179 srenew(sb->ival, sb->nr);
180 sb->ival_alloc = sb->nr;
183 case xdr_datatype_large_int:
184 if (sb->nr > sb->lval_alloc)
186 srenew(sb->lval, sb->nr);
187 sb->lval_alloc = sb->nr;
190 case xdr_datatype_char:
191 if (sb->nr > sb->cval_alloc)
193 srenew(sb->cval, sb->nr);
194 sb->cval_alloc = sb->nr;
197 case xdr_datatype_string:
198 if (sb->nr > sb->sval_alloc)
202 srenew(sb->sval, sb->nr);
203 for (i = sb->sval_alloc; i < sb->nr; i++)
207 sb->sval_alloc = sb->nr;
211 gmx_incons("Unknown block type: this file is corrupted or from the future");
215 static void enxblock_init(t_enxblock *eb)
223 static void enxblock_free(t_enxblock *eb)
225 if (eb->nsub_alloc > 0)
228 for (i = 0; i < eb->nsub_alloc; i++)
230 enxsubblock_free(&(eb->sub[i]));
238 void init_enxframe(t_enxframe *fr)
249 fr->nblock_alloc = 0;
254 void free_enxframe(t_enxframe *fr)
262 for (b = 0; b < fr->nblock_alloc; b++)
264 enxblock_free(&(fr->block[b]));
269 void add_blocks_enxframe(t_enxframe *fr, int n)
272 if (n > fr->nblock_alloc)
276 srenew(fr->block, n);
277 for (b = fr->nblock_alloc; b < fr->nblock; b++)
279 enxblock_init(&(fr->block[b]));
281 fr->nblock_alloc = n;
285 t_enxblock *find_block_id_enxframe(t_enxframe *ef, int id, t_enxblock *prev)
287 gmx_off_t starti = 0;
292 starti = (prev - ef->block) + 1;
294 for (i = starti; i < ef->nblock; i++)
296 if (ef->block[i].id == id)
298 return &(ef->block[i]);
304 void add_subblocks_enxblock(t_enxblock *eb, int n)
307 if (eb->nsub > eb->nsub_alloc)
312 for (b = eb->nsub_alloc; b < n; b++)
314 enxsubblock_init(&(eb->sub[b]));
320 static void enx_warning(const char *msg)
322 if (getenv("GMX_ENX_NO_FATAL") != NULL)
328 gmx_fatal(FARGS, "%s\n%s",
330 "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");
334 static void edr_strings(XDR *xdr, gmx_bool bRead, int file_version,
335 int n, gmx_enxnm_t **nms)
344 for (i = 0; i < n; i++)
360 if (!xdr_string(xdr, &(nm->name), STRLEN))
362 gmx_file("Cannot write energy names to file; maybe you are out of disk space?");
364 if (file_version >= 2)
366 if (!xdr_string(xdr, &(nm->unit), STRLEN))
368 gmx_file("Cannot write energy names to file; maybe you are out of disk space?");
373 nm->unit = strdup("kJ/mol");
378 void do_enxnms(ener_file_t ef, int *nre, gmx_enxnm_t **nms)
382 gmx_bool bRead = gmx_fio_getread(ef->fio);
386 gmx_fio_checktype(ef->fio);
388 xdr = gmx_fio_getxdr(ef->fio);
390 if (!xdr_int(xdr, &magic))
394 gmx_file("Cannot write energy names to file; maybe you are out of disk space?");
401 /* Assume this is an old edr format */
404 ef->eo.bOldFileOpen = TRUE;
405 ef->eo.bReadFirstStep = FALSE;
406 srenew(ef->eo.ener_prev, *nre);
410 ef->eo.bOldFileOpen = FALSE;
414 gmx_fatal(FARGS, "Energy names magic number mismatch, this is not a GROMACS edr file");
416 file_version = enx_version;
417 xdr_int(xdr, &file_version);
418 if (file_version > enx_version)
420 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program", gmx_fio_getname(ef->fio), file_version, enx_version);
424 if (file_version != enx_version)
426 fprintf(stderr, "Note: enx file_version %d, software version %d\n",
427 file_version, enx_version);
430 edr_strings(xdr, bRead, file_version, *nre, nms);
433 static gmx_bool do_eheader(ener_file_t ef, int *file_version, t_enxframe *fr,
434 int nre_test, gmx_bool *bWrongPrecision, gmx_bool *bOK)
436 int magic = -7777777;
437 real first_real_to_check;
438 int b, i, zero = 0, dum = 0;
439 gmx_bool bRead = gmx_fio_getread(ef->fio);
444 xdr_datatype dtreal = xdr_datatype_float;
446 xdr_datatype dtreal = xdr_datatype_double;
451 *bWrongPrecision = FALSE;
455 /* The original energy frame started with a real,
456 * so we have to use a real for compatibility.
457 * This is VERY DIRTY code, since do_eheader can be called
458 * with the wrong precision set and then we could read r > -1e10,
459 * while actually the intention was r < -1e10.
460 * When nre_test >= 0, do_eheader should therefore terminate
461 * before the number of i/o calls starts depending on what has been read
462 * (which is the case for for instance the block sizes for variable
463 * number of blocks, where this number is read before).
465 first_real_to_check = -2e10;
466 if (!gmx_fio_do_real(ef->fio, first_real_to_check))
470 if (first_real_to_check > -1e10)
472 /* Assume we are reading an old format */
474 fr->t = first_real_to_check;
475 if (!gmx_fio_do_int(ef->fio, dum))
483 if (!gmx_fio_do_int(ef->fio, magic))
487 if (magic != -7777777)
489 enx_warning("Energy header magic number mismatch, this is not a GROMACS edr file");
493 *file_version = enx_version;
494 if (!gmx_fio_do_int(ef->fio, *file_version))
498 if (*bOK && *file_version > enx_version)
500 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program", gmx_fio_getname(ef->fio), file_version, enx_version);
502 if (!gmx_fio_do_double(ef->fio, fr->t))
506 if (!gmx_fio_do_gmx_large_int(ef->fio, fr->step))
510 if (!bRead && fr->nsum == 1)
512 /* Do not store sums of length 1,
513 * since this does not add information.
515 if (!gmx_fio_do_int(ef->fio, zero))
522 if (!gmx_fio_do_int(ef->fio, fr->nsum))
527 if (*file_version >= 3)
529 if (!gmx_fio_do_gmx_large_int(ef->fio, fr->nsteps))
536 fr->nsteps = max(1, fr->nsum);
538 if (*file_version >= 5)
540 if (!gmx_fio_do_double(ef->fio, fr->dt))
550 if (!gmx_fio_do_int(ef->fio, fr->nre))
554 if (*file_version < 4)
556 if (!gmx_fio_do_int(ef->fio, ndisre))
563 /* now reserved for possible future use */
564 if (!gmx_fio_do_int(ef->fio, dum))
570 if (!gmx_fio_do_int(ef->fio, fr->nblock))
581 if (*file_version >= 4)
583 enx_warning("Distance restraint blocks in old style in new style file");
591 /* Frames could have nre=0, so we can not rely only on the fr->nre check */
592 if (bRead && nre_test >= 0 &&
593 ((fr->nre > 0 && fr->nre != nre_test) ||
594 fr->nre < 0 || ndisre < 0 || fr->nblock < 0))
596 *bWrongPrecision = TRUE;
600 /* we now know what these should be, or we've already bailed out because
601 of wrong precision */
602 if (*file_version == 1 && (fr->t < 0 || fr->t > 1e20 || fr->step < 0 ) )
604 enx_warning("edr file with negative step number or unreasonable time (and without version number).");
612 add_blocks_enxframe(fr, fr->nblock);
618 /* sub[0] is the instantaneous data, sub[1] is time averaged */
619 add_subblocks_enxblock(&(fr->block[0]), 2);
620 fr->block[0].id = enxDISRE;
621 fr->block[0].sub[0].nr = ndisre;
622 fr->block[0].sub[1].nr = ndisre;
623 fr->block[0].sub[0].type = dtreal;
624 fr->block[0].sub[1].type = dtreal;
628 /* read block header info */
629 for (b = startb; b < fr->nblock; b++)
631 if (*file_version < 4)
633 /* blocks in old version files always have 1 subblock that
634 consists of reals. */
639 add_subblocks_enxblock(&(fr->block[b]), 1);
643 if (fr->block[b].nsub != 1)
645 gmx_incons("Writing an old version .edr file with too many subblocks");
647 if (fr->block[b].sub[0].type != dtreal)
649 gmx_incons("Writing an old version .edr file the wrong subblock type");
652 nrint = fr->block[b].sub[0].nr;
654 if (!gmx_fio_do_int(ef->fio, nrint))
658 fr->block[b].id = b - startb;
659 fr->block[b].sub[0].nr = nrint;
660 fr->block[b].sub[0].type = dtreal;
665 /* in the new version files, the block header only contains
666 the ID and the number of subblocks */
667 int nsub = fr->block[b].nsub;
668 *bOK = *bOK && gmx_fio_do_int(ef->fio, fr->block[b].id);
669 *bOK = *bOK && gmx_fio_do_int(ef->fio, nsub);
671 fr->block[b].nsub = nsub;
674 add_subblocks_enxblock(&(fr->block[b]), nsub);
677 /* read/write type & size for each subblock */
678 for (i = 0; i < nsub; i++)
680 t_enxsubblock *sub = &(fr->block[b].sub[i]); /* shortcut */
681 int typenr = sub->type;
683 *bOK = *bOK && gmx_fio_do_int(ef->fio, typenr);
684 *bOK = *bOK && gmx_fio_do_int(ef->fio, sub->nr);
686 sub->type = (xdr_datatype)typenr;
690 if (!gmx_fio_do_int(ef->fio, fr->e_size))
695 /* now reserved for possible future use */
696 if (!gmx_fio_do_int(ef->fio, dum))
701 /* Do a dummy int to keep the format compatible with the old code */
702 if (!gmx_fio_do_int(ef->fio, dum))
707 if (*bOK && *file_version == 1 && nre_test < 0)
710 if (fp >= ener_old_nalloc)
712 gmx_incons("Problem with reading old format energy files");
716 if (!ef->eo.bReadFirstStep)
718 ef->eo.bReadFirstStep = TRUE;
719 ef->eo.first_step = fr->step;
720 ef->eo.step_prev = fr->step;
721 ef->eo.nsum_prev = 0;
724 fr->nsum = fr->step - ef->eo.first_step + 1;
725 fr->nsteps = fr->step - ef->eo.step_prev;
732 void free_enxnms(int n, gmx_enxnm_t *nms)
736 for (i = 0; i < n; i++)
745 void close_enx(ener_file_t ef)
747 if (gmx_fio_close(ef->fio) != 0)
749 gmx_file("Cannot close energy file; it might be corrupt, or maybe you are out of disk space?");
753 static gmx_bool empty_file(const char *fn)
760 fp = gmx_fio_fopen(fn, "r");
761 ret = fread(&dum, sizeof(dum), 1, fp);
769 ener_file_t open_enx(const char *fn, const char *mode)
772 gmx_enxnm_t *nms = NULL;
773 int file_version = -1;
775 gmx_bool bWrongPrecision, bOK = TRUE;
776 struct ener_file *ef;
782 ef->fio = gmx_fio_open(fn, mode);
783 gmx_fio_checktype(ef->fio);
784 gmx_fio_setprecision(ef->fio, FALSE);
785 do_enxnms(ef, &nre, &nms);
787 do_eheader(ef, &file_version, fr, nre, &bWrongPrecision, &bOK);
790 gmx_file("Cannot read energy file header. Corrupt file?");
793 /* Now check whether this file is in single precision */
794 if (!bWrongPrecision &&
795 ((fr->e_size && (fr->nre == nre) &&
796 (nre*4*(long int)sizeof(float) == fr->e_size)) ) )
798 fprintf(stderr, "Opened %s as single precision energy file\n", fn);
799 free_enxnms(nre, nms);
803 gmx_fio_rewind(ef->fio);
804 gmx_fio_checktype(ef->fio);
805 gmx_fio_setprecision(ef->fio, TRUE);
806 do_enxnms(ef, &nre, &nms);
807 do_eheader(ef, &file_version, fr, nre, &bWrongPrecision, &bOK);
810 gmx_file("Cannot write energy file header; maybe you are out of disk space?");
813 if (((fr->e_size && (fr->nre == nre) &&
814 (nre*4*(long int)sizeof(double) == fr->e_size)) ))
816 fprintf(stderr, "Opened %s as double precision energy file\n",
823 gmx_fatal(FARGS, "File %s is empty", fn);
827 gmx_fatal(FARGS, "Energy file %s not recognized, maybe different CPU?",
831 free_enxnms(nre, nms);
835 gmx_fio_rewind(ef->fio);
839 ef->fio = gmx_fio_open(fn, mode);
847 t_fileio *enx_file_pointer(const ener_file_t ef)
852 static void convert_full_sums(ener_old_t *ener_old, t_enxframe *fr)
856 double esum_all, eav_all;
862 for (i = 0; i < fr->nre; i++)
864 if (fr->ener[i].e != 0)
868 if (fr->ener[i].esum != 0)
873 if (ne > 0 && ns == 0)
875 /* We do not have all energy sums */
880 /* Convert old full simulation sums to sums between energy frames */
881 nstep_all = fr->step - ener_old->first_step + 1;
882 if (fr->nsum > 1 && fr->nsum == nstep_all && ener_old->nsum_prev > 0)
884 /* Set the new sum length: the frame step difference */
885 fr->nsum = fr->step - ener_old->step_prev;
886 for (i = 0; i < fr->nre; i++)
888 esum_all = fr->ener[i].esum;
889 eav_all = fr->ener[i].eav;
890 fr->ener[i].esum = esum_all - ener_old->ener_prev[i].esum;
891 fr->ener[i].eav = eav_all - ener_old->ener_prev[i].eav
892 - dsqr(ener_old->ener_prev[i].esum/(nstep_all - fr->nsum)
893 - esum_all/nstep_all)*
894 (nstep_all - fr->nsum)*nstep_all/(double)fr->nsum;
895 ener_old->ener_prev[i].esum = esum_all;
896 ener_old->ener_prev[i].eav = eav_all;
898 ener_old->nsum_prev = nstep_all;
900 else if (fr->nsum > 0)
902 if (fr->nsum != nstep_all)
904 fprintf(stderr, "\nWARNING: something is wrong with the energy sums, will not use exact averages\n");
905 ener_old->nsum_prev = 0;
909 ener_old->nsum_prev = nstep_all;
911 /* Copy all sums to ener_prev */
912 for (i = 0; i < fr->nre; i++)
914 ener_old->ener_prev[i].esum = fr->ener[i].esum;
915 ener_old->ener_prev[i].eav = fr->ener[i].eav;
919 ener_old->step_prev = fr->step;
922 gmx_bool do_enx(ener_file_t ef, t_enxframe *fr)
924 int file_version = -1;
926 gmx_bool bRead, bOK, bOK1, bSane;
927 real tmp1, tmp2, rdum;
932 bRead = gmx_fio_getread(ef->fio);
935 fr->e_size = fr->nre*sizeof(fr->ener[0].e)*4;
936 /*d_size = fr->ndisre*(sizeof(real)*2);*/
938 gmx_fio_checktype(ef->fio);
940 if (!do_eheader(ef, &file_version, fr, -1, NULL, &bOK))
944 fprintf(stderr, "\rLast energy frame read %d time %8.3f ",
945 ef->framenr-1, ef->frametime);
949 "\nWARNING: Incomplete energy frame: nr %d time %8.3f\n",
955 gmx_file("Cannot write energy file header; maybe you are out of disk space?");
961 if ((ef->framenr < 20 || ef->framenr % 10 == 0) &&
962 (ef->framenr < 200 || ef->framenr % 100 == 0) &&
963 (ef->framenr < 2000 || ef->framenr % 1000 == 0))
965 fprintf(stderr, "\rReading energy frame %6d time %8.3f ",
969 ef->frametime = fr->t;
971 /* Check sanity of this header */
973 for (b = 0; b < fr->nblock; b++)
975 bSane = bSane || (fr->block[b].nsub > 0);
977 if (!((fr->step >= 0) && bSane))
979 fprintf(stderr, "\nWARNING: there may be something wrong with energy file %s\n",
980 gmx_fio_getname(ef->fio));
981 fprintf(stderr, "Found: step=%s, nre=%d, nblock=%d, time=%g.\n"
982 "Trying to skip frame expect a crash though\n",
983 gmx_step_str(fr->step, buf), fr->nre, fr->nblock, fr->t);
985 if (bRead && fr->nre > fr->e_alloc)
987 srenew(fr->ener, fr->nre);
988 for (i = fr->e_alloc; (i < fr->nre); i++)
992 fr->ener[i].esum = 0;
994 fr->e_alloc = fr->nre;
997 for (i = 0; i < fr->nre; i++)
999 bOK = bOK && gmx_fio_do_real(ef->fio, fr->ener[i].e);
1001 /* Do not store sums of length 1,
1002 * since this does not add information.
1004 if (file_version == 1 ||
1005 (bRead && fr->nsum > 0) || fr->nsum > 1)
1007 tmp1 = fr->ener[i].eav;
1008 bOK = bOK && gmx_fio_do_real(ef->fio, tmp1);
1011 fr->ener[i].eav = tmp1;
1014 /* This is to save only in single precision (unless compiled in DP) */
1015 tmp2 = fr->ener[i].esum;
1016 bOK = bOK && gmx_fio_do_real(ef->fio, tmp2);
1019 fr->ener[i].esum = tmp2;
1022 if (file_version == 1)
1024 /* Old, unused real */
1026 bOK = bOK && gmx_fio_do_real(ef->fio, rdum);
1031 /* Here we can not check for file_version==1, since one could have
1032 * continued an old format simulation with a new one with mdrun -append.
1034 if (bRead && ef->eo.bOldFileOpen)
1036 /* Convert old full simulation sums to sums between energy frames */
1037 convert_full_sums(&(ef->eo), fr);
1039 /* read the blocks */
1040 for (b = 0; b < fr->nblock; b++)
1042 /* now read the subblocks. */
1043 int nsub = fr->block[b].nsub; /* shortcut */
1046 for (i = 0; i < nsub; i++)
1048 t_enxsubblock *sub = &(fr->block[b].sub[i]); /* shortcut */
1052 enxsubblock_alloc(sub);
1055 /* read/write data */
1059 case xdr_datatype_float:
1060 bOK1 = gmx_fio_ndo_float(ef->fio, sub->fval, sub->nr);
1062 case xdr_datatype_double:
1063 bOK1 = gmx_fio_ndo_double(ef->fio, sub->dval, sub->nr);
1065 case xdr_datatype_int:
1066 bOK1 = gmx_fio_ndo_int(ef->fio, sub->ival, sub->nr);
1068 case xdr_datatype_large_int:
1069 bOK1 = gmx_fio_ndo_gmx_large_int(ef->fio, sub->lval, sub->nr);
1071 case xdr_datatype_char:
1072 bOK1 = gmx_fio_ndo_uchar(ef->fio, sub->cval, sub->nr);
1074 case xdr_datatype_string:
1075 bOK1 = gmx_fio_ndo_string(ef->fio, sub->sval, sub->nr);
1078 gmx_incons("Reading unknown block data type: this file is corrupted or from the future");
1086 if (gmx_fio_flush(ef->fio) != 0)
1088 gmx_file("Cannot write energy file; maybe you are out of disk space?");
1096 fprintf(stderr, "\nLast energy frame read %d",
1098 fprintf(stderr, "\nWARNING: Incomplete energy frame: nr %d time %8.3f\n",
1099 ef->framenr, fr->t);
1103 gmx_fatal(FARGS, "could not write energies");
1111 static real find_energy(const char *name, int nre, gmx_enxnm_t *enm,
1116 for (i = 0; i < nre; i++)
1118 if (strcmp(enm[i].name, name) == 0)
1120 return fr->ener[i].e;
1124 gmx_fatal(FARGS, "Could not find energy term named '%s'", name);
1130 void get_enx_state(const char *fn, real t, gmx_groups_t *groups, t_inputrec *ir,
1133 /* Should match the names in mdebin.c */
1134 static const char *boxvel_nm[] = {
1135 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
1136 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
1139 static const char *pcouplmu_nm[] = {
1140 "Pcoupl-Mu-XX", "Pcoupl-Mu-YY", "Pcoupl-Mu-ZZ",
1141 "Pcoupl-Mu-YX", "Pcoupl-Mu-ZX", "Pcoupl-Mu-ZY"
1143 static const char *baro_nm[] = {
1148 int ind0[] = { XX, YY, ZZ, YY, ZZ, ZZ };
1149 int ind1[] = { XX, YY, ZZ, XX, XX, YY };
1150 int nre, nfr, i, j, ni, npcoupl;
1153 gmx_enxnm_t *enm = NULL;
1157 in = open_enx(fn, "r");
1158 do_enxnms(in, &nre, &enm);
1161 while ((nfr == 0 || fr->t != t) && do_enx(in, fr))
1166 fprintf(stderr, "\n");
1168 if (nfr == 0 || fr->t != t)
1170 gmx_fatal(FARGS, "Could not find frame with time %f in '%s'", t, fn);
1173 npcoupl = TRICLINIC(ir->compress) ? 6 : 3;
1174 if (ir->epc == epcPARRINELLORAHMAN)
1176 clear_mat(state->boxv);
1177 for (i = 0; i < npcoupl; i++)
1179 state->boxv[ind0[i]][ind1[i]] =
1180 find_energy(boxvel_nm[i], nre, enm, fr);
1182 fprintf(stderr, "\nREAD %d BOX VELOCITIES FROM %s\n\n", npcoupl, fn);
1185 if (ir->etc == etcNOSEHOOVER)
1191 for (i = 0; i < state->ngtc; i++)
1193 ni = groups->grps[egcTC].nm_ind[i];
1194 bufi = *(groups->grpname[ni]);
1195 for (j = 0; (j < state->nhchainlength); j++)
1197 if (IR_NVT_TROTTER(ir))
1199 sprintf(cns, "-%d", j);
1201 sprintf(buf, "Xi%s-%s", cns, bufi);
1202 state->nosehoover_xi[i] = find_energy(buf, nre, enm, fr);
1203 sprintf(buf, "vXi%s-%s", cns, bufi);
1204 state->nosehoover_vxi[i] = find_energy(buf, nre, enm, fr);
1208 fprintf(stderr, "\nREAD %d NOSE-HOOVER Xi chains FROM %s\n\n", state->ngtc, fn);
1210 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir))
1212 for (i = 0; i < state->nnhpres; i++)
1214 bufi = baro_nm[0]; /* All barostat DOF's together for now */
1215 for (j = 0; (j < state->nhchainlength); j++)
1217 sprintf(buf, "Xi-%d-%s", j, bufi);
1218 state->nhpres_xi[i] = find_energy(buf, nre, enm, fr);
1219 sprintf(buf, "vXi-%d-%s", j, bufi);
1220 state->nhpres_vxi[i] = find_energy(buf, nre, enm, fr);
1223 fprintf(stderr, "\nREAD %d NOSE-HOOVER BAROSTAT Xi chains FROM %s\n\n", state->nnhpres, fn);
1227 free_enxnms(nre, enm);