Merge release-4-6 into release-5-0
[gromacs.git] / src / gromacs / gmxlib / checkpoint.c
blob825a8968dad1cbb59ef6e5b94ce8b750b8019395
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2008,2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 /* The source code in this file should be thread-safe.
37 Please keep it that way. */
40 #ifdef HAVE_CONFIG_H
41 #include <config.h>
42 #endif
44 #include <string.h>
45 #include <time.h>
47 #ifdef HAVE_SYS_TIME_H
48 #include <sys/time.h>
49 #endif
51 #ifdef HAVE_UNISTD_H
52 #include <unistd.h>
53 #endif
55 #ifdef GMX_NATIVE_WINDOWS
56 /* _chsize_s */
57 #include <io.h>
58 #include <sys/locking.h>
59 #endif
61 #include "copyrite.h"
62 #include "names.h"
63 #include "typedefs.h"
64 #include "types/commrec.h"
65 #include "gromacs/utility/smalloc.h"
66 #include "txtdump.h"
67 #include "vec.h"
68 #include "network.h"
69 #include "checkpoint.h"
70 #include "main.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include <fcntl.h>
74 #include "gromacs/fileio/filenm.h"
75 #include "gromacs/fileio/futil.h"
76 #include "gromacs/fileio/gmxfio.h"
77 #include "gromacs/fileio/xdrf.h"
78 #include "gromacs/fileio/xdr_datatype.h"
79 #include "gromacs/utility/baseversion.h"
80 #include "gmx_fatal.h"
82 #include "buildinfo.h"
84 #ifdef GMX_FAHCORE
85 #include "corewrap.h"
86 #endif
88 #define CPT_MAGIC1 171817
89 #define CPT_MAGIC2 171819
90 #define CPTSTRLEN 1024
92 #ifdef GMX_DOUBLE
93 #define GMX_CPT_BUILD_DP 1
94 #else
95 #define GMX_CPT_BUILD_DP 0
96 #endif
98 /* cpt_version should normally only be changed
99 * when the header of footer format changes.
100 * The state data format itself is backward and forward compatible.
101 * But old code can not read a new entry that is present in the file
102 * (but can read a new format when new entries are not present).
104 static const int cpt_version = 16;
107 const char *est_names[estNR] =
109 "FE-lambda",
110 "box", "box-rel", "box-v", "pres_prev",
111 "nosehoover-xi", "thermostat-integral",
112 "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
113 "disre_initf", "disre_rm3tav",
114 "orire_initf", "orire_Dtav",
115 "svir_prev", "nosehoover-vxi", "v_eta", "vol0", "nhpres_xi", "nhpres_vxi", "fvir_prev", "fep_state", "MC-rng", "MC-rng-i"
118 enum {
119 eeksEKIN_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksEKINF, eeksEKINO, eeksEKINSCALEF, eeksEKINSCALEH, eeksVSCALE, eeksEKINTOTAL, eeksNR
122 const char *eeks_names[eeksNR] =
124 "Ekin_n", "Ekinh", "dEkindlambda", "mv_cos",
125 "Ekinf", "Ekinh_old", "EkinScaleF_NHC", "EkinScaleH_NHC", "Vscale_NHC", "Ekin_Total"
128 enum {
129 eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
130 eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
131 eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM,
132 eenhENERGY_DELTA_H_NN,
133 eenhENERGY_DELTA_H_LIST,
134 eenhENERGY_DELTA_H_STARTTIME,
135 eenhENERGY_DELTA_H_STARTLAMBDA,
136 eenhNR
139 const char *eenh_names[eenhNR] =
141 "energy_n", "energy_aver", "energy_sum", "energy_nsum",
142 "energy_sum_sim", "energy_nsum_sim",
143 "energy_nsteps", "energy_nsteps_sim",
144 "energy_delta_h_nn",
145 "energy_delta_h_list",
146 "energy_delta_h_start_time",
147 "energy_delta_h_start_lambda"
150 /* free energy history variables -- need to be preserved over checkpoint */
151 enum {
152 edfhBEQUIL, edfhNATLAMBDA, edfhWLHISTO, edfhWLDELTA, edfhSUMWEIGHTS, edfhSUMDG, edfhSUMMINVAR, edfhSUMVAR,
153 edfhACCUMP, edfhACCUMM, edfhACCUMP2, edfhACCUMM2, edfhTIJ, edfhTIJEMP, edfhNR
155 /* free energy history variable names */
156 const char *edfh_names[edfhNR] =
158 "bEquilibrated", "N_at_state", "Wang-Landau Histogram", "Wang-Landau Delta", "Weights", "Free Energies", "minvar", "variance",
159 "accumulated_plus", "accumulated_minus", "accumulated_plus_2", "accumulated_minus_2", "Tij", "Tij_empirical"
162 #ifdef GMX_NATIVE_WINDOWS
163 static int
164 gmx_wintruncate(const char *filename, __int64 size)
166 #ifdef GMX_FAHCORE
167 /*we do this elsewhere*/
168 return 0;
169 #else
170 FILE *fp;
171 int rc;
173 fp = fopen(filename, "rb+");
175 if (fp == NULL)
177 return -1;
180 return _chsize_s( fileno(fp), size);
181 #endif
183 #endif
186 enum {
187 ecprREAL, ecprRVEC, ecprMATRIX
190 enum {
191 cptpEST, cptpEEKS, cptpEENH, cptpEDFH
193 /* enums for the different components of checkpoint variables, replacing the hard coded ones.
194 cptpEST - state variables.
195 cptpEEKS - Kinetic energy state variables.
196 cptpEENH - Energy history state variables.
197 cptpEDFH - free energy history variables.
201 static const char *st_names(int cptp, int ecpt)
203 switch (cptp)
205 case cptpEST: return est_names [ecpt]; break;
206 case cptpEEKS: return eeks_names[ecpt]; break;
207 case cptpEENH: return eenh_names[ecpt]; break;
208 case cptpEDFH: return edfh_names[ecpt]; break;
211 return NULL;
214 static void cp_warning(FILE *fp)
216 fprintf(fp, "\nWARNING: Checkpoint file is corrupted or truncated\n\n");
219 static void cp_error()
221 gmx_fatal(FARGS, "Checkpoint file corrupted/truncated, or maybe you are out of disk space?");
224 static void do_cpt_string_err(XDR *xd, gmx_bool bRead, const char *desc, char **s, FILE *list)
226 bool_t res = 0;
228 if (bRead)
230 snew(*s, CPTSTRLEN);
232 res = xdr_string(xd, s, CPTSTRLEN);
233 if (res == 0)
235 cp_error();
237 if (list)
239 fprintf(list, "%s = %s\n", desc, *s);
240 sfree(*s);
244 static int do_cpt_int(XDR *xd, const char *desc, int *i, FILE *list)
246 bool_t res = 0;
248 res = xdr_int(xd, i);
249 if (res == 0)
251 return -1;
253 if (list)
255 fprintf(list, "%s = %d\n", desc, *i);
257 return 0;
260 static int do_cpt_u_chars(XDR *xd, const char *desc, int n, unsigned char *i, FILE *list)
262 bool_t res = 1;
263 int j;
264 if (list)
266 fprintf(list, "%s = ", desc);
268 for (j = 0; j < n && res; j++)
270 res &= xdr_u_char(xd, &i[j]);
271 if (list)
273 fprintf(list, "%02x", i[j]);
276 if (list)
278 fprintf(list, "\n");
280 if (res == 0)
282 return -1;
285 return 0;
288 static void do_cpt_int_err(XDR *xd, const char *desc, int *i, FILE *list)
290 if (do_cpt_int(xd, desc, i, list) < 0)
292 cp_error();
296 static void do_cpt_step_err(XDR *xd, const char *desc, gmx_int64_t *i, FILE *list)
298 bool_t res = 0;
299 char buf[STEPSTRSIZE];
301 res = xdr_int64(xd, i);
302 if (res == 0)
304 cp_error();
306 if (list)
308 fprintf(list, "%s = %s\n", desc, gmx_step_str(*i, buf));
312 static void do_cpt_double_err(XDR *xd, const char *desc, double *f, FILE *list)
314 bool_t res = 0;
316 res = xdr_double(xd, f);
317 if (res == 0)
319 cp_error();
321 if (list)
323 fprintf(list, "%s = %f\n", desc, *f);
327 static void do_cpt_real_err(XDR *xd, real *f)
329 bool_t res = 0;
331 #ifdef GMX_DOUBLE
332 res = xdr_double(xd, f);
333 #else
334 res = xdr_float(xd, f);
335 #endif
336 if (res == 0)
338 cp_error();
342 static void do_cpt_n_rvecs_err(XDR *xd, const char *desc, int n, rvec f[], FILE *list)
344 int i, j;
346 for (i = 0; i < n; i++)
348 for (j = 0; j < DIM; j++)
350 do_cpt_real_err(xd, &f[i][j]);
354 if (list)
356 pr_rvecs(list, 0, desc, f, n);
360 /* If nval >= 0, nval is used; on read this should match the passed value.
361 * If nval n<0, *nptr is used; on read the value is stored in nptr
363 static int do_cpte_reals_low(XDR *xd, int cptp, int ecpt, int sflags,
364 int nval, int *nptr, real **v,
365 FILE *list, int erealtype)
367 bool_t res = 0;
368 #ifndef GMX_DOUBLE
369 int dtc = xdr_datatype_float;
370 #else
371 int dtc = xdr_datatype_double;
372 #endif
373 real *vp, *va = NULL;
374 float *vf;
375 double *vd;
376 int nf, dt, i;
378 if (list == NULL)
380 if (nval >= 0)
382 nf = nval;
384 else
386 if (nptr == NULL)
388 gmx_incons("*ntpr=NULL in do_cpte_reals_low");
390 nf = *nptr;
393 res = xdr_int(xd, &nf);
394 if (res == 0)
396 return -1;
398 if (list == NULL)
400 if (nval >= 0)
402 if (nf != nval)
404 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), nval, nf);
407 else
409 *nptr = nf;
412 dt = dtc;
413 res = xdr_int(xd, &dt);
414 if (res == 0)
416 return -1;
418 if (dt != dtc)
420 fprintf(stderr, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
421 st_names(cptp, ecpt), xdr_datatype_names[dtc],
422 xdr_datatype_names[dt]);
424 if (list || !(sflags & (1<<ecpt)))
426 snew(va, nf);
427 vp = va;
429 else
431 if (*v == NULL)
433 snew(*v, nf);
435 vp = *v;
437 if (dt == xdr_datatype_float)
439 if (dtc == xdr_datatype_float)
441 vf = (float *)vp;
443 else
445 snew(vf, nf);
447 res = xdr_vector(xd, (char *)vf, nf,
448 (unsigned int)sizeof(float), (xdrproc_t)xdr_float);
449 if (res == 0)
451 return -1;
453 if (dtc != xdr_datatype_float)
455 for (i = 0; i < nf; i++)
457 vp[i] = vf[i];
459 sfree(vf);
462 else
464 if (dtc == xdr_datatype_double)
466 vd = (double *)vp;
468 else
470 snew(vd, nf);
472 res = xdr_vector(xd, (char *)vd, nf,
473 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
474 if (res == 0)
476 return -1;
478 if (dtc != xdr_datatype_double)
480 for (i = 0; i < nf; i++)
482 vp[i] = vd[i];
484 sfree(vd);
488 if (list)
490 switch (erealtype)
492 case ecprREAL:
493 pr_reals(list, 0, st_names(cptp, ecpt), vp, nf);
494 break;
495 case ecprRVEC:
496 pr_rvecs(list, 0, st_names(cptp, ecpt), (rvec *)vp, nf/3);
497 break;
498 default:
499 gmx_incons("Unknown checkpoint real type");
502 if (va)
504 sfree(va);
507 return 0;
511 /* This function stores n along with the reals for reading,
512 * but on reading it assumes that n matches the value in the checkpoint file,
513 * a fatal error is generated when this is not the case.
515 static int do_cpte_reals(XDR *xd, int cptp, int ecpt, int sflags,
516 int n, real **v, FILE *list)
518 return do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, v, list, ecprREAL);
521 /* This function does the same as do_cpte_reals,
522 * except that on reading it ignores the passed value of *n
523 * and stored the value read from the checkpoint file in *n.
525 static int do_cpte_n_reals(XDR *xd, int cptp, int ecpt, int sflags,
526 int *n, real **v, FILE *list)
528 return do_cpte_reals_low(xd, cptp, ecpt, sflags, -1, n, v, list, ecprREAL);
531 static int do_cpte_real(XDR *xd, int cptp, int ecpt, int sflags,
532 real *r, FILE *list)
534 int n;
536 return do_cpte_reals_low(xd, cptp, ecpt, sflags, 1, NULL, &r, list, ecprREAL);
539 static int do_cpte_ints(XDR *xd, int cptp, int ecpt, int sflags,
540 int n, int **v, FILE *list)
542 bool_t res = 0;
543 int dtc = xdr_datatype_int;
544 int *vp, *va = NULL;
545 int nf, dt, i;
547 nf = n;
548 res = xdr_int(xd, &nf);
549 if (res == 0)
551 return -1;
553 if (list == NULL && v != NULL && nf != n)
555 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
557 dt = dtc;
558 res = xdr_int(xd, &dt);
559 if (res == 0)
561 return -1;
563 if (dt != dtc)
565 gmx_fatal(FARGS, "Type mismatch for state entry %s, code type is %s, file type is %s\n",
566 st_names(cptp, ecpt), xdr_datatype_names[dtc],
567 xdr_datatype_names[dt]);
569 if (list || !(sflags & (1<<ecpt)) || v == NULL)
571 snew(va, nf);
572 vp = va;
574 else
576 if (*v == NULL)
578 snew(*v, nf);
580 vp = *v;
582 res = xdr_vector(xd, (char *)vp, nf,
583 (unsigned int)sizeof(int), (xdrproc_t)xdr_int);
584 if (res == 0)
586 return -1;
588 if (list)
590 pr_ivec(list, 0, st_names(cptp, ecpt), vp, nf, TRUE);
592 if (va)
594 sfree(va);
597 return 0;
600 static int do_cpte_int(XDR *xd, int cptp, int ecpt, int sflags,
601 int *i, FILE *list)
603 return do_cpte_ints(xd, cptp, ecpt, sflags, 1, &i, list);
606 static int do_cpte_doubles(XDR *xd, int cptp, int ecpt, int sflags,
607 int n, double **v, FILE *list)
609 bool_t res = 0;
610 int dtc = xdr_datatype_double;
611 double *vp, *va = NULL;
612 int nf, dt, i;
614 nf = n;
615 res = xdr_int(xd, &nf);
616 if (res == 0)
618 return -1;
620 if (list == NULL && nf != n)
622 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
624 dt = dtc;
625 res = xdr_int(xd, &dt);
626 if (res == 0)
628 return -1;
630 if (dt != dtc)
632 gmx_fatal(FARGS, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
633 st_names(cptp, ecpt), xdr_datatype_names[dtc],
634 xdr_datatype_names[dt]);
636 if (list || !(sflags & (1<<ecpt)))
638 snew(va, nf);
639 vp = va;
641 else
643 if (*v == NULL)
645 snew(*v, nf);
647 vp = *v;
649 res = xdr_vector(xd, (char *)vp, nf,
650 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
651 if (res == 0)
653 return -1;
655 if (list)
657 pr_doubles(list, 0, st_names(cptp, ecpt), vp, nf);
659 if (va)
661 sfree(va);
664 return 0;
667 static int do_cpte_double(XDR *xd, int cptp, int ecpt, int sflags,
668 double *r, FILE *list)
670 return do_cpte_doubles(xd, cptp, ecpt, sflags, 1, &r, list);
674 static int do_cpte_rvecs(XDR *xd, int cptp, int ecpt, int sflags,
675 int n, rvec **v, FILE *list)
677 int n3;
679 return do_cpte_reals_low(xd, cptp, ecpt, sflags,
680 n*DIM, NULL, (real **)v, list, ecprRVEC);
683 static int do_cpte_matrix(XDR *xd, int cptp, int ecpt, int sflags,
684 matrix v, FILE *list)
686 real *vr;
687 real ret;
689 vr = (real *)&(v[0][0]);
690 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
691 DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
693 if (list && ret == 0)
695 pr_rvecs(list, 0, st_names(cptp, ecpt), v, DIM);
698 return ret;
702 static int do_cpte_nmatrix(XDR *xd, int cptp, int ecpt, int sflags,
703 int n, real **v, FILE *list)
705 int i;
706 real *vr;
707 real ret, reti;
708 char name[CPTSTRLEN];
710 ret = 0;
711 if (v == NULL)
713 snew(v, n);
715 for (i = 0; i < n; i++)
717 reti = 0;
718 vr = v[i];
719 reti = do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, &(v[i]), NULL, ecprREAL);
720 if (list && reti == 0)
722 sprintf(name, "%s[%d]", st_names(cptp, ecpt), i);
723 pr_reals(list, 0, name, v[i], n);
725 if (reti == 0)
727 ret = 0;
730 return ret;
733 static int do_cpte_matrices(XDR *xd, int cptp, int ecpt, int sflags,
734 int n, matrix **v, FILE *list)
736 bool_t res = 0;
737 matrix *vp, *va = NULL;
738 real *vr;
739 int nf, i, j, k;
740 int ret;
742 nf = n;
743 res = xdr_int(xd, &nf);
744 if (res == 0)
746 return -1;
748 if (list == NULL && nf != n)
750 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
752 if (list || !(sflags & (1<<ecpt)))
754 snew(va, nf);
755 vp = va;
757 else
759 if (*v == NULL)
761 snew(*v, nf);
763 vp = *v;
765 snew(vr, nf*DIM*DIM);
766 for (i = 0; i < nf; i++)
768 for (j = 0; j < DIM; j++)
770 for (k = 0; k < DIM; k++)
772 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
776 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
777 nf*DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
778 for (i = 0; i < nf; i++)
780 for (j = 0; j < DIM; j++)
782 for (k = 0; k < DIM; k++)
784 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
788 sfree(vr);
790 if (list && ret == 0)
792 for (i = 0; i < nf; i++)
794 pr_rvecs(list, 0, st_names(cptp, ecpt), vp[i], DIM);
797 if (va)
799 sfree(va);
802 return ret;
805 static void do_cpt_header(XDR *xd, gmx_bool bRead, int *file_version,
806 char **version, char **btime, char **buser, char **bhost,
807 int *double_prec,
808 char **fprog, char **ftime,
809 int *eIntegrator, int *simulation_part,
810 gmx_int64_t *step, double *t,
811 int *nnodes, int *dd_nc, int *npme,
812 int *natoms, int *ngtc, int *nnhpres, int *nhchainlength,
813 int *nlambda, int *flags_state,
814 int *flags_eks, int *flags_enh, int *flags_dfh,
815 int *nED, int *eSwapCoords,
816 FILE *list)
818 bool_t res = 0;
819 int magic;
820 int idum = 0;
821 int i;
822 char *fhost;
824 if (bRead)
826 magic = -1;
828 else
830 magic = CPT_MAGIC1;
832 res = xdr_int(xd, &magic);
833 if (res == 0)
835 gmx_fatal(FARGS, "The checkpoint file is empty/corrupted, or maybe you are out of disk space?");
837 if (magic != CPT_MAGIC1)
839 gmx_fatal(FARGS, "Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
840 "The checkpoint file is corrupted or not a checkpoint file",
841 magic, CPT_MAGIC1);
843 if (!bRead)
845 snew(fhost, 255);
846 gmx_gethostname(fhost, 255);
848 do_cpt_string_err(xd, bRead, "GROMACS version", version, list);
849 do_cpt_string_err(xd, bRead, "GROMACS build time", btime, list);
850 do_cpt_string_err(xd, bRead, "GROMACS build user", buser, list);
851 do_cpt_string_err(xd, bRead, "GROMACS build host", bhost, list);
852 do_cpt_string_err(xd, bRead, "generating program", fprog, list);
853 do_cpt_string_err(xd, bRead, "generation time", ftime, list);
854 *file_version = cpt_version;
855 do_cpt_int_err(xd, "checkpoint file version", file_version, list);
856 if (*file_version > cpt_version)
858 gmx_fatal(FARGS, "Attempting to read a checkpoint file of version %d with code of version %d\n", *file_version, cpt_version);
860 if (*file_version >= 13)
862 do_cpt_int_err(xd, "GROMACS double precision", double_prec, list);
864 else
866 *double_prec = -1;
868 if (*file_version >= 12)
870 do_cpt_string_err(xd, bRead, "generating host", &fhost, list);
871 if (list == NULL)
873 sfree(fhost);
876 do_cpt_int_err(xd, "#atoms", natoms, list);
877 do_cpt_int_err(xd, "#T-coupling groups", ngtc, list);
878 if (*file_version >= 10)
880 do_cpt_int_err(xd, "#Nose-Hoover T-chains", nhchainlength, list);
882 else
884 *nhchainlength = 1;
886 if (*file_version >= 11)
888 do_cpt_int_err(xd, "#Nose-Hoover T-chains for barostat ", nnhpres, list);
890 else
892 *nnhpres = 0;
894 if (*file_version >= 14)
896 do_cpt_int_err(xd, "# of total lambda states ", nlambda, list);
898 else
900 *nlambda = 0;
902 do_cpt_int_err(xd, "integrator", eIntegrator, list);
903 if (*file_version >= 3)
905 do_cpt_int_err(xd, "simulation part #", simulation_part, list);
907 else
909 *simulation_part = 1;
911 if (*file_version >= 5)
913 do_cpt_step_err(xd, "step", step, list);
915 else
917 do_cpt_int_err(xd, "step", &idum, list);
918 *step = idum;
920 do_cpt_double_err(xd, "t", t, list);
921 do_cpt_int_err(xd, "#PP-nodes", nnodes, list);
922 idum = 1;
923 do_cpt_int_err(xd, "dd_nc[x]", dd_nc ? &(dd_nc[0]) : &idum, list);
924 do_cpt_int_err(xd, "dd_nc[y]", dd_nc ? &(dd_nc[1]) : &idum, list);
925 do_cpt_int_err(xd, "dd_nc[z]", dd_nc ? &(dd_nc[2]) : &idum, list);
926 do_cpt_int_err(xd, "#PME-only nodes", npme, list);
927 do_cpt_int_err(xd, "state flags", flags_state, list);
928 if (*file_version >= 4)
930 do_cpt_int_err(xd, "ekin data flags", flags_eks, list);
931 do_cpt_int_err(xd, "energy history flags", flags_enh, list);
933 else
935 *flags_eks = 0;
936 *flags_enh = (*flags_state >> (estORIRE_DTAV+1));
937 *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
938 (1<<(estORIRE_DTAV+2)) |
939 (1<<(estORIRE_DTAV+3))));
941 if (*file_version >= 14)
943 do_cpt_int_err(xd, "df history flags", flags_dfh, list);
945 else
947 *flags_dfh = 0;
950 if (*file_version >= 15)
952 do_cpt_int_err(xd, "ED data sets", nED, list);
954 else
956 *nED = 0;
958 if (*file_version >= 16)
960 do_cpt_int_err(xd, "swap", eSwapCoords, list);
964 static int do_cpt_footer(XDR *xd, int file_version)
966 bool_t res = 0;
967 int magic;
969 if (file_version >= 2)
971 magic = CPT_MAGIC2;
972 res = xdr_int(xd, &magic);
973 if (res == 0)
975 cp_error();
977 if (magic != CPT_MAGIC2)
979 return -1;
983 return 0;
986 static int do_cpt_state(XDR *xd, gmx_bool bRead,
987 int fflags, t_state *state,
988 FILE *list)
990 int sflags;
991 int i;
992 int ret;
993 int nnht, nnhtp;
995 ret = 0;
997 nnht = state->nhchainlength*state->ngtc;
998 nnhtp = state->nhchainlength*state->nnhpres;
1000 if (bRead) /* we need to allocate space for dfhist if we are reading */
1002 init_df_history(&state->dfhist, state->dfhist.nlambda);
1005 sflags = state->flags;
1006 for (i = 0; (i < estNR && ret == 0); i++)
1008 if (fflags & (1<<i))
1010 switch (i)
1012 case estLAMBDA: ret = do_cpte_reals(xd, cptpEST, i, sflags, efptNR, &(state->lambda), list); break;
1013 case estFEPSTATE: ret = do_cpte_int (xd, cptpEST, i, sflags, &state->fep_state, list); break;
1014 case estBOX: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box, list); break;
1015 case estBOX_REL: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box_rel, list); break;
1016 case estBOXV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->boxv, list); break;
1017 case estPRES_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->pres_prev, list); break;
1018 case estSVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->svir_prev, list); break;
1019 case estFVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->fvir_prev, list); break;
1020 case estNH_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_xi, list); break;
1021 case estNH_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_vxi, list); break;
1022 case estNHPRES_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_xi, list); break;
1023 case estNHPRES_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_vxi, list); break;
1024 case estTC_INT: ret = do_cpte_doubles(xd, cptpEST, i, sflags, state->ngtc, &state->therm_integral, list); break;
1025 case estVETA: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->veta, list); break;
1026 case estVOL0: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->vol0, list); break;
1027 case estX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->x, list); break;
1028 case estV: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->v, list); break;
1029 case estSDX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->sd_X, list); break;
1030 /* The RNG entries are no longer written,
1031 * the next 4 lines are only for reading old files.
1033 case estLD_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1034 case estLD_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1035 case estMC_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1036 case estMC_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1037 case estDISRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.disre_initf, list); break;
1038 case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.ndisrepairs, &state->hist.disre_rm3tav, list); break;
1039 case estORIRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.orire_initf, list); break;
1040 case estORIRE_DTAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.norire_Dtav, &state->hist.orire_Dtav, list); break;
1041 default:
1042 gmx_fatal(FARGS, "Unknown state entry %d\n"
1043 "You are probably reading a new checkpoint file with old code", i);
1048 return ret;
1051 static int do_cpt_ekinstate(XDR *xd, int fflags, ekinstate_t *ekins,
1052 FILE *list)
1054 int i;
1055 int ret;
1057 ret = 0;
1059 for (i = 0; (i < eeksNR && ret == 0); i++)
1061 if (fflags & (1<<i))
1063 switch (i)
1066 case eeksEKIN_N: ret = do_cpte_int(xd, cptpEEKS, i, fflags, &ekins->ekin_n, list); break;
1067 case eeksEKINH: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh, list); break;
1068 case eeksEKINF: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinf, list); break;
1069 case eeksEKINO: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh_old, list); break;
1070 case eeksEKINTOTAL: ret = do_cpte_matrix(xd, cptpEEKS, i, fflags, ekins->ekin_total, list); break;
1071 case eeksEKINSCALEF: ret = do_cpte_doubles(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinscalef_nhc, list); break;
1072 case eeksVSCALE: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->vscale_nhc, list); break;
1073 case eeksEKINSCALEH: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->ekinscaleh_nhc, list); break;
1074 case eeksDEKINDL: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->dekindl, list); break;
1075 case eeksMVCOS: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->mvcos, list); break;
1076 default:
1077 gmx_fatal(FARGS, "Unknown ekin data state entry %d\n"
1078 "You are probably reading a new checkpoint file with old code", i);
1083 return ret;
1087 static int do_cpt_swapstate(XDR *xd, gmx_bool bRead, swapstate_t *swapstate, FILE *list)
1089 int ii, ic, j;
1090 int ret = 0;
1091 int swap_cpt_version = 1;
1094 if (eswapNO == swapstate->eSwapCoords)
1096 return ret;
1099 swapstate->bFromCpt = bRead;
1101 do_cpt_int_err(xd, "swap checkpoint version", &swap_cpt_version, list);
1102 do_cpt_int_err(xd, "swap coupling steps", &swapstate->nAverage, list);
1104 /* When reading, init_swapcoords has not been called yet,
1105 * so we have to allocate memory first. */
1107 for (ic = 0; ic < eCompNR; ic++)
1109 for (ii = 0; ii < eIonNR; ii++)
1111 if (bRead)
1113 do_cpt_int_err(xd, "swap requested atoms", &swapstate->nat_req[ic][ii], list);
1115 else
1117 do_cpt_int_err(xd, "swap requested atoms p", swapstate->nat_req_p[ic][ii], list);
1120 if (bRead)
1122 do_cpt_int_err(xd, "swap influx netto", &swapstate->inflow_netto[ic][ii], list);
1124 else
1126 do_cpt_int_err(xd, "swap influx netto p", swapstate->inflow_netto_p[ic][ii], list);
1129 if (bRead && (NULL == swapstate->nat_past[ic][ii]) )
1131 snew(swapstate->nat_past[ic][ii], swapstate->nAverage);
1134 for (j = 0; j < swapstate->nAverage; j++)
1136 if (bRead)
1138 do_cpt_int_err(xd, "swap past atom counts", &swapstate->nat_past[ic][ii][j], list);
1140 else
1142 do_cpt_int_err(xd, "swap past atom counts p", &swapstate->nat_past_p[ic][ii][j], list);
1148 /* Ion flux per channel */
1149 for (ic = 0; ic < eChanNR; ic++)
1151 for (ii = 0; ii < eIonNR; ii++)
1153 if (bRead)
1155 do_cpt_int_err(xd, "channel flux", &swapstate->fluxfromAtoB[ic][ii], list);
1157 else
1159 do_cpt_int_err(xd, "channel flux p", swapstate->fluxfromAtoB_p[ic][ii], list);
1164 /* Ion flux leakage */
1165 if (bRead)
1167 snew(swapstate->fluxleak, 1);
1169 do_cpt_int_err(xd, "flux leakage", swapstate->fluxleak, list);
1171 /* Ion history */
1172 do_cpt_int_err(xd, "number of ions", &swapstate->nions, list);
1174 if (bRead)
1176 snew(swapstate->channel_label, swapstate->nions);
1177 snew(swapstate->comp_from, swapstate->nions);
1180 do_cpt_u_chars(xd, "channel history", swapstate->nions, swapstate->channel_label, list);
1181 do_cpt_u_chars(xd, "domain history", swapstate->nions, swapstate->comp_from, list);
1183 /* Save the last known whole positions to checkpoint
1184 * file to be able to also make multimeric channels whole in PBC */
1185 do_cpt_int_err(xd, "Ch0 atoms", &swapstate->nat[eChan0], list);
1186 do_cpt_int_err(xd, "Ch1 atoms", &swapstate->nat[eChan1], list);
1187 if (bRead)
1189 snew(swapstate->xc_old_whole[eChan0], swapstate->nat[eChan0]);
1190 snew(swapstate->xc_old_whole[eChan1], swapstate->nat[eChan1]);
1191 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], swapstate->xc_old_whole[eChan0], list);
1192 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], swapstate->xc_old_whole[eChan1], list);
1194 else
1196 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], *swapstate->xc_old_whole_p[eChan0], list);
1197 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], *swapstate->xc_old_whole_p[eChan1], list);
1200 return ret;
1204 static int do_cpt_enerhist(XDR *xd, gmx_bool bRead,
1205 int fflags, energyhistory_t *enerhist,
1206 FILE *list)
1208 int i;
1209 int j;
1210 int ret;
1212 ret = 0;
1214 if (bRead)
1216 enerhist->nsteps = 0;
1217 enerhist->nsum = 0;
1218 enerhist->nsteps_sim = 0;
1219 enerhist->nsum_sim = 0;
1220 enerhist->dht = NULL;
1222 if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
1224 snew(enerhist->dht, 1);
1225 enerhist->dht->ndh = NULL;
1226 enerhist->dht->dh = NULL;
1227 enerhist->dht->start_lambda_set = FALSE;
1231 for (i = 0; (i < eenhNR && ret == 0); i++)
1233 if (fflags & (1<<i))
1235 switch (i)
1237 case eenhENERGY_N: ret = do_cpte_int(xd, cptpEENH, i, fflags, &enerhist->nener, list); break;
1238 case eenhENERGY_AVER: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_ave, list); break;
1239 case eenhENERGY_SUM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum, list); break;
1240 case eenhENERGY_NSUM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum, list); break;
1241 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum_sim, list); break;
1242 case eenhENERGY_NSUM_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum_sim, list); break;
1243 case eenhENERGY_NSTEPS: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps, list); break;
1244 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps_sim, list); break;
1245 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd, eenh_names[i], &(enerhist->dht->nndh), list);
1246 if (bRead) /* now allocate memory for it */
1248 snew(enerhist->dht->dh, enerhist->dht->nndh);
1249 snew(enerhist->dht->ndh, enerhist->dht->nndh);
1250 for (j = 0; j < enerhist->dht->nndh; j++)
1252 enerhist->dht->ndh[j] = 0;
1253 enerhist->dht->dh[j] = NULL;
1256 break;
1257 case eenhENERGY_DELTA_H_LIST:
1258 for (j = 0; j < enerhist->dht->nndh; j++)
1260 ret = do_cpte_n_reals(xd, cptpEENH, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list);
1262 break;
1263 case eenhENERGY_DELTA_H_STARTTIME:
1264 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_time), list); break;
1265 case eenhENERGY_DELTA_H_STARTLAMBDA:
1266 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_lambda), list); break;
1267 default:
1268 gmx_fatal(FARGS, "Unknown energy history entry %d\n"
1269 "You are probably reading a new checkpoint file with old code", i);
1274 if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
1276 /* Assume we have an old file format and copy sum to sum_sim */
1277 srenew(enerhist->ener_sum_sim, enerhist->nener);
1278 for (i = 0; i < enerhist->nener; i++)
1280 enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1282 fflags |= (1<<eenhENERGY_SUM_SIM);
1285 if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1286 !(fflags & (1<<eenhENERGY_NSTEPS)))
1288 /* Assume we have an old file format and copy nsum to nsteps */
1289 enerhist->nsteps = enerhist->nsum;
1290 fflags |= (1<<eenhENERGY_NSTEPS);
1292 if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1293 !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1295 /* Assume we have an old file format and copy nsum to nsteps */
1296 enerhist->nsteps_sim = enerhist->nsum_sim;
1297 fflags |= (1<<eenhENERGY_NSTEPS_SIM);
1300 return ret;
1303 static int do_cpt_df_hist(XDR *xd, int fflags, df_history_t *dfhist, FILE *list)
1305 int i, nlambda;
1306 int ret;
1308 nlambda = dfhist->nlambda;
1309 ret = 0;
1311 for (i = 0; (i < edfhNR && ret == 0); i++)
1313 if (fflags & (1<<i))
1315 switch (i)
1317 case edfhBEQUIL: ret = do_cpte_int(xd, cptpEDFH, i, fflags, &dfhist->bEquil, list); break;
1318 case edfhNATLAMBDA: ret = do_cpte_ints(xd, cptpEDFH, i, fflags, nlambda, &dfhist->n_at_lam, list); break;
1319 case edfhWLHISTO: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->wl_histo, list); break;
1320 case edfhWLDELTA: ret = do_cpte_real(xd, cptpEDFH, i, fflags, &dfhist->wl_delta, list); break;
1321 case edfhSUMWEIGHTS: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_weights, list); break;
1322 case edfhSUMDG: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_dg, list); break;
1323 case edfhSUMMINVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_minvar, list); break;
1324 case edfhSUMVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_variance, list); break;
1325 case edfhACCUMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p, list); break;
1326 case edfhACCUMM: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m, list); break;
1327 case edfhACCUMP2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p2, list); break;
1328 case edfhACCUMM2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m2, list); break;
1329 case edfhTIJ: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij, list); break;
1330 case edfhTIJEMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij_empirical, list); break;
1332 default:
1333 gmx_fatal(FARGS, "Unknown df history entry %d\n"
1334 "You are probably reading a new checkpoint file with old code", i);
1339 return ret;
1343 /* This function stores the last whole configuration of the reference and
1344 * average structure in the .cpt file
1346 static int do_cpt_EDstate(XDR *xd, gmx_bool bRead,
1347 edsamstate_t *EDstate, FILE *list)
1349 int i, j;
1350 int ret = 0;
1351 char buf[STRLEN];
1354 EDstate->bFromCpt = bRead;
1356 if (EDstate->nED <= 0)
1358 return ret;
1361 /* When reading, init_edsam has not been called yet,
1362 * so we have to allocate memory first. */
1363 if (bRead)
1365 snew(EDstate->nref, EDstate->nED);
1366 snew(EDstate->old_sref, EDstate->nED);
1367 snew(EDstate->nav, EDstate->nED);
1368 snew(EDstate->old_sav, EDstate->nED);
1371 /* Read/write the last whole conformation of SREF and SAV for each ED dataset (usually only one) */
1372 for (i = 0; i < EDstate->nED; i++)
1374 /* Reference structure SREF */
1375 sprintf(buf, "ED%d # of atoms in reference structure", i+1);
1376 do_cpt_int_err(xd, buf, &EDstate->nref[i], list);
1377 sprintf(buf, "ED%d x_ref", i+1);
1378 if (bRead)
1380 snew(EDstate->old_sref[i], EDstate->nref[i]);
1381 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref[i], list);
1383 else
1385 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref_p[i], list);
1388 /* Average structure SAV */
1389 sprintf(buf, "ED%d # of atoms in average structure", i+1);
1390 do_cpt_int_err(xd, buf, &EDstate->nav[i], list);
1391 sprintf(buf, "ED%d x_av", i+1);
1392 if (bRead)
1394 snew(EDstate->old_sav[i], EDstate->nav[i]);
1395 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav[i], list);
1397 else
1399 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav_p[i], list);
1403 return ret;
1407 static int do_cpt_files(XDR *xd, gmx_bool bRead,
1408 gmx_file_position_t **p_outputfiles, int *nfiles,
1409 FILE *list, int file_version)
1411 int i, j;
1412 gmx_off_t offset;
1413 gmx_off_t mask = 0xFFFFFFFFL;
1414 int offset_high, offset_low;
1415 char *buf;
1416 gmx_file_position_t *outputfiles;
1418 if (do_cpt_int(xd, "number of output files", nfiles, list) != 0)
1420 return -1;
1423 if (bRead)
1425 snew(*p_outputfiles, *nfiles);
1428 outputfiles = *p_outputfiles;
1430 for (i = 0; i < *nfiles; i++)
1432 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1433 if (bRead)
1435 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1436 strncpy(outputfiles[i].filename, buf, CPTSTRLEN-1);
1437 if (list == NULL)
1439 sfree(buf);
1442 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1444 return -1;
1446 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1448 return -1;
1450 outputfiles[i].offset = ( ((gmx_off_t) offset_high) << 32 ) | ( (gmx_off_t) offset_low & mask );
1452 else
1454 buf = outputfiles[i].filename;
1455 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1456 /* writing */
1457 offset = outputfiles[i].offset;
1458 if (offset == -1)
1460 offset_low = -1;
1461 offset_high = -1;
1463 else
1465 offset_low = (int) (offset & mask);
1466 offset_high = (int) ((offset >> 32) & mask);
1468 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1470 return -1;
1472 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1474 return -1;
1477 if (file_version >= 8)
1479 if (do_cpt_int(xd, "file_checksum_size", &(outputfiles[i].chksum_size),
1480 list) != 0)
1482 return -1;
1484 if (do_cpt_u_chars(xd, "file_checksum", 16, outputfiles[i].chksum, list) != 0)
1486 return -1;
1489 else
1491 outputfiles[i].chksum_size = -1;
1494 return 0;
1498 void write_checkpoint(const char *fn, gmx_bool bNumberAndKeep,
1499 FILE *fplog, t_commrec *cr,
1500 int eIntegrator, int simulation_part,
1501 gmx_bool bExpanded, int elamstats,
1502 gmx_int64_t step, double t, t_state *state)
1504 t_fileio *fp;
1505 int file_version;
1506 char *version;
1507 char *btime;
1508 char *buser;
1509 char *bhost;
1510 int double_prec;
1511 char *fprog;
1512 char *fntemp; /* the temporary checkpoint file name */
1513 time_t now;
1514 char timebuf[STRLEN];
1515 int nppnodes, npmenodes, flag_64bit;
1516 char buf[1024], suffix[5+STEPSTRSIZE], sbuf[STEPSTRSIZE];
1517 gmx_file_position_t *outputfiles;
1518 int noutputfiles;
1519 char *ftime;
1520 int flags_eks, flags_enh, flags_dfh, i;
1521 t_fileio *ret;
1523 if (DOMAINDECOMP(cr))
1525 nppnodes = cr->dd->nnodes;
1526 npmenodes = cr->npmenodes;
1528 else
1530 nppnodes = 1;
1531 npmenodes = 0;
1534 #ifndef GMX_NO_RENAME
1535 /* make the new temporary filename */
1536 snew(fntemp, strlen(fn)+5+STEPSTRSIZE);
1537 strcpy(fntemp, fn);
1538 fntemp[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1539 sprintf(suffix, "_%s%s", "step", gmx_step_str(step, sbuf));
1540 strcat(fntemp, suffix);
1541 strcat(fntemp, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1542 #else
1543 /* if we can't rename, we just overwrite the cpt file.
1544 * dangerous if interrupted.
1546 snew(fntemp, strlen(fn));
1547 strcpy(fntemp, fn);
1548 #endif
1549 time(&now);
1550 gmx_ctime_r(&now, timebuf, STRLEN);
1552 if (fplog)
1554 fprintf(fplog, "Writing checkpoint, step %s at %s\n\n",
1555 gmx_step_str(step, buf), timebuf);
1558 /* Get offsets for open files */
1559 gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1561 fp = gmx_fio_open(fntemp, "w");
1563 if (state->ekinstate.bUpToDate)
1565 flags_eks =
1566 ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) |
1567 (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) |
1568 (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1570 else
1572 flags_eks = 0;
1575 flags_enh = 0;
1576 if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1578 flags_enh |= (1<<eenhENERGY_N);
1579 if (state->enerhist.nsum > 0)
1581 flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1582 (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSUM));
1584 if (state->enerhist.nsum_sim > 0)
1586 flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSTEPS_SIM) |
1587 (1<<eenhENERGY_NSUM_SIM));
1589 if (state->enerhist.dht)
1591 flags_enh |= ( (1<< eenhENERGY_DELTA_H_NN) |
1592 (1<< eenhENERGY_DELTA_H_LIST) |
1593 (1<< eenhENERGY_DELTA_H_STARTTIME) |
1594 (1<< eenhENERGY_DELTA_H_STARTLAMBDA) );
1598 if (bExpanded)
1600 flags_dfh = ((1<<edfhBEQUIL) | (1<<edfhNATLAMBDA) | (1<<edfhSUMWEIGHTS) | (1<<edfhSUMDG) |
1601 (1<<edfhTIJ) | (1<<edfhTIJEMP));
1602 if (EWL(elamstats))
1604 flags_dfh |= ((1<<edfhWLDELTA) | (1<<edfhWLHISTO));
1606 if ((elamstats == elamstatsMINVAR) || (elamstats == elamstatsBARKER) || (elamstats == elamstatsMETROPOLIS))
1608 flags_dfh |= ((1<<edfhACCUMP) | (1<<edfhACCUMM) | (1<<edfhACCUMP2) | (1<<edfhACCUMM2)
1609 | (1<<edfhSUMMINVAR) | (1<<edfhSUMVAR));
1612 else
1614 flags_dfh = 0;
1617 /* We can check many more things now (CPU, acceleration, etc), but
1618 * it is highly unlikely to have two separate builds with exactly
1619 * the same version, user, time, and build host!
1622 version = gmx_strdup(gmx_version());
1623 btime = gmx_strdup(BUILD_TIME);
1624 buser = gmx_strdup(BUILD_USER);
1625 bhost = gmx_strdup(BUILD_HOST);
1627 double_prec = GMX_CPT_BUILD_DP;
1628 fprog = gmx_strdup(Program());
1630 ftime = &(timebuf[0]);
1632 do_cpt_header(gmx_fio_getxdr(fp), FALSE, &file_version,
1633 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1634 &eIntegrator, &simulation_part, &step, &t, &nppnodes,
1635 DOMAINDECOMP(cr) ? cr->dd->nc : NULL, &npmenodes,
1636 &state->natoms, &state->ngtc, &state->nnhpres,
1637 &state->nhchainlength, &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
1638 &state->edsamstate.nED, &state->swapstate.eSwapCoords,
1639 NULL);
1641 sfree(version);
1642 sfree(btime);
1643 sfree(buser);
1644 sfree(bhost);
1645 sfree(fprog);
1647 if ((do_cpt_state(gmx_fio_getxdr(fp), FALSE, state->flags, state, NULL) < 0) ||
1648 (do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL) < 0) ||
1649 (do_cpt_enerhist(gmx_fio_getxdr(fp), FALSE, flags_enh, &state->enerhist, NULL) < 0) ||
1650 (do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL) < 0) ||
1651 (do_cpt_EDstate(gmx_fio_getxdr(fp), FALSE, &state->edsamstate, NULL) < 0) ||
1652 (do_cpt_swapstate(gmx_fio_getxdr(fp), FALSE, &state->swapstate, NULL) < 0) ||
1653 (do_cpt_files(gmx_fio_getxdr(fp), FALSE, &outputfiles, &noutputfiles, NULL,
1654 file_version) < 0))
1656 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1659 do_cpt_footer(gmx_fio_getxdr(fp), file_version);
1661 /* we really, REALLY, want to make sure to physically write the checkpoint,
1662 and all the files it depends on, out to disk. Because we've
1663 opened the checkpoint with gmx_fio_open(), it's in our list
1664 of open files. */
1665 ret = gmx_fio_all_output_fsync();
1667 if (ret)
1669 char buf[STRLEN];
1670 sprintf(buf,
1671 "Cannot fsync '%s'; maybe you are out of disk space?",
1672 gmx_fio_getname(ret));
1674 if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV) == NULL)
1676 gmx_file(buf);
1678 else
1680 gmx_warning(buf);
1684 if (gmx_fio_close(fp) != 0)
1686 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1689 /* we don't move the checkpoint if the user specified they didn't want it,
1690 or if the fsyncs failed */
1691 #ifndef GMX_NO_RENAME
1692 if (!bNumberAndKeep && !ret)
1694 if (gmx_fexist(fn))
1696 /* Rename the previous checkpoint file */
1697 strcpy(buf, fn);
1698 buf[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1699 strcat(buf, "_prev");
1700 strcat(buf, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1701 #ifndef GMX_FAHCORE
1702 /* we copy here so that if something goes wrong between now and
1703 * the rename below, there's always a state.cpt.
1704 * If renames are atomic (such as in POSIX systems),
1705 * this copying should be unneccesary.
1707 gmx_file_copy(fn, buf, FALSE);
1708 /* We don't really care if this fails:
1709 * there's already a new checkpoint.
1711 #else
1712 gmx_file_rename(fn, buf);
1713 #endif
1715 if (gmx_file_rename(fntemp, fn) != 0)
1717 gmx_file("Cannot rename checkpoint file; maybe you are out of disk space?");
1720 #endif /* GMX_NO_RENAME */
1722 sfree(outputfiles);
1723 sfree(fntemp);
1725 #ifdef GMX_FAHCORE
1726 /*code for alternate checkpointing scheme. moved from top of loop over
1727 steps */
1728 fcRequestCheckPoint();
1729 if (fcCheckPointParallel( cr->nodeid, NULL, 0) == 0)
1731 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", step );
1733 #endif /* end GMX_FAHCORE block */
1736 static void print_flag_mismatch(FILE *fplog, int sflags, int fflags)
1738 int i;
1740 fprintf(fplog, "\nState entry mismatch between the simulation and the checkpoint file\n");
1741 fprintf(fplog, "Entries which are not present in the checkpoint file will not be updated\n");
1742 fprintf(fplog, " %24s %11s %11s\n", "", "simulation", "checkpoint");
1743 for (i = 0; i < estNR; i++)
1745 if ((sflags & (1<<i)) || (fflags & (1<<i)))
1747 fprintf(fplog, " %24s %11s %11s\n",
1748 est_names[i],
1749 (sflags & (1<<i)) ? " present " : "not present",
1750 (fflags & (1<<i)) ? " present " : "not present");
1755 static void check_int(FILE *fplog, const char *type, int p, int f, gmx_bool *mm)
1757 FILE *fp = fplog ? fplog : stderr;
1759 if (p != f)
1761 fprintf(fp, " %s mismatch,\n", type);
1762 fprintf(fp, " current program: %d\n", p);
1763 fprintf(fp, " checkpoint file: %d\n", f);
1764 fprintf(fp, "\n");
1765 *mm = TRUE;
1769 static void check_string(FILE *fplog, const char *type, const char *p,
1770 const char *f, gmx_bool *mm)
1772 FILE *fp = fplog ? fplog : stderr;
1774 if (strcmp(p, f) != 0)
1776 fprintf(fp, " %s mismatch,\n", type);
1777 fprintf(fp, " current program: %s\n", p);
1778 fprintf(fp, " checkpoint file: %s\n", f);
1779 fprintf(fp, "\n");
1780 *mm = TRUE;
1784 static void check_match(FILE *fplog,
1785 char *version,
1786 char *btime, char *buser, char *bhost, int double_prec,
1787 char *fprog,
1788 t_commrec *cr, int npp_f, int npme_f,
1789 ivec dd_nc, ivec dd_nc_f)
1791 int npp;
1792 gmx_bool mm = FALSE;
1793 gmx_bool patchlevel_differs = FALSE;
1794 gmx_bool version_differs = FALSE;
1796 check_string(fplog, "Version", gmx_version(), version, &mm);
1797 patchlevel_differs = mm;
1799 if (patchlevel_differs)
1801 /* Gromacs should be able to continue from checkpoints between
1802 * different patch level versions, but we do not guarantee
1803 * compatibility between different major/minor versions - check this.
1805 int gmx_major, gmx_minor;
1806 int cpt_major, cpt_minor;
1807 sscanf(gmx_version(), "VERSION %d.%d", &gmx_major, &gmx_minor);
1808 sscanf(version, "VERSION %d.%d", &cpt_major, &cpt_minor);
1809 version_differs = (gmx_major != cpt_major || gmx_minor != cpt_minor);
1812 check_string(fplog, "Build time", BUILD_TIME, btime, &mm);
1813 check_string(fplog, "Build user", BUILD_USER, buser, &mm);
1814 check_string(fplog, "Build host", BUILD_HOST, bhost, &mm);
1815 check_int (fplog, "Double prec.", GMX_CPT_BUILD_DP, double_prec, &mm);
1816 check_string(fplog, "Program name", Program(), fprog, &mm);
1818 check_int (fplog, "#nodes", cr->nnodes, npp_f+npme_f, &mm);
1819 if (cr->nnodes > 1)
1821 check_int (fplog, "#PME-nodes", cr->npmenodes, npme_f, &mm);
1823 npp = cr->nnodes;
1824 if (cr->npmenodes >= 0)
1826 npp -= cr->npmenodes;
1828 if (npp == npp_f)
1830 check_int (fplog, "#DD-cells[x]", dd_nc[XX], dd_nc_f[XX], &mm);
1831 check_int (fplog, "#DD-cells[y]", dd_nc[YY], dd_nc_f[YY], &mm);
1832 check_int (fplog, "#DD-cells[z]", dd_nc[ZZ], dd_nc_f[ZZ], &mm);
1836 if (mm)
1838 const char msg_version_difference[] =
1839 "The current Gromacs major & minor version are not identical to those that\n"
1840 "generated the checkpoint file. In principle Gromacs does not support\n"
1841 "continuation from checkpoints between different versions, so we advise\n"
1842 "against this. If you still want to try your luck we recommend that you use\n"
1843 "the -noappend flag to keep your output files from the two versions separate.\n"
1844 "This might also work around errors where the output fields in the energy\n"
1845 "file have changed between the different major & minor versions.\n";
1847 const char msg_mismatch_notice[] =
1848 "Gromacs patchlevel, binary or parallel settings differ from previous run.\n"
1849 "Continuation is exact, but not guaranteed to be binary identical.\n";
1851 const char msg_logdetails[] =
1852 "See the log file for details.\n";
1854 if (version_differs)
1856 fprintf(stderr, "%s%s\n", msg_version_difference, fplog ? msg_logdetails : "");
1858 if (fplog)
1860 fprintf(fplog, "%s\n", msg_version_difference);
1863 else
1865 /* Major & minor versions match at least, but something is different. */
1866 fprintf(stderr, "%s%s\n", msg_mismatch_notice, fplog ? msg_logdetails : "");
1867 if (fplog)
1869 fprintf(fplog, "%s\n", msg_mismatch_notice);
1875 static void read_checkpoint(const char *fn, FILE **pfplog,
1876 t_commrec *cr, ivec dd_nc,
1877 int eIntegrator, int *init_fep_state, gmx_int64_t *step, double *t,
1878 t_state *state, gmx_bool *bReadEkin,
1879 int *simulation_part,
1880 gmx_bool bAppendOutputFiles, gmx_bool bForceAppend)
1882 t_fileio *fp;
1883 int i, j, rc;
1884 int file_version;
1885 char *version, *btime, *buser, *bhost, *fprog, *ftime;
1886 int double_prec;
1887 char filename[STRLEN], buf[STEPSTRSIZE];
1888 int nppnodes, eIntegrator_f, nppnodes_f, npmenodes_f;
1889 ivec dd_nc_f;
1890 int natoms, ngtc, nnhpres, nhchainlength, nlambda, fflags, flags_eks, flags_enh, flags_dfh;
1891 int d;
1892 int ret;
1893 gmx_file_position_t *outputfiles;
1894 int nfiles;
1895 t_fileio *chksum_file;
1896 FILE * fplog = *pfplog;
1897 unsigned char digest[16];
1898 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1899 struct flock fl; /* don't initialize here: the struct order is OS
1900 dependent! */
1901 #endif
1903 const char *int_warn =
1904 "WARNING: The checkpoint file was generated with integrator %s,\n"
1905 " while the simulation uses integrator %s\n\n";
1907 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1908 fl.l_type = F_WRLCK;
1909 fl.l_whence = SEEK_SET;
1910 fl.l_start = 0;
1911 fl.l_len = 0;
1912 fl.l_pid = 0;
1913 #endif
1915 fp = gmx_fio_open(fn, "r");
1916 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
1917 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1918 &eIntegrator_f, simulation_part, step, t,
1919 &nppnodes_f, dd_nc_f, &npmenodes_f,
1920 &natoms, &ngtc, &nnhpres, &nhchainlength, &nlambda,
1921 &fflags, &flags_eks, &flags_enh, &flags_dfh,
1922 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
1924 if (bAppendOutputFiles &&
1925 file_version >= 13 && double_prec != GMX_CPT_BUILD_DP)
1927 gmx_fatal(FARGS, "Output file appending requested, but the code and checkpoint file precision (single/double) don't match");
1930 if (cr == NULL || MASTER(cr))
1932 fprintf(stderr, "\nReading checkpoint file %s generated: %s\n\n",
1933 fn, ftime);
1936 /* This will not be written if we do appending, since fplog is still NULL then */
1937 if (fplog)
1939 fprintf(fplog, "\n");
1940 fprintf(fplog, "Reading checkpoint file %s\n", fn);
1941 fprintf(fplog, " file generated by: %s\n", fprog);
1942 fprintf(fplog, " file generated at: %s\n", ftime);
1943 fprintf(fplog, " GROMACS build time: %s\n", btime);
1944 fprintf(fplog, " GROMACS build user: %s\n", buser);
1945 fprintf(fplog, " GROMACS build host: %s\n", bhost);
1946 fprintf(fplog, " GROMACS double prec.: %d\n", double_prec);
1947 fprintf(fplog, " simulation part #: %d\n", *simulation_part);
1948 fprintf(fplog, " step: %s\n", gmx_step_str(*step, buf));
1949 fprintf(fplog, " time: %f\n", *t);
1950 fprintf(fplog, "\n");
1953 if (natoms != state->natoms)
1955 gmx_fatal(FARGS, "Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms", natoms, state->natoms);
1957 if (ngtc != state->ngtc)
1959 gmx_fatal(FARGS, "Checkpoint file is for a system of %d T-coupling groups, while the current system consists of %d T-coupling groups", ngtc, state->ngtc);
1961 if (nnhpres != state->nnhpres)
1963 gmx_fatal(FARGS, "Checkpoint file is for a system of %d NH-pressure-coupling variables, while the current system consists of %d NH-pressure-coupling variables", nnhpres, state->nnhpres);
1966 if (nlambda != state->dfhist.nlambda)
1968 gmx_fatal(FARGS, "Checkpoint file is for a system with %d lambda states, while the current system consists of %d lambda states", nlambda, state->dfhist.nlambda);
1971 init_gtc_state(state, state->ngtc, state->nnhpres, nhchainlength); /* need to keep this here to keep the tpr format working */
1972 /* write over whatever was read; we use the number of Nose-Hoover chains from the checkpoint */
1974 if (eIntegrator_f != eIntegrator)
1976 if (MASTER(cr))
1978 fprintf(stderr, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1980 if (bAppendOutputFiles)
1982 gmx_fatal(FARGS,
1983 "Output file appending requested, but input/checkpoint integrators do not match.\n"
1984 "Stopping the run to prevent you from ruining all your data...\n"
1985 "If you _really_ know what you are doing, try with the -noappend option.\n");
1987 if (fplog)
1989 fprintf(fplog, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1993 if (!PAR(cr))
1995 nppnodes = 1;
1996 cr->npmenodes = 0;
1998 else if (cr->nnodes == nppnodes_f + npmenodes_f)
2000 if (cr->npmenodes < 0)
2002 cr->npmenodes = npmenodes_f;
2004 nppnodes = cr->nnodes - cr->npmenodes;
2005 if (nppnodes == nppnodes_f)
2007 for (d = 0; d < DIM; d++)
2009 if (dd_nc[d] == 0)
2011 dd_nc[d] = dd_nc_f[d];
2016 else
2018 /* The number of PP nodes has not been set yet */
2019 nppnodes = -1;
2022 if (fflags != state->flags)
2025 if (MASTER(cr))
2027 if (bAppendOutputFiles)
2029 gmx_fatal(FARGS,
2030 "Output file appending requested, but input and checkpoint states are not identical.\n"
2031 "Stopping the run to prevent you from ruining all your data...\n"
2032 "You can try with the -noappend option, and get more info in the log file.\n");
2035 if (getenv("GMX_ALLOW_CPT_MISMATCH") == NULL)
2037 gmx_fatal(FARGS, "You seem to have switched ensemble, integrator, T and/or P-coupling algorithm between the cpt and tpr file. The recommended way of doing this is passing the cpt file to grompp (with option -t) instead of to mdrun. If you know what you are doing, you can override this error by setting the env.var. GMX_ALLOW_CPT_MISMATCH");
2039 else
2041 fprintf(stderr,
2042 "WARNING: The checkpoint state entries do not match the simulation,\n"
2043 " see the log file for details\n\n");
2047 if (fplog)
2049 print_flag_mismatch(fplog, state->flags, fflags);
2052 else
2054 if (MASTER(cr))
2056 check_match(fplog, version, btime, buser, bhost, double_prec, fprog,
2057 cr, nppnodes_f, npmenodes_f, dd_nc, dd_nc_f);
2060 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, fflags, state, NULL);
2061 *init_fep_state = state->fep_state; /* there should be a better way to do this than setting it here.
2062 Investigate for 5.0. */
2063 if (ret)
2065 cp_error();
2067 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2068 if (ret)
2070 cp_error();
2072 *bReadEkin = ((flags_eks & (1<<eeksEKINH)) || (flags_eks & (1<<eeksEKINF)) || (flags_eks & (1<<eeksEKINO)) ||
2073 ((flags_eks & (1<<eeksEKINSCALEF)) | (flags_eks & (1<<eeksEKINSCALEH)) | (flags_eks & (1<<eeksVSCALE))));
2075 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2076 flags_enh, &state->enerhist, NULL);
2077 if (ret)
2079 cp_error();
2082 if (file_version < 6)
2084 const char *warn = "Reading checkpoint file in old format, assuming that the run that generated this file started at step 0, if this is not the case the averages stored in the energy file will be incorrect.";
2086 fprintf(stderr, "\nWARNING: %s\n\n", warn);
2087 if (fplog)
2089 fprintf(fplog, "\nWARNING: %s\n\n", warn);
2091 state->enerhist.nsum = *step;
2092 state->enerhist.nsum_sim = *step;
2095 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2096 if (ret)
2098 cp_error();
2101 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2102 if (ret)
2104 cp_error();
2107 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2108 if (ret)
2110 cp_error();
2113 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, NULL, file_version);
2114 if (ret)
2116 cp_error();
2119 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2120 if (ret)
2122 cp_error();
2124 if (gmx_fio_close(fp) != 0)
2126 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2129 sfree(fprog);
2130 sfree(ftime);
2131 sfree(btime);
2132 sfree(buser);
2133 sfree(bhost);
2135 /* If the user wants to append to output files,
2136 * we use the file pointer positions of the output files stored
2137 * in the checkpoint file and truncate the files such that any frames
2138 * written after the checkpoint time are removed.
2139 * All files are md5sum checked such that we can be sure that
2140 * we do not truncate other (maybe imprortant) files.
2142 if (bAppendOutputFiles)
2144 if (fn2ftp(outputfiles[0].filename) != efLOG)
2146 /* make sure first file is log file so that it is OK to use it for
2147 * locking
2149 gmx_fatal(FARGS, "The first output file should always be the log "
2150 "file but instead is: %s. Cannot do appending because of this condition.", outputfiles[0].filename);
2152 for (i = 0; i < nfiles; i++)
2154 if (outputfiles[i].offset < 0)
2156 gmx_fatal(FARGS, "The original run wrote a file called '%s' which "
2157 "is larger than 2 GB, but mdrun did not support large file"
2158 " offsets. Can not append. Run mdrun with -noappend",
2159 outputfiles[i].filename);
2161 #ifdef GMX_FAHCORE
2162 chksum_file = gmx_fio_open(outputfiles[i].filename, "a");
2164 #else
2165 chksum_file = gmx_fio_open(outputfiles[i].filename, "r+");
2167 /* lock log file */
2168 if (i == 0)
2170 /* Note that there are systems where the lock operation
2171 * will succeed, but a second process can also lock the file.
2172 * We should probably try to detect this.
2174 #if defined __native_client__
2175 errno = ENOSYS;
2176 if (1)
2178 #elif defined GMX_NATIVE_WINDOWS
2179 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX) == -1)
2180 #else
2181 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl) == -1)
2182 #endif
2184 if (errno == ENOSYS)
2186 if (!bForceAppend)
2188 gmx_fatal(FARGS, "File locking is not supported on this system. Use -noappend or specify -append explicitly to append anyhow.");
2190 else
2192 fprintf(stderr, "\nNOTE: File locking is not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2193 if (fplog)
2195 fprintf(fplog, "\nNOTE: File locking not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2199 else if (errno == EACCES || errno == EAGAIN)
2201 gmx_fatal(FARGS, "Failed to lock: %s. Already running "
2202 "simulation?", outputfiles[i].filename);
2204 else
2206 gmx_fatal(FARGS, "Failed to lock: %s. %s.",
2207 outputfiles[i].filename, strerror(errno));
2212 /* compute md5 chksum */
2213 if (outputfiles[i].chksum_size != -1)
2215 if (gmx_fio_get_file_md5(chksum_file, outputfiles[i].offset,
2216 digest) != outputfiles[i].chksum_size) /*at the end of the call the file position is at the end of the file*/
2218 gmx_fatal(FARGS, "Can't read %d bytes of '%s' to compute checksum. The file has been replaced or its contents have been modified. Cannot do appending because of this condition.",
2219 outputfiles[i].chksum_size,
2220 outputfiles[i].filename);
2223 if (i == 0) /*log file needs to be seeked in case we need to truncate (other files are truncated below)*/
2225 if (gmx_fio_seek(chksum_file, outputfiles[i].offset))
2227 gmx_fatal(FARGS, "Seek error! Failed to truncate log-file: %s.", strerror(errno));
2230 #endif
2232 if (i == 0) /*open log file here - so that lock is never lifted
2233 after chksum is calculated */
2235 *pfplog = gmx_fio_getfp(chksum_file);
2237 else
2239 gmx_fio_close(chksum_file);
2241 #ifndef GMX_FAHCORE
2242 /* compare md5 chksum */
2243 if (outputfiles[i].chksum_size != -1 &&
2244 memcmp(digest, outputfiles[i].chksum, 16) != 0)
2246 if (debug)
2248 fprintf(debug, "chksum for %s: ", outputfiles[i].filename);
2249 for (j = 0; j < 16; j++)
2251 fprintf(debug, "%02x", digest[j]);
2253 fprintf(debug, "\n");
2255 gmx_fatal(FARGS, "Checksum wrong for '%s'. The file has been replaced or its contents have been modified. Cannot do appending because of this condition.",
2256 outputfiles[i].filename);
2258 #endif
2261 if (i != 0) /*log file is already seeked to correct position */
2263 #ifdef GMX_NATIVE_WINDOWS
2264 rc = gmx_wintruncate(outputfiles[i].filename, outputfiles[i].offset);
2265 #else
2266 rc = truncate(outputfiles[i].filename, outputfiles[i].offset);
2267 #endif
2268 if (rc != 0)
2270 gmx_fatal(FARGS, "Truncation of file %s failed. Cannot do appending because of this failure.", outputfiles[i].filename);
2276 sfree(outputfiles);
2280 void load_checkpoint(const char *fn, FILE **fplog,
2281 t_commrec *cr, ivec dd_nc,
2282 t_inputrec *ir, t_state *state,
2283 gmx_bool *bReadEkin,
2284 gmx_bool bAppend, gmx_bool bForceAppend)
2286 gmx_int64_t step;
2287 double t;
2289 if (SIMMASTER(cr))
2291 /* Read the state from the checkpoint file */
2292 read_checkpoint(fn, fplog,
2293 cr, dd_nc,
2294 ir->eI, &(ir->fepvals->init_fep_state), &step, &t, state, bReadEkin,
2295 &ir->simulation_part, bAppend, bForceAppend);
2297 if (PAR(cr))
2299 gmx_bcast(sizeof(cr->npmenodes), &cr->npmenodes, cr);
2300 gmx_bcast(DIM*sizeof(dd_nc[0]), dd_nc, cr);
2301 gmx_bcast(sizeof(step), &step, cr);
2302 gmx_bcast(sizeof(*bReadEkin), bReadEkin, cr);
2304 ir->bContinuation = TRUE;
2305 if (ir->nsteps >= 0)
2307 ir->nsteps += ir->init_step - step;
2309 ir->init_step = step;
2310 ir->simulation_part += 1;
2313 static void read_checkpoint_data(t_fileio *fp, int *simulation_part,
2314 gmx_int64_t *step, double *t, t_state *state,
2315 int *nfiles, gmx_file_position_t **outputfiles)
2317 int file_version;
2318 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2319 int double_prec;
2320 int eIntegrator;
2321 int nppnodes, npme;
2322 ivec dd_nc;
2323 int flags_eks, flags_enh, flags_dfh;
2324 int nfiles_loc;
2325 gmx_file_position_t *files_loc = NULL;
2326 int ret;
2328 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2329 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2330 &eIntegrator, simulation_part, step, t, &nppnodes, dd_nc, &npme,
2331 &state->natoms, &state->ngtc, &state->nnhpres, &state->nhchainlength,
2332 &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
2333 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
2334 ret =
2335 do_cpt_state(gmx_fio_getxdr(fp), TRUE, state->flags, state, NULL);
2336 if (ret)
2338 cp_error();
2340 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2341 if (ret)
2343 cp_error();
2345 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2346 flags_enh, &state->enerhist, NULL);
2347 if (ret)
2349 cp_error();
2351 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2352 if (ret)
2354 cp_error();
2357 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2358 if (ret)
2360 cp_error();
2363 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2364 if (ret)
2366 cp_error();
2369 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE,
2370 outputfiles != NULL ? outputfiles : &files_loc,
2371 outputfiles != NULL ? nfiles : &nfiles_loc,
2372 NULL, file_version);
2373 if (files_loc != NULL)
2375 sfree(files_loc);
2378 if (ret)
2380 cp_error();
2383 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2384 if (ret)
2386 cp_error();
2389 sfree(fprog);
2390 sfree(ftime);
2391 sfree(btime);
2392 sfree(buser);
2393 sfree(bhost);
2396 void
2397 read_checkpoint_state(const char *fn, int *simulation_part,
2398 gmx_int64_t *step, double *t, t_state *state)
2400 t_fileio *fp;
2402 fp = gmx_fio_open(fn, "r");
2403 read_checkpoint_data(fp, simulation_part, step, t, state, NULL, NULL);
2404 if (gmx_fio_close(fp) != 0)
2406 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2410 void read_checkpoint_trxframe(t_fileio *fp, t_trxframe *fr)
2412 /* This next line is nasty because the sub-structures of t_state
2413 * cannot be assumed to be zeroed (or even initialized in ways the
2414 * rest of the code might assume). Using snew would be better, but
2415 * this will all go away for 5.0. */
2416 t_state state;
2417 int simulation_part;
2418 gmx_int64_t step;
2419 double t;
2421 init_state(&state, 0, 0, 0, 0, 0);
2423 read_checkpoint_data(fp, &simulation_part, &step, &t, &state, NULL, NULL);
2425 fr->natoms = state.natoms;
2426 fr->bTitle = FALSE;
2427 fr->bStep = TRUE;
2428 fr->step = gmx_int64_to_int(step,
2429 "conversion of checkpoint to trajectory");
2430 fr->bTime = TRUE;
2431 fr->time = t;
2432 fr->bLambda = TRUE;
2433 fr->lambda = state.lambda[efptFEP];
2434 fr->fep_state = state.fep_state;
2435 fr->bAtoms = FALSE;
2436 fr->bX = (state.flags & (1<<estX));
2437 if (fr->bX)
2439 fr->x = state.x;
2440 state.x = NULL;
2442 fr->bV = (state.flags & (1<<estV));
2443 if (fr->bV)
2445 fr->v = state.v;
2446 state.v = NULL;
2448 fr->bF = FALSE;
2449 fr->bBox = (state.flags & (1<<estBOX));
2450 if (fr->bBox)
2452 copy_mat(state.box, fr->box);
2454 done_state(&state);
2457 void list_checkpoint(const char *fn, FILE *out)
2459 t_fileio *fp;
2460 int file_version;
2461 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2462 int double_prec;
2463 int eIntegrator, simulation_part, nppnodes, npme;
2464 gmx_int64_t step;
2465 double t;
2466 ivec dd_nc;
2467 t_state state;
2468 int flags_eks, flags_enh, flags_dfh;
2469 int indent;
2470 int i, j;
2471 int ret;
2472 gmx_file_position_t *outputfiles;
2473 int nfiles;
2475 init_state(&state, -1, -1, -1, -1, 0);
2477 fp = gmx_fio_open(fn, "r");
2478 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2479 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2480 &eIntegrator, &simulation_part, &step, &t, &nppnodes, dd_nc, &npme,
2481 &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
2482 &(state.dfhist.nlambda), &state.flags,
2483 &flags_eks, &flags_enh, &flags_dfh, &state.edsamstate.nED,
2484 &state.swapstate.eSwapCoords, out);
2485 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, state.flags, &state, out);
2486 if (ret)
2488 cp_error();
2490 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state.ekinstate, out);
2491 if (ret)
2493 cp_error();
2495 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2496 flags_enh, &state.enerhist, out);
2498 if (ret == 0)
2500 ret = do_cpt_df_hist(gmx_fio_getxdr(fp),
2501 flags_dfh, &state.dfhist, out);
2504 if (ret == 0)
2506 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state.edsamstate, out);
2509 if (ret == 0)
2511 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state.swapstate, out);
2514 if (ret == 0)
2516 do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, out, file_version);
2519 if (ret == 0)
2521 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2524 if (ret)
2526 cp_warning(out);
2528 if (gmx_fio_close(fp) != 0)
2530 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2533 done_state(&state);
2537 static gmx_bool exist_output_file(const char *fnm_cp, int nfile, const t_filenm fnm[])
2539 int i;
2541 /* Check if the output file name stored in the checkpoint file
2542 * is one of the output file names of mdrun.
2544 i = 0;
2545 while (i < nfile &&
2546 !(is_output(&fnm[i]) && strcmp(fnm_cp, fnm[i].fns[0]) == 0))
2548 i++;
2551 return (i < nfile && gmx_fexist(fnm_cp));
2554 /* This routine cannot print tons of data, since it is called before the log file is opened. */
2555 gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
2556 gmx_int64_t *cpt_step, t_commrec *cr,
2557 gmx_bool bAppendReq,
2558 int nfile, const t_filenm fnm[],
2559 const char *part_suffix, gmx_bool *bAddPart)
2561 t_fileio *fp;
2562 gmx_int64_t step = 0;
2563 double t;
2564 /* This next line is nasty because the sub-structures of t_state
2565 * cannot be assumed to be zeroed (or even initialized in ways the
2566 * rest of the code might assume). Using snew would be better, but
2567 * this will all go away for 5.0. */
2568 t_state state;
2569 int nfiles;
2570 gmx_file_position_t *outputfiles;
2571 int nexist, f;
2572 gmx_bool bAppend;
2573 char *fn, suf_up[STRLEN];
2575 bAppend = FALSE;
2577 if (SIMMASTER(cr))
2579 if (!gmx_fexist(filename) || (!(fp = gmx_fio_open(filename, "r")) ))
2581 *simulation_part = 0;
2583 else
2585 init_state(&state, 0, 0, 0, 0, 0);
2587 read_checkpoint_data(fp, simulation_part, &step, &t, &state,
2588 &nfiles, &outputfiles);
2589 if (gmx_fio_close(fp) != 0)
2591 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2593 done_state(&state);
2595 if (bAppendReq)
2597 nexist = 0;
2598 for (f = 0; f < nfiles; f++)
2600 if (exist_output_file(outputfiles[f].filename, nfile, fnm))
2602 nexist++;
2605 if (nexist == nfiles)
2607 bAppend = bAppendReq;
2609 else if (nexist > 0)
2611 fprintf(stderr,
2612 "Output file appending has been requested,\n"
2613 "but some output files listed in the checkpoint file %s\n"
2614 "are not present or are named differently by the current program:\n",
2615 filename);
2616 fprintf(stderr, "output files present:");
2617 for (f = 0; f < nfiles; f++)
2619 if (exist_output_file(outputfiles[f].filename,
2620 nfile, fnm))
2622 fprintf(stderr, " %s", outputfiles[f].filename);
2625 fprintf(stderr, "\n");
2626 fprintf(stderr, "output files not present or named differently:");
2627 for (f = 0; f < nfiles; f++)
2629 if (!exist_output_file(outputfiles[f].filename,
2630 nfile, fnm))
2632 fprintf(stderr, " %s", outputfiles[f].filename);
2635 fprintf(stderr, "\n");
2637 gmx_fatal(FARGS, "File appending requested, but %d of the %d output files are not present or are named differently", nfiles-nexist, nfiles);
2641 if (bAppend)
2643 if (nfiles == 0)
2645 gmx_fatal(FARGS, "File appending requested, but no output file information is stored in the checkpoint file");
2647 fn = outputfiles[0].filename;
2648 if (strlen(fn) < 4 ||
2649 gmx_strcasecmp(fn+strlen(fn)-4, ftp2ext(efLOG)) == 0)
2651 gmx_fatal(FARGS, "File appending requested, but the log file is not the first file listed in the checkpoint file");
2653 /* Set bAddPart to whether the suffix string '.part' is present
2654 * in the log file name.
2656 strcpy(suf_up, part_suffix);
2657 upstring(suf_up);
2658 *bAddPart = (strstr(fn, part_suffix) != NULL ||
2659 strstr(fn, suf_up) != NULL);
2662 sfree(outputfiles);
2665 if (PAR(cr))
2667 gmx_bcast(sizeof(*simulation_part), simulation_part, cr);
2669 if (*simulation_part > 0 && bAppendReq)
2671 gmx_bcast(sizeof(bAppend), &bAppend, cr);
2672 gmx_bcast(sizeof(*bAddPart), bAddPart, cr);
2675 if (NULL != cpt_step)
2677 *cpt_step = step;
2680 return bAppend;