2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2008,2009,2010,2011,2012,2013,2014,2015, 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 #include "gromacs/legacyheaders/checkpoint.h"
49 #ifdef GMX_NATIVE_WINDOWS
51 #include <sys/locking.h>
54 #include "buildinfo.h"
55 #include "gromacs/fileio/filenm.h"
56 #include "gromacs/fileio/gmxfio.h"
57 #include "gromacs/fileio/xdr_datatype.h"
58 #include "gromacs/fileio/xdrf.h"
59 #include "gromacs/legacyheaders/copyrite.h"
60 #include "gromacs/legacyheaders/names.h"
61 #include "gromacs/legacyheaders/network.h"
62 #include "gromacs/legacyheaders/txtdump.h"
63 #include "gromacs/legacyheaders/typedefs.h"
64 #include "gromacs/legacyheaders/types/commrec.h"
65 #include "gromacs/math/vec.h"
66 #include "gromacs/utility/baseversion.h"
67 #include "gromacs/utility/cstringutil.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/futil.h"
70 #include "gromacs/utility/smalloc.h"
71 #include "gromacs/utility/sysinfo.h"
77 #define CPT_MAGIC1 171817
78 #define CPT_MAGIC2 171819
79 #define CPTSTRLEN 1024
82 #define GMX_CPT_BUILD_DP 1
84 #define GMX_CPT_BUILD_DP 0
87 /* cpt_version should normally only be changed
88 * when the header of footer format changes.
89 * The state data format itself is backward and forward compatible.
90 * But old code can not read a new entry that is present in the file
91 * (but can read a new format when new entries are not present).
93 static const int cpt_version = 16;
96 const char *est_names[estNR] =
99 "box", "box-rel", "box-v", "pres_prev",
100 "nosehoover-xi", "thermostat-integral",
101 "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
102 "disre_initf", "disre_rm3tav",
103 "orire_initf", "orire_Dtav",
104 "svir_prev", "nosehoover-vxi", "v_eta", "vol0", "nhpres_xi", "nhpres_vxi", "fvir_prev", "fep_state", "MC-rng", "MC-rng-i"
108 eeksEKIN_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksEKINF, eeksEKINO, eeksEKINSCALEF, eeksEKINSCALEH, eeksVSCALE, eeksEKINTOTAL, eeksNR
111 const char *eeks_names[eeksNR] =
113 "Ekin_n", "Ekinh", "dEkindlambda", "mv_cos",
114 "Ekinf", "Ekinh_old", "EkinScaleF_NHC", "EkinScaleH_NHC", "Vscale_NHC", "Ekin_Total"
118 eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
119 eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
120 eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM,
121 eenhENERGY_DELTA_H_NN,
122 eenhENERGY_DELTA_H_LIST,
123 eenhENERGY_DELTA_H_STARTTIME,
124 eenhENERGY_DELTA_H_STARTLAMBDA,
128 const char *eenh_names[eenhNR] =
130 "energy_n", "energy_aver", "energy_sum", "energy_nsum",
131 "energy_sum_sim", "energy_nsum_sim",
132 "energy_nsteps", "energy_nsteps_sim",
134 "energy_delta_h_list",
135 "energy_delta_h_start_time",
136 "energy_delta_h_start_lambda"
139 /* free energy history variables -- need to be preserved over checkpoint */
141 edfhBEQUIL, edfhNATLAMBDA, edfhWLHISTO, edfhWLDELTA, edfhSUMWEIGHTS, edfhSUMDG, edfhSUMMINVAR, edfhSUMVAR,
142 edfhACCUMP, edfhACCUMM, edfhACCUMP2, edfhACCUMM2, edfhTIJ, edfhTIJEMP, edfhNR
144 /* free energy history variable names */
145 const char *edfh_names[edfhNR] =
147 "bEquilibrated", "N_at_state", "Wang-Landau Histogram", "Wang-Landau Delta", "Weights", "Free Energies", "minvar", "variance",
148 "accumulated_plus", "accumulated_minus", "accumulated_plus_2", "accumulated_minus_2", "Tij", "Tij_empirical"
152 ecprREAL, ecprRVEC, ecprMATRIX
156 cptpEST, cptpEEKS, cptpEENH, cptpEDFH
158 /* enums for the different components of checkpoint variables, replacing the hard coded ones.
159 cptpEST - state variables.
160 cptpEEKS - Kinetic energy state variables.
161 cptpEENH - Energy history state variables.
162 cptpEDFH - free energy history variables.
166 static const char *st_names(int cptp, int ecpt)
170 case cptpEST: return est_names [ecpt]; break;
171 case cptpEEKS: return eeks_names[ecpt]; break;
172 case cptpEENH: return eenh_names[ecpt]; break;
173 case cptpEDFH: return edfh_names[ecpt]; break;
179 static void cp_warning(FILE *fp)
181 fprintf(fp, "\nWARNING: Checkpoint file is corrupted or truncated\n\n");
184 static void cp_error()
186 gmx_fatal(FARGS, "Checkpoint file corrupted/truncated, or maybe you are out of disk space?");
189 static void do_cpt_string_err(XDR *xd, gmx_bool bRead, const char *desc, char **s, FILE *list)
197 res = xdr_string(xd, s, CPTSTRLEN);
204 fprintf(list, "%s = %s\n", desc, *s);
209 static int do_cpt_int(XDR *xd, const char *desc, int *i, FILE *list)
213 res = xdr_int(xd, i);
220 fprintf(list, "%s = %d\n", desc, *i);
225 static int do_cpt_u_chars(XDR *xd, const char *desc, int n, unsigned char *i, FILE *list)
231 fprintf(list, "%s = ", desc);
233 for (j = 0; j < n && res; j++)
235 res &= xdr_u_char(xd, &i[j]);
238 fprintf(list, "%02x", i[j]);
253 static void do_cpt_int_err(XDR *xd, const char *desc, int *i, FILE *list)
255 if (do_cpt_int(xd, desc, i, list) < 0)
261 static void do_cpt_step_err(XDR *xd, const char *desc, gmx_int64_t *i, FILE *list)
264 char buf[STEPSTRSIZE];
266 res = xdr_int64(xd, i);
273 fprintf(list, "%s = %s\n", desc, gmx_step_str(*i, buf));
277 static void do_cpt_double_err(XDR *xd, const char *desc, double *f, FILE *list)
281 res = xdr_double(xd, f);
288 fprintf(list, "%s = %f\n", desc, *f);
292 static void do_cpt_real_err(XDR *xd, real *f)
297 res = xdr_double(xd, f);
299 res = xdr_float(xd, f);
307 static void do_cpt_n_rvecs_err(XDR *xd, const char *desc, int n, rvec f[], FILE *list)
311 for (i = 0; i < n; i++)
313 for (j = 0; j < DIM; j++)
315 do_cpt_real_err(xd, &f[i][j]);
321 pr_rvecs(list, 0, desc, f, n);
325 /* If nval >= 0, nval is used; on read this should match the passed value.
326 * If nval n<0, *nptr is used; on read the value is stored in nptr
328 static int do_cpte_reals_low(XDR *xd, int cptp, int ecpt, int sflags,
329 int nval, int *nptr, real **v,
330 FILE *list, int erealtype)
334 int dtc = xdr_datatype_float;
336 int dtc = xdr_datatype_double;
338 real *vp, *va = NULL;
353 gmx_incons("*ntpr=NULL in do_cpte_reals_low");
358 res = xdr_int(xd, &nf);
369 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), nval, nf);
378 res = xdr_int(xd, &dt);
385 fprintf(stderr, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
386 st_names(cptp, ecpt), xdr_datatype_names[dtc],
387 xdr_datatype_names[dt]);
389 if (list || !(sflags & (1<<ecpt)))
402 if (dt == xdr_datatype_float)
404 if (dtc == xdr_datatype_float)
412 res = xdr_vector(xd, (char *)vf, nf,
413 (unsigned int)sizeof(float), (xdrproc_t)xdr_float);
418 if (dtc != xdr_datatype_float)
420 for (i = 0; i < nf; i++)
429 if (dtc == xdr_datatype_double)
431 /* cppcheck-suppress invalidPointerCast
432 * Only executed if real is anyhow double */
439 res = xdr_vector(xd, (char *)vd, nf,
440 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
445 if (dtc != xdr_datatype_double)
447 for (i = 0; i < nf; i++)
460 pr_reals(list, 0, st_names(cptp, ecpt), vp, nf);
463 pr_rvecs(list, 0, st_names(cptp, ecpt), (rvec *)vp, nf/3);
466 gmx_incons("Unknown checkpoint real type");
478 /* This function stores n along with the reals for reading,
479 * but on reading it assumes that n matches the value in the checkpoint file,
480 * a fatal error is generated when this is not the case.
482 static int do_cpte_reals(XDR *xd, int cptp, int ecpt, int sflags,
483 int n, real **v, FILE *list)
485 return do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, v, list, ecprREAL);
488 /* This function does the same as do_cpte_reals,
489 * except that on reading it ignores the passed value of *n
490 * and stored the value read from the checkpoint file in *n.
492 static int do_cpte_n_reals(XDR *xd, int cptp, int ecpt, int sflags,
493 int *n, real **v, FILE *list)
495 return do_cpte_reals_low(xd, cptp, ecpt, sflags, -1, n, v, list, ecprREAL);
498 static int do_cpte_real(XDR *xd, int cptp, int ecpt, int sflags,
501 return do_cpte_reals_low(xd, cptp, ecpt, sflags, 1, NULL, &r, list, ecprREAL);
504 static int do_cpte_ints(XDR *xd, int cptp, int ecpt, int sflags,
505 int n, int **v, FILE *list)
508 int dtc = xdr_datatype_int;
513 res = xdr_int(xd, &nf);
518 if (list == NULL && v != NULL && nf != n)
520 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
523 res = xdr_int(xd, &dt);
530 gmx_fatal(FARGS, "Type mismatch for state entry %s, code type is %s, file type is %s\n",
531 st_names(cptp, ecpt), xdr_datatype_names[dtc],
532 xdr_datatype_names[dt]);
534 if (list || !(sflags & (1<<ecpt)) || v == NULL)
547 res = xdr_vector(xd, (char *)vp, nf,
548 (unsigned int)sizeof(int), (xdrproc_t)xdr_int);
555 pr_ivec(list, 0, st_names(cptp, ecpt), vp, nf, TRUE);
565 static int do_cpte_int(XDR *xd, int cptp, int ecpt, int sflags,
568 return do_cpte_ints(xd, cptp, ecpt, sflags, 1, &i, list);
571 static int do_cpte_doubles(XDR *xd, int cptp, int ecpt, int sflags,
572 int n, double **v, FILE *list)
575 int dtc = xdr_datatype_double;
576 double *vp, *va = NULL;
580 res = xdr_int(xd, &nf);
585 if (list == NULL && nf != n)
587 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
590 res = xdr_int(xd, &dt);
597 gmx_fatal(FARGS, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
598 st_names(cptp, ecpt), xdr_datatype_names[dtc],
599 xdr_datatype_names[dt]);
601 if (list || !(sflags & (1<<ecpt)))
614 res = xdr_vector(xd, (char *)vp, nf,
615 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
622 pr_doubles(list, 0, st_names(cptp, ecpt), vp, nf);
632 static int do_cpte_double(XDR *xd, int cptp, int ecpt, int sflags,
633 double *r, FILE *list)
635 return do_cpte_doubles(xd, cptp, ecpt, sflags, 1, &r, list);
639 static int do_cpte_rvecs(XDR *xd, int cptp, int ecpt, int sflags,
640 int n, rvec **v, FILE *list)
642 return do_cpte_reals_low(xd, cptp, ecpt, sflags,
643 n*DIM, NULL, (real **)v, list, ecprRVEC);
646 static int do_cpte_matrix(XDR *xd, int cptp, int ecpt, int sflags,
647 matrix v, FILE *list)
652 vr = (real *)&(v[0][0]);
653 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
654 DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
656 if (list && ret == 0)
658 pr_rvecs(list, 0, st_names(cptp, ecpt), v, DIM);
665 static int do_cpte_nmatrix(XDR *xd, int cptp, int ecpt, int sflags,
666 int n, real **v, FILE *list)
670 char name[CPTSTRLEN];
677 for (i = 0; i < n; i++)
679 reti = do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, &(v[i]), NULL, ecprREAL);
680 if (list && reti == 0)
682 sprintf(name, "%s[%d]", st_names(cptp, ecpt), i);
683 pr_reals(list, 0, name, v[i], n);
693 static int do_cpte_matrices(XDR *xd, int cptp, int ecpt, int sflags,
694 int n, matrix **v, FILE *list)
697 matrix *vp, *va = NULL;
703 res = xdr_int(xd, &nf);
708 if (list == NULL && nf != n)
710 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
712 if (list || !(sflags & (1<<ecpt)))
725 snew(vr, nf*DIM*DIM);
726 for (i = 0; i < nf; i++)
728 for (j = 0; j < DIM; j++)
730 for (k = 0; k < DIM; k++)
732 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
736 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
737 nf*DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
738 for (i = 0; i < nf; i++)
740 for (j = 0; j < DIM; j++)
742 for (k = 0; k < DIM; k++)
744 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
750 if (list && ret == 0)
752 for (i = 0; i < nf; i++)
754 pr_rvecs(list, 0, st_names(cptp, ecpt), vp[i], DIM);
765 static void do_cpt_header(XDR *xd, gmx_bool bRead, int *file_version,
766 char **version, char **btime, char **buser, char **bhost,
768 char **fprog, char **ftime,
769 int *eIntegrator, int *simulation_part,
770 gmx_int64_t *step, double *t,
771 int *nnodes, int *dd_nc, int *npme,
772 int *natoms, int *ngtc, int *nnhpres, int *nhchainlength,
773 int *nlambda, int *flags_state,
774 int *flags_eks, int *flags_enh, int *flags_dfh,
775 int *nED, int *eSwapCoords,
791 res = xdr_int(xd, &magic);
794 gmx_fatal(FARGS, "The checkpoint file is empty/corrupted, or maybe you are out of disk space?");
796 if (magic != CPT_MAGIC1)
798 gmx_fatal(FARGS, "Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
799 "The checkpoint file is corrupted or not a checkpoint file",
805 gmx_gethostname(fhost, 255);
807 do_cpt_string_err(xd, bRead, "GROMACS version", version, list);
808 do_cpt_string_err(xd, bRead, "GROMACS build time", btime, list);
809 do_cpt_string_err(xd, bRead, "GROMACS build user", buser, list);
810 do_cpt_string_err(xd, bRead, "GROMACS build host", bhost, list);
811 do_cpt_string_err(xd, bRead, "generating program", fprog, list);
812 do_cpt_string_err(xd, bRead, "generation time", ftime, list);
813 *file_version = cpt_version;
814 do_cpt_int_err(xd, "checkpoint file version", file_version, list);
815 if (*file_version > cpt_version)
817 gmx_fatal(FARGS, "Attempting to read a checkpoint file of version %d with code of version %d\n", *file_version, cpt_version);
819 if (*file_version >= 13)
821 do_cpt_int_err(xd, "GROMACS double precision", double_prec, list);
827 if (*file_version >= 12)
829 do_cpt_string_err(xd, bRead, "generating host", &fhost, list);
835 do_cpt_int_err(xd, "#atoms", natoms, list);
836 do_cpt_int_err(xd, "#T-coupling groups", ngtc, list);
837 if (*file_version >= 10)
839 do_cpt_int_err(xd, "#Nose-Hoover T-chains", nhchainlength, list);
845 if (*file_version >= 11)
847 do_cpt_int_err(xd, "#Nose-Hoover T-chains for barostat ", nnhpres, list);
853 if (*file_version >= 14)
855 do_cpt_int_err(xd, "# of total lambda states ", nlambda, list);
861 do_cpt_int_err(xd, "integrator", eIntegrator, list);
862 if (*file_version >= 3)
864 do_cpt_int_err(xd, "simulation part #", simulation_part, list);
868 *simulation_part = 1;
870 if (*file_version >= 5)
872 do_cpt_step_err(xd, "step", step, list);
876 do_cpt_int_err(xd, "step", &idum, list);
879 do_cpt_double_err(xd, "t", t, list);
880 do_cpt_int_err(xd, "#PP-ranks", nnodes, list);
882 do_cpt_int_err(xd, "dd_nc[x]", dd_nc ? &(dd_nc[0]) : &idum, list);
883 do_cpt_int_err(xd, "dd_nc[y]", dd_nc ? &(dd_nc[1]) : &idum, list);
884 do_cpt_int_err(xd, "dd_nc[z]", dd_nc ? &(dd_nc[2]) : &idum, list);
885 do_cpt_int_err(xd, "#PME-only ranks", npme, list);
886 do_cpt_int_err(xd, "state flags", flags_state, list);
887 if (*file_version >= 4)
889 do_cpt_int_err(xd, "ekin data flags", flags_eks, list);
890 do_cpt_int_err(xd, "energy history flags", flags_enh, list);
895 *flags_enh = (*flags_state >> (estORIRE_DTAV+1));
896 *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
897 (1<<(estORIRE_DTAV+2)) |
898 (1<<(estORIRE_DTAV+3))));
900 if (*file_version >= 14)
902 do_cpt_int_err(xd, "df history flags", flags_dfh, list);
909 if (*file_version >= 15)
911 do_cpt_int_err(xd, "ED data sets", nED, list);
917 if (*file_version >= 16)
919 do_cpt_int_err(xd, "swap", eSwapCoords, list);
923 static int do_cpt_footer(XDR *xd, int file_version)
928 if (file_version >= 2)
931 res = xdr_int(xd, &magic);
936 if (magic != CPT_MAGIC2)
945 static int do_cpt_state(XDR *xd, gmx_bool bRead,
946 int fflags, t_state *state,
956 nnht = state->nhchainlength*state->ngtc;
957 nnhtp = state->nhchainlength*state->nnhpres;
959 if (bRead) /* we need to allocate space for dfhist if we are reading */
961 init_df_history(&state->dfhist, state->dfhist.nlambda);
964 sflags = state->flags;
965 for (i = 0; (i < estNR && ret == 0); i++)
971 case estLAMBDA: ret = do_cpte_reals(xd, cptpEST, i, sflags, efptNR, &(state->lambda), list); break;
972 case estFEPSTATE: ret = do_cpte_int (xd, cptpEST, i, sflags, &state->fep_state, list); break;
973 case estBOX: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box, list); break;
974 case estBOX_REL: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box_rel, list); break;
975 case estBOXV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->boxv, list); break;
976 case estPRES_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->pres_prev, list); break;
977 case estSVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->svir_prev, list); break;
978 case estFVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->fvir_prev, list); break;
979 case estNH_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_xi, list); break;
980 case estNH_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_vxi, list); break;
981 case estNHPRES_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_xi, list); break;
982 case estNHPRES_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_vxi, list); break;
983 case estTC_INT: ret = do_cpte_doubles(xd, cptpEST, i, sflags, state->ngtc, &state->therm_integral, list); break;
984 case estVETA: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->veta, list); break;
985 case estVOL0: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->vol0, list); break;
986 case estX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->x, list); break;
987 case estV: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->v, list); break;
988 case estSDX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->sd_X, list); break;
989 /* The RNG entries are no longer written,
990 * the next 4 lines are only for reading old files.
992 case estLD_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
993 case estLD_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
994 case estMC_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
995 case estMC_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
996 case estDISRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.disre_initf, list); break;
997 case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.ndisrepairs, &state->hist.disre_rm3tav, list); break;
998 case estORIRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.orire_initf, list); break;
999 case estORIRE_DTAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.norire_Dtav, &state->hist.orire_Dtav, list); break;
1001 gmx_fatal(FARGS, "Unknown state entry %d\n"
1002 "You are probably reading a new checkpoint file with old code", i);
1010 static int do_cpt_ekinstate(XDR *xd, int fflags, ekinstate_t *ekins,
1018 for (i = 0; (i < eeksNR && ret == 0); i++)
1020 if (fflags & (1<<i))
1025 case eeksEKIN_N: ret = do_cpte_int(xd, cptpEEKS, i, fflags, &ekins->ekin_n, list); break;
1026 case eeksEKINH: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh, list); break;
1027 case eeksEKINF: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinf, list); break;
1028 case eeksEKINO: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh_old, list); break;
1029 case eeksEKINTOTAL: ret = do_cpte_matrix(xd, cptpEEKS, i, fflags, ekins->ekin_total, list); break;
1030 case eeksEKINSCALEF: ret = do_cpte_doubles(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinscalef_nhc, list); break;
1031 case eeksVSCALE: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->vscale_nhc, list); break;
1032 case eeksEKINSCALEH: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->ekinscaleh_nhc, list); break;
1033 case eeksDEKINDL: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->dekindl, list); break;
1034 case eeksMVCOS: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->mvcos, list); break;
1036 gmx_fatal(FARGS, "Unknown ekin data state entry %d\n"
1037 "You are probably reading a new checkpoint file with old code", i);
1046 static int do_cpt_swapstate(XDR *xd, gmx_bool bRead, swapstate_t *swapstate, FILE *list)
1050 int swap_cpt_version = 1;
1053 if (eswapNO == swapstate->eSwapCoords)
1058 swapstate->bFromCpt = bRead;
1060 do_cpt_int_err(xd, "swap checkpoint version", &swap_cpt_version, list);
1061 do_cpt_int_err(xd, "swap coupling steps", &swapstate->nAverage, list);
1063 /* When reading, init_swapcoords has not been called yet,
1064 * so we have to allocate memory first. */
1066 for (ic = 0; ic < eCompNR; ic++)
1068 for (ii = 0; ii < eIonNR; ii++)
1072 do_cpt_int_err(xd, "swap requested atoms", &swapstate->nat_req[ic][ii], list);
1076 do_cpt_int_err(xd, "swap requested atoms p", swapstate->nat_req_p[ic][ii], list);
1081 do_cpt_int_err(xd, "swap influx netto", &swapstate->inflow_netto[ic][ii], list);
1085 do_cpt_int_err(xd, "swap influx netto p", swapstate->inflow_netto_p[ic][ii], list);
1088 if (bRead && (NULL == swapstate->nat_past[ic][ii]) )
1090 snew(swapstate->nat_past[ic][ii], swapstate->nAverage);
1093 for (j = 0; j < swapstate->nAverage; j++)
1097 do_cpt_int_err(xd, "swap past atom counts", &swapstate->nat_past[ic][ii][j], list);
1101 do_cpt_int_err(xd, "swap past atom counts p", &swapstate->nat_past_p[ic][ii][j], list);
1107 /* Ion flux per channel */
1108 for (ic = 0; ic < eChanNR; ic++)
1110 for (ii = 0; ii < eIonNR; ii++)
1114 do_cpt_int_err(xd, "channel flux", &swapstate->fluxfromAtoB[ic][ii], list);
1118 do_cpt_int_err(xd, "channel flux p", swapstate->fluxfromAtoB_p[ic][ii], list);
1123 /* Ion flux leakage */
1126 snew(swapstate->fluxleak, 1);
1128 do_cpt_int_err(xd, "flux leakage", swapstate->fluxleak, list);
1131 do_cpt_int_err(xd, "number of ions", &swapstate->nions, list);
1135 snew(swapstate->channel_label, swapstate->nions);
1136 snew(swapstate->comp_from, swapstate->nions);
1139 do_cpt_u_chars(xd, "channel history", swapstate->nions, swapstate->channel_label, list);
1140 do_cpt_u_chars(xd, "domain history", swapstate->nions, swapstate->comp_from, list);
1142 /* Save the last known whole positions to checkpoint
1143 * file to be able to also make multimeric channels whole in PBC */
1144 do_cpt_int_err(xd, "Ch0 atoms", &swapstate->nat[eChan0], list);
1145 do_cpt_int_err(xd, "Ch1 atoms", &swapstate->nat[eChan1], list);
1148 snew(swapstate->xc_old_whole[eChan0], swapstate->nat[eChan0]);
1149 snew(swapstate->xc_old_whole[eChan1], swapstate->nat[eChan1]);
1150 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], swapstate->xc_old_whole[eChan0], list);
1151 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], swapstate->xc_old_whole[eChan1], list);
1155 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], *swapstate->xc_old_whole_p[eChan0], list);
1156 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], *swapstate->xc_old_whole_p[eChan1], list);
1163 static int do_cpt_enerhist(XDR *xd, gmx_bool bRead,
1164 int fflags, energyhistory_t *enerhist,
1175 enerhist->nsteps = 0;
1177 enerhist->nsteps_sim = 0;
1178 enerhist->nsum_sim = 0;
1179 enerhist->dht = NULL;
1181 if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
1183 snew(enerhist->dht, 1);
1184 enerhist->dht->ndh = NULL;
1185 enerhist->dht->dh = NULL;
1186 enerhist->dht->start_lambda_set = FALSE;
1190 for (i = 0; (i < eenhNR && ret == 0); i++)
1192 if (fflags & (1<<i))
1196 case eenhENERGY_N: ret = do_cpte_int(xd, cptpEENH, i, fflags, &enerhist->nener, list); break;
1197 case eenhENERGY_AVER: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_ave, list); break;
1198 case eenhENERGY_SUM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum, list); break;
1199 case eenhENERGY_NSUM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum, list); break;
1200 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum_sim, list); break;
1201 case eenhENERGY_NSUM_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum_sim, list); break;
1202 case eenhENERGY_NSTEPS: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps, list); break;
1203 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps_sim, list); break;
1204 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd, eenh_names[i], &(enerhist->dht->nndh), list);
1205 if (bRead) /* now allocate memory for it */
1207 snew(enerhist->dht->dh, enerhist->dht->nndh);
1208 snew(enerhist->dht->ndh, enerhist->dht->nndh);
1209 for (j = 0; j < enerhist->dht->nndh; j++)
1211 enerhist->dht->ndh[j] = 0;
1212 enerhist->dht->dh[j] = NULL;
1216 case eenhENERGY_DELTA_H_LIST:
1217 for (j = 0; j < enerhist->dht->nndh; j++)
1219 ret = do_cpte_n_reals(xd, cptpEENH, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list);
1222 case eenhENERGY_DELTA_H_STARTTIME:
1223 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_time), list); break;
1224 case eenhENERGY_DELTA_H_STARTLAMBDA:
1225 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_lambda), list); break;
1227 gmx_fatal(FARGS, "Unknown energy history entry %d\n"
1228 "You are probably reading a new checkpoint file with old code", i);
1233 if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
1235 /* Assume we have an old file format and copy sum to sum_sim */
1236 srenew(enerhist->ener_sum_sim, enerhist->nener);
1237 for (i = 0; i < enerhist->nener; i++)
1239 enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1243 if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1244 !(fflags & (1<<eenhENERGY_NSTEPS)))
1246 /* Assume we have an old file format and copy nsum to nsteps */
1247 enerhist->nsteps = enerhist->nsum;
1249 if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1250 !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1252 /* Assume we have an old file format and copy nsum to nsteps */
1253 enerhist->nsteps_sim = enerhist->nsum_sim;
1259 static int do_cpt_df_hist(XDR *xd, int fflags, df_history_t *dfhist, FILE *list)
1264 nlambda = dfhist->nlambda;
1267 for (i = 0; (i < edfhNR && ret == 0); i++)
1269 if (fflags & (1<<i))
1273 case edfhBEQUIL: ret = do_cpte_int(xd, cptpEDFH, i, fflags, &dfhist->bEquil, list); break;
1274 case edfhNATLAMBDA: ret = do_cpte_ints(xd, cptpEDFH, i, fflags, nlambda, &dfhist->n_at_lam, list); break;
1275 case edfhWLHISTO: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->wl_histo, list); break;
1276 case edfhWLDELTA: ret = do_cpte_real(xd, cptpEDFH, i, fflags, &dfhist->wl_delta, list); break;
1277 case edfhSUMWEIGHTS: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_weights, list); break;
1278 case edfhSUMDG: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_dg, list); break;
1279 case edfhSUMMINVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_minvar, list); break;
1280 case edfhSUMVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_variance, list); break;
1281 case edfhACCUMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p, list); break;
1282 case edfhACCUMM: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m, list); break;
1283 case edfhACCUMP2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p2, list); break;
1284 case edfhACCUMM2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m2, list); break;
1285 case edfhTIJ: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij, list); break;
1286 case edfhTIJEMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij_empirical, list); break;
1289 gmx_fatal(FARGS, "Unknown df history entry %d\n"
1290 "You are probably reading a new checkpoint file with old code", i);
1299 /* This function stores the last whole configuration of the reference and
1300 * average structure in the .cpt file
1302 static int do_cpt_EDstate(XDR *xd, gmx_bool bRead,
1303 edsamstate_t *EDstate, FILE *list)
1310 EDstate->bFromCpt = bRead;
1312 if (EDstate->nED <= 0)
1317 /* When reading, init_edsam has not been called yet,
1318 * so we have to allocate memory first. */
1321 snew(EDstate->nref, EDstate->nED);
1322 snew(EDstate->old_sref, EDstate->nED);
1323 snew(EDstate->nav, EDstate->nED);
1324 snew(EDstate->old_sav, EDstate->nED);
1327 /* Read/write the last whole conformation of SREF and SAV for each ED dataset (usually only one) */
1328 for (i = 0; i < EDstate->nED; i++)
1330 /* Reference structure SREF */
1331 sprintf(buf, "ED%d # of atoms in reference structure", i+1);
1332 do_cpt_int_err(xd, buf, &EDstate->nref[i], list);
1333 sprintf(buf, "ED%d x_ref", i+1);
1336 snew(EDstate->old_sref[i], EDstate->nref[i]);
1337 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref[i], list);
1341 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref_p[i], list);
1344 /* Average structure SAV */
1345 sprintf(buf, "ED%d # of atoms in average structure", i+1);
1346 do_cpt_int_err(xd, buf, &EDstate->nav[i], list);
1347 sprintf(buf, "ED%d x_av", i+1);
1350 snew(EDstate->old_sav[i], EDstate->nav[i]);
1351 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav[i], list);
1355 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav_p[i], list);
1363 static int do_cpt_files(XDR *xd, gmx_bool bRead,
1364 gmx_file_position_t **p_outputfiles, int *nfiles,
1365 FILE *list, int file_version)
1369 gmx_off_t mask = 0xFFFFFFFFL;
1370 int offset_high, offset_low;
1372 gmx_file_position_t *outputfiles;
1374 if (do_cpt_int(xd, "number of output files", nfiles, list) != 0)
1381 snew(*p_outputfiles, *nfiles);
1384 outputfiles = *p_outputfiles;
1386 for (i = 0; i < *nfiles; i++)
1388 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1391 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1392 strncpy(outputfiles[i].filename, buf, CPTSTRLEN-1);
1398 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1402 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1406 outputfiles[i].offset = ( ((gmx_off_t) offset_high) << 32 ) | ( (gmx_off_t) offset_low & mask );
1410 buf = outputfiles[i].filename;
1411 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1413 offset = outputfiles[i].offset;
1421 offset_low = (int) (offset & mask);
1422 offset_high = (int) ((offset >> 32) & mask);
1424 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1428 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1433 if (file_version >= 8)
1435 if (do_cpt_int(xd, "file_checksum_size", &(outputfiles[i].chksum_size),
1440 if (do_cpt_u_chars(xd, "file_checksum", 16, outputfiles[i].chksum, list) != 0)
1447 outputfiles[i].chksum_size = -1;
1454 void write_checkpoint(const char *fn, gmx_bool bNumberAndKeep,
1455 FILE *fplog, t_commrec *cr,
1456 int eIntegrator, int simulation_part,
1457 gmx_bool bExpanded, int elamstats,
1458 gmx_int64_t step, double t, t_state *state)
1468 char *fntemp; /* the temporary checkpoint file name */
1469 char timebuf[STRLEN];
1470 int nppnodes, npmenodes;
1471 char buf[1024], suffix[5+STEPSTRSIZE], sbuf[STEPSTRSIZE];
1472 gmx_file_position_t *outputfiles;
1475 int flags_eks, flags_enh, flags_dfh;
1478 if (DOMAINDECOMP(cr))
1480 nppnodes = cr->dd->nnodes;
1481 npmenodes = cr->npmenodes;
1489 #ifndef GMX_NO_RENAME
1490 /* make the new temporary filename */
1491 snew(fntemp, strlen(fn)+5+STEPSTRSIZE);
1493 fntemp[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1494 sprintf(suffix, "_%s%s", "step", gmx_step_str(step, sbuf));
1495 strcat(fntemp, suffix);
1496 strcat(fntemp, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1498 /* if we can't rename, we just overwrite the cpt file.
1499 * dangerous if interrupted.
1501 snew(fntemp, strlen(fn));
1504 gmx_format_current_time(timebuf, STRLEN);
1508 fprintf(fplog, "Writing checkpoint, step %s at %s\n\n",
1509 gmx_step_str(step, buf), timebuf);
1512 /* Get offsets for open files */
1513 gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1515 fp = gmx_fio_open(fntemp, "w");
1517 if (state->ekinstate.bUpToDate)
1520 ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) |
1521 (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) |
1522 (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1530 if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1532 flags_enh |= (1<<eenhENERGY_N);
1533 if (state->enerhist.nsum > 0)
1535 flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1536 (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSUM));
1538 if (state->enerhist.nsum_sim > 0)
1540 flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSTEPS_SIM) |
1541 (1<<eenhENERGY_NSUM_SIM));
1543 if (state->enerhist.dht)
1545 flags_enh |= ( (1<< eenhENERGY_DELTA_H_NN) |
1546 (1<< eenhENERGY_DELTA_H_LIST) |
1547 (1<< eenhENERGY_DELTA_H_STARTTIME) |
1548 (1<< eenhENERGY_DELTA_H_STARTLAMBDA) );
1554 flags_dfh = ((1<<edfhBEQUIL) | (1<<edfhNATLAMBDA) | (1<<edfhSUMWEIGHTS) | (1<<edfhSUMDG) |
1555 (1<<edfhTIJ) | (1<<edfhTIJEMP));
1558 flags_dfh |= ((1<<edfhWLDELTA) | (1<<edfhWLHISTO));
1560 if ((elamstats == elamstatsMINVAR) || (elamstats == elamstatsBARKER) || (elamstats == elamstatsMETROPOLIS))
1562 flags_dfh |= ((1<<edfhACCUMP) | (1<<edfhACCUMM) | (1<<edfhACCUMP2) | (1<<edfhACCUMM2)
1563 | (1<<edfhSUMMINVAR) | (1<<edfhSUMVAR));
1571 /* We can check many more things now (CPU, acceleration, etc), but
1572 * it is highly unlikely to have two separate builds with exactly
1573 * the same version, user, time, and build host!
1576 version = gmx_strdup(gmx_version());
1577 btime = gmx_strdup(BUILD_TIME);
1578 buser = gmx_strdup(BUILD_USER);
1579 bhost = gmx_strdup(BUILD_HOST);
1581 double_prec = GMX_CPT_BUILD_DP;
1582 fprog = gmx_strdup(Program());
1584 ftime = &(timebuf[0]);
1586 do_cpt_header(gmx_fio_getxdr(fp), FALSE, &file_version,
1587 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1588 &eIntegrator, &simulation_part, &step, &t, &nppnodes,
1589 DOMAINDECOMP(cr) ? cr->dd->nc : NULL, &npmenodes,
1590 &state->natoms, &state->ngtc, &state->nnhpres,
1591 &state->nhchainlength, &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
1592 &state->edsamstate.nED, &state->swapstate.eSwapCoords,
1601 if ((do_cpt_state(gmx_fio_getxdr(fp), FALSE, state->flags, state, NULL) < 0) ||
1602 (do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL) < 0) ||
1603 (do_cpt_enerhist(gmx_fio_getxdr(fp), FALSE, flags_enh, &state->enerhist, NULL) < 0) ||
1604 (do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL) < 0) ||
1605 (do_cpt_EDstate(gmx_fio_getxdr(fp), FALSE, &state->edsamstate, NULL) < 0) ||
1606 (do_cpt_swapstate(gmx_fio_getxdr(fp), FALSE, &state->swapstate, NULL) < 0) ||
1607 (do_cpt_files(gmx_fio_getxdr(fp), FALSE, &outputfiles, &noutputfiles, NULL,
1610 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1613 do_cpt_footer(gmx_fio_getxdr(fp), file_version);
1615 /* we really, REALLY, want to make sure to physically write the checkpoint,
1616 and all the files it depends on, out to disk. Because we've
1617 opened the checkpoint with gmx_fio_open(), it's in our list
1619 ret = gmx_fio_all_output_fsync();
1625 "Cannot fsync '%s'; maybe you are out of disk space?",
1626 gmx_fio_getname(ret));
1628 if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV) == NULL)
1638 if (gmx_fio_close(fp) != 0)
1640 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1643 /* we don't move the checkpoint if the user specified they didn't want it,
1644 or if the fsyncs failed */
1645 #ifndef GMX_NO_RENAME
1646 if (!bNumberAndKeep && !ret)
1650 /* Rename the previous checkpoint file */
1652 buf[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1653 strcat(buf, "_prev");
1654 strcat(buf, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1656 /* we copy here so that if something goes wrong between now and
1657 * the rename below, there's always a state.cpt.
1658 * If renames are atomic (such as in POSIX systems),
1659 * this copying should be unneccesary.
1661 gmx_file_copy(fn, buf, FALSE);
1662 /* We don't really care if this fails:
1663 * there's already a new checkpoint.
1666 gmx_file_rename(fn, buf);
1669 if (gmx_file_rename(fntemp, fn) != 0)
1671 gmx_file("Cannot rename checkpoint file; maybe you are out of disk space?");
1674 #endif /* GMX_NO_RENAME */
1680 /*code for alternate checkpointing scheme. moved from top of loop over
1682 fcRequestCheckPoint();
1683 if (fcCheckPointParallel( cr->nodeid, NULL, 0) == 0)
1685 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", step );
1687 #endif /* end GMX_FAHCORE block */
1690 static void print_flag_mismatch(FILE *fplog, int sflags, int fflags)
1694 fprintf(fplog, "\nState entry mismatch between the simulation and the checkpoint file\n");
1695 fprintf(fplog, "Entries which are not present in the checkpoint file will not be updated\n");
1696 fprintf(fplog, " %24s %11s %11s\n", "", "simulation", "checkpoint");
1697 for (i = 0; i < estNR; i++)
1699 if ((sflags & (1<<i)) || (fflags & (1<<i)))
1701 fprintf(fplog, " %24s %11s %11s\n",
1703 (sflags & (1<<i)) ? " present " : "not present",
1704 (fflags & (1<<i)) ? " present " : "not present");
1709 static void check_int(FILE *fplog, const char *type, int p, int f, gmx_bool *mm)
1711 FILE *fp = fplog ? fplog : stderr;
1715 fprintf(fp, " %s mismatch,\n", type);
1716 fprintf(fp, " current program: %d\n", p);
1717 fprintf(fp, " checkpoint file: %d\n", f);
1723 static void check_string(FILE *fplog, const char *type, const char *p,
1724 const char *f, gmx_bool *mm)
1726 FILE *fp = fplog ? fplog : stderr;
1728 if (strcmp(p, f) != 0)
1730 fprintf(fp, " %s mismatch,\n", type);
1731 fprintf(fp, " current program: %s\n", p);
1732 fprintf(fp, " checkpoint file: %s\n", f);
1738 static void check_match(FILE *fplog,
1740 char *btime, char *buser, char *bhost, int double_prec,
1742 t_commrec *cr, int npp_f, int npme_f,
1743 ivec dd_nc, ivec dd_nc_f)
1746 gmx_bool mm = FALSE;
1747 gmx_bool patchlevel_differs = FALSE;
1748 gmx_bool version_differs = FALSE;
1750 check_string(fplog, "Version", gmx_version(), version, &mm);
1751 patchlevel_differs = mm;
1753 if (patchlevel_differs)
1755 /* Gromacs should be able to continue from checkpoints between
1756 * different patch level versions, but we do not guarantee
1757 * compatibility between different major/minor versions - check this.
1759 int gmx_major, gmx_minor;
1760 int cpt_major, cpt_minor;
1761 sscanf(gmx_version(), "VERSION %5d.%5d", &gmx_major, &gmx_minor);
1762 int ret = sscanf(version, "VERSION %5d.%5d", &cpt_major, &cpt_minor);
1763 version_differs = (ret < 2 || gmx_major != cpt_major ||
1764 gmx_minor != cpt_minor);
1767 check_string(fplog, "Build time", BUILD_TIME, btime, &mm);
1768 check_string(fplog, "Build user", BUILD_USER, buser, &mm);
1769 check_string(fplog, "Build host", BUILD_HOST, bhost, &mm);
1770 check_int (fplog, "Double prec.", GMX_CPT_BUILD_DP, double_prec, &mm);
1771 check_string(fplog, "Program name", Program(), fprog, &mm);
1773 check_int (fplog, "#ranks", cr->nnodes, npp_f+npme_f, &mm);
1776 check_int (fplog, "#PME-ranks", cr->npmenodes, npme_f, &mm);
1779 if (cr->npmenodes >= 0)
1781 npp -= cr->npmenodes;
1785 check_int (fplog, "#DD-cells[x]", dd_nc[XX], dd_nc_f[XX], &mm);
1786 check_int (fplog, "#DD-cells[y]", dd_nc[YY], dd_nc_f[YY], &mm);
1787 check_int (fplog, "#DD-cells[z]", dd_nc[ZZ], dd_nc_f[ZZ], &mm);
1793 const char msg_version_difference[] =
1794 "The current GROMACS major & minor version are not identical to those that\n"
1795 "generated the checkpoint file. In principle GROMACS does not support\n"
1796 "continuation from checkpoints between different versions, so we advise\n"
1797 "against this. If you still want to try your luck we recommend that you use\n"
1798 "the -noappend flag to keep your output files from the two versions separate.\n"
1799 "This might also work around errors where the output fields in the energy\n"
1800 "file have changed between the different major & minor versions.\n";
1802 const char msg_mismatch_notice[] =
1803 "GROMACS patchlevel, binary or parallel settings differ from previous run.\n"
1804 "Continuation is exact, but not guaranteed to be binary identical.\n";
1806 const char msg_logdetails[] =
1807 "See the log file for details.\n";
1809 if (version_differs)
1811 fprintf(stderr, "%s%s\n", msg_version_difference, fplog ? msg_logdetails : "");
1815 fprintf(fplog, "%s\n", msg_version_difference);
1820 /* Major & minor versions match at least, but something is different. */
1821 fprintf(stderr, "%s%s\n", msg_mismatch_notice, fplog ? msg_logdetails : "");
1824 fprintf(fplog, "%s\n", msg_mismatch_notice);
1830 static void read_checkpoint(const char *fn, FILE **pfplog,
1831 t_commrec *cr, ivec dd_nc,
1832 int eIntegrator, int *init_fep_state, gmx_int64_t *step, double *t,
1833 t_state *state, gmx_bool *bReadEkin,
1834 int *simulation_part,
1835 gmx_bool bAppendOutputFiles, gmx_bool bForceAppend)
1840 char *version, *btime, *buser, *bhost, *fprog, *ftime;
1842 char buf[STEPSTRSIZE];
1843 int eIntegrator_f, nppnodes_f, npmenodes_f;
1845 int natoms, ngtc, nnhpres, nhchainlength, nlambda, fflags, flags_eks, flags_enh, flags_dfh;
1848 gmx_file_position_t *outputfiles;
1850 t_fileio *chksum_file;
1851 FILE * fplog = *pfplog;
1852 unsigned char digest[16];
1853 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1854 struct flock fl; /* don't initialize here: the struct order is OS
1858 const char *int_warn =
1859 "WARNING: The checkpoint file was generated with integrator %s,\n"
1860 " while the simulation uses integrator %s\n\n";
1862 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1863 fl.l_type = F_WRLCK;
1864 fl.l_whence = SEEK_SET;
1870 fp = gmx_fio_open(fn, "r");
1871 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
1872 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1873 &eIntegrator_f, simulation_part, step, t,
1874 &nppnodes_f, dd_nc_f, &npmenodes_f,
1875 &natoms, &ngtc, &nnhpres, &nhchainlength, &nlambda,
1876 &fflags, &flags_eks, &flags_enh, &flags_dfh,
1877 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
1879 if (bAppendOutputFiles &&
1880 file_version >= 13 && double_prec != GMX_CPT_BUILD_DP)
1882 gmx_fatal(FARGS, "Output file appending requested, but the code and checkpoint file precision (single/double) don't match");
1885 if (cr == NULL || MASTER(cr))
1887 fprintf(stderr, "\nReading checkpoint file %s generated: %s\n\n",
1891 /* This will not be written if we do appending, since fplog is still NULL then */
1894 fprintf(fplog, "\n");
1895 fprintf(fplog, "Reading checkpoint file %s\n", fn);
1896 fprintf(fplog, " file generated by: %s\n", fprog);
1897 fprintf(fplog, " file generated at: %s\n", ftime);
1898 fprintf(fplog, " GROMACS build time: %s\n", btime);
1899 fprintf(fplog, " GROMACS build user: %s\n", buser);
1900 fprintf(fplog, " GROMACS build host: %s\n", bhost);
1901 fprintf(fplog, " GROMACS double prec.: %d\n", double_prec);
1902 fprintf(fplog, " simulation part #: %d\n", *simulation_part);
1903 fprintf(fplog, " step: %s\n", gmx_step_str(*step, buf));
1904 fprintf(fplog, " time: %f\n", *t);
1905 fprintf(fplog, "\n");
1908 if (natoms != state->natoms)
1910 gmx_fatal(FARGS, "Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms", natoms, state->natoms);
1912 if (ngtc != state->ngtc)
1914 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);
1916 if (nnhpres != state->nnhpres)
1918 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);
1921 if (nlambda != state->dfhist.nlambda)
1923 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);
1926 init_gtc_state(state, state->ngtc, state->nnhpres, nhchainlength); /* need to keep this here to keep the tpr format working */
1927 /* write over whatever was read; we use the number of Nose-Hoover chains from the checkpoint */
1929 if (eIntegrator_f != eIntegrator)
1933 fprintf(stderr, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1935 if (bAppendOutputFiles)
1938 "Output file appending requested, but input/checkpoint integrators do not match.\n"
1939 "Stopping the run to prevent you from ruining all your data...\n"
1940 "If you _really_ know what you are doing, try with the -noappend option.\n");
1944 fprintf(fplog, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1952 else if (cr->nnodes == nppnodes_f + npmenodes_f)
1954 if (cr->npmenodes < 0)
1956 cr->npmenodes = npmenodes_f;
1958 int nppnodes = cr->nnodes - cr->npmenodes;
1959 if (nppnodes == nppnodes_f)
1961 for (d = 0; d < DIM; d++)
1965 dd_nc[d] = dd_nc_f[d];
1971 if (fflags != state->flags)
1976 if (bAppendOutputFiles)
1979 "Output file appending requested, but input and checkpoint states are not identical.\n"
1980 "Stopping the run to prevent you from ruining all your data...\n"
1981 "You can try with the -noappend option, and get more info in the log file.\n");
1984 if (getenv("GMX_ALLOW_CPT_MISMATCH") == NULL)
1986 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");
1991 "WARNING: The checkpoint state entries do not match the simulation,\n"
1992 " see the log file for details\n\n");
1998 print_flag_mismatch(fplog, state->flags, fflags);
2005 check_match(fplog, version, btime, buser, bhost, double_prec, fprog,
2006 cr, nppnodes_f, npmenodes_f, dd_nc, dd_nc_f);
2009 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, fflags, state, NULL);
2010 *init_fep_state = state->fep_state; /* there should be a better way to do this than setting it here.
2011 Investigate for 5.0. */
2016 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2021 *bReadEkin = ((flags_eks & (1<<eeksEKINH)) || (flags_eks & (1<<eeksEKINF)) || (flags_eks & (1<<eeksEKINO)) ||
2022 ((flags_eks & (1<<eeksEKINSCALEF)) | (flags_eks & (1<<eeksEKINSCALEH)) | (flags_eks & (1<<eeksVSCALE))));
2024 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2025 flags_enh, &state->enerhist, NULL);
2031 if (file_version < 6)
2033 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.";
2035 fprintf(stderr, "\nWARNING: %s\n\n", warn);
2038 fprintf(fplog, "\nWARNING: %s\n\n", warn);
2040 state->enerhist.nsum = *step;
2041 state->enerhist.nsum_sim = *step;
2044 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2050 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2056 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2062 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, NULL, file_version);
2068 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2073 if (gmx_fio_close(fp) != 0)
2075 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2084 /* If the user wants to append to output files,
2085 * we use the file pointer positions of the output files stored
2086 * in the checkpoint file and truncate the files such that any frames
2087 * written after the checkpoint time are removed.
2088 * All files are md5sum checked such that we can be sure that
2089 * we do not truncate other (maybe imprortant) files.
2091 if (bAppendOutputFiles)
2093 if (fn2ftp(outputfiles[0].filename) != efLOG)
2095 /* make sure first file is log file so that it is OK to use it for
2098 gmx_fatal(FARGS, "The first output file should always be the log "
2099 "file but instead is: %s. Cannot do appending because of this condition.", outputfiles[0].filename);
2101 for (i = 0; i < nfiles; i++)
2103 if (outputfiles[i].offset < 0)
2105 gmx_fatal(FARGS, "The original run wrote a file called '%s' which "
2106 "is larger than 2 GB, but mdrun did not support large file"
2107 " offsets. Can not append. Run mdrun with -noappend",
2108 outputfiles[i].filename);
2111 chksum_file = gmx_fio_open(outputfiles[i].filename, "a");
2114 chksum_file = gmx_fio_open(outputfiles[i].filename, "r+");
2119 /* Note that there are systems where the lock operation
2120 * will succeed, but a second process can also lock the file.
2121 * We should probably try to detect this.
2123 #if defined __native_client__
2127 #elif defined GMX_NATIVE_WINDOWS
2128 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX) == -1)
2130 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl) == -1)
2133 if (errno == ENOSYS)
2137 gmx_fatal(FARGS, "File locking is not supported on this system. Use -noappend or specify -append explicitly to append anyhow.");
2141 fprintf(stderr, "\nNOTE: File locking is not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2144 fprintf(fplog, "\nNOTE: File locking not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2148 else if (errno == EACCES || errno == EAGAIN)
2150 gmx_fatal(FARGS, "Failed to lock: %s. Already running "
2151 "simulation?", outputfiles[i].filename);
2155 gmx_fatal(FARGS, "Failed to lock: %s. %s.",
2156 outputfiles[i].filename, strerror(errno));
2161 /* compute md5 chksum */
2162 if (outputfiles[i].chksum_size != -1)
2164 if (gmx_fio_get_file_md5(chksum_file, outputfiles[i].offset,
2165 digest) != outputfiles[i].chksum_size) /*at the end of the call the file position is at the end of the file*/
2167 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.",
2168 outputfiles[i].chksum_size,
2169 outputfiles[i].filename);
2172 if (i == 0) /*log file needs to be seeked in case we need to truncate (other files are truncated below)*/
2174 if (gmx_fio_seek(chksum_file, outputfiles[i].offset))
2176 gmx_fatal(FARGS, "Seek error! Failed to truncate log-file: %s.", strerror(errno));
2181 if (i == 0) /*open log file here - so that lock is never lifted
2182 after chksum is calculated */
2184 *pfplog = gmx_fio_getfp(chksum_file);
2188 gmx_fio_close(chksum_file);
2191 /* compare md5 chksum */
2192 if (outputfiles[i].chksum_size != -1 &&
2193 memcmp(digest, outputfiles[i].chksum, 16) != 0)
2197 fprintf(debug, "chksum for %s: ", outputfiles[i].filename);
2198 for (j = 0; j < 16; j++)
2200 fprintf(debug, "%02x", digest[j]);
2202 fprintf(debug, "\n");
2204 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.",
2205 outputfiles[i].filename);
2210 if (i != 0) /*log file is already seeked to correct position */
2212 #if !defined(GMX_NATIVE_WINDOWS) || !defined(GMX_FAHCORE)
2213 /* For FAHCORE, we do this elsewhere*/
2214 rc = gmx_truncate(outputfiles[i].filename, outputfiles[i].offset);
2217 gmx_fatal(FARGS, "Truncation of file %s failed. Cannot do appending because of this failure.", outputfiles[i].filename);
2228 void load_checkpoint(const char *fn, FILE **fplog,
2229 t_commrec *cr, ivec dd_nc,
2230 t_inputrec *ir, t_state *state,
2231 gmx_bool *bReadEkin,
2232 gmx_bool bAppend, gmx_bool bForceAppend)
2239 /* Read the state from the checkpoint file */
2240 read_checkpoint(fn, fplog,
2242 ir->eI, &(ir->fepvals->init_fep_state), &step, &t, state, bReadEkin,
2243 &ir->simulation_part, bAppend, bForceAppend);
2247 gmx_bcast(sizeof(cr->npmenodes), &cr->npmenodes, cr);
2248 gmx_bcast(DIM*sizeof(dd_nc[0]), dd_nc, cr);
2249 gmx_bcast(sizeof(step), &step, cr);
2250 gmx_bcast(sizeof(*bReadEkin), bReadEkin, cr);
2252 ir->bContinuation = TRUE;
2253 if (ir->nsteps >= 0)
2255 ir->nsteps += ir->init_step - step;
2257 ir->init_step = step;
2258 ir->simulation_part += 1;
2261 void read_checkpoint_part_and_step(const char *filename,
2262 int *simulation_part,
2266 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2271 int flags_eks, flags_enh, flags_dfh;
2276 if (filename == NULL ||
2277 !gmx_fexist(filename) ||
2278 (!(fp = gmx_fio_open(filename, "r"))))
2280 *simulation_part = 0;
2285 /* Not calling initializing state before use is nasty, but all we
2286 do is read into its member variables and throw the struct away
2287 again immediately. */
2289 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2290 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2291 &eIntegrator, simulation_part, step, &t, &nppnodes, dd_nc, &npme,
2292 &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
2293 &(state.dfhist.nlambda), &state.flags, &flags_eks, &flags_enh, &flags_dfh,
2294 &state.edsamstate.nED, &state.swapstate.eSwapCoords, NULL);
2299 static void read_checkpoint_data(t_fileio *fp, int *simulation_part,
2300 gmx_int64_t *step, double *t, t_state *state,
2301 int *nfiles, gmx_file_position_t **outputfiles)
2304 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2309 int flags_eks, flags_enh, flags_dfh;
2311 gmx_file_position_t *files_loc = NULL;
2314 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2315 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2316 &eIntegrator, simulation_part, step, t, &nppnodes, dd_nc, &npme,
2317 &state->natoms, &state->ngtc, &state->nnhpres, &state->nhchainlength,
2318 &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
2319 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
2321 do_cpt_state(gmx_fio_getxdr(fp), TRUE, state->flags, state, NULL);
2326 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2331 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2332 flags_enh, &state->enerhist, NULL);
2337 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2343 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2349 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2355 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE,
2356 outputfiles != NULL ? outputfiles : &files_loc,
2357 outputfiles != NULL ? nfiles : &nfiles_loc,
2358 NULL, file_version);
2359 if (files_loc != NULL)
2369 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2383 read_checkpoint_state(const char *fn, int *simulation_part,
2384 gmx_int64_t *step, double *t, t_state *state)
2388 fp = gmx_fio_open(fn, "r");
2389 read_checkpoint_data(fp, simulation_part, step, t, state, NULL, NULL);
2390 if (gmx_fio_close(fp) != 0)
2392 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2396 void read_checkpoint_trxframe(t_fileio *fp, t_trxframe *fr)
2398 /* This next line is nasty because the sub-structures of t_state
2399 * cannot be assumed to be zeroed (or even initialized in ways the
2400 * rest of the code might assume). Using snew would be better, but
2401 * this will all go away for 5.0. */
2403 int simulation_part;
2407 init_state(&state, 0, 0, 0, 0, 0);
2409 read_checkpoint_data(fp, &simulation_part, &step, &t, &state, NULL, NULL);
2411 fr->natoms = state.natoms;
2414 fr->step = gmx_int64_to_int(step,
2415 "conversion of checkpoint to trajectory");
2419 fr->lambda = state.lambda[efptFEP];
2420 fr->fep_state = state.fep_state;
2422 fr->bX = (state.flags & (1<<estX));
2428 fr->bV = (state.flags & (1<<estV));
2435 fr->bBox = (state.flags & (1<<estBOX));
2438 copy_mat(state.box, fr->box);
2443 void list_checkpoint(const char *fn, FILE *out)
2447 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2449 int eIntegrator, simulation_part, nppnodes, npme;
2454 int flags_eks, flags_enh, flags_dfh;
2456 gmx_file_position_t *outputfiles;
2459 init_state(&state, -1, -1, -1, -1, 0);
2461 fp = gmx_fio_open(fn, "r");
2462 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2463 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2464 &eIntegrator, &simulation_part, &step, &t, &nppnodes, dd_nc, &npme,
2465 &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
2466 &(state.dfhist.nlambda), &state.flags,
2467 &flags_eks, &flags_enh, &flags_dfh, &state.edsamstate.nED,
2468 &state.swapstate.eSwapCoords, out);
2469 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, state.flags, &state, out);
2474 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state.ekinstate, out);
2479 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2480 flags_enh, &state.enerhist, out);
2484 ret = do_cpt_df_hist(gmx_fio_getxdr(fp),
2485 flags_dfh, &state.dfhist, out);
2490 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state.edsamstate, out);
2495 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state.swapstate, out);
2500 do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, out, file_version);
2505 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2512 if (gmx_fio_close(fp) != 0)
2514 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2520 /* This routine cannot print tons of data, since it is called before the log file is opened. */
2522 read_checkpoint_simulation_part_and_filenames(t_fileio *fp,
2523 int *simulation_part,
2525 gmx_file_position_t **outputfiles)
2527 gmx_int64_t step = 0;
2531 init_state(&state, 0, 0, 0, 0, 0);
2533 read_checkpoint_data(fp, simulation_part, &step, &t, &state,
2534 nfiles, outputfiles);
2535 if (gmx_fio_close(fp) != 0)
2537 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");