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/gmxfio-xdr.h"
58 #include "gromacs/fileio/xdr_datatype.h"
59 #include "gromacs/fileio/xdrf.h"
60 #include "gromacs/legacyheaders/copyrite.h"
61 #include "gromacs/legacyheaders/names.h"
62 #include "gromacs/legacyheaders/network.h"
63 #include "gromacs/legacyheaders/txtdump.h"
64 #include "gromacs/legacyheaders/typedefs.h"
65 #include "gromacs/legacyheaders/types/commrec.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/utility/baseversion.h"
68 #include "gromacs/utility/cstringutil.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/futil.h"
71 #include "gromacs/utility/smalloc.h"
72 #include "gromacs/utility/sysinfo.h"
78 #define CPT_MAGIC1 171817
79 #define CPT_MAGIC2 171819
80 #define CPTSTRLEN 1024
83 #define GMX_CPT_BUILD_DP 1
85 #define GMX_CPT_BUILD_DP 0
88 /* cpt_version should normally only be changed
89 * when the header of footer format changes.
90 * The state data format itself is backward and forward compatible.
91 * But old code can not read a new entry that is present in the file
92 * (but can read a new format when new entries are not present).
94 static const int cpt_version = 16;
97 const char *est_names[estNR] =
100 "box", "box-rel", "box-v", "pres_prev",
101 "nosehoover-xi", "thermostat-integral",
102 "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
103 "disre_initf", "disre_rm3tav",
104 "orire_initf", "orire_Dtav",
105 "svir_prev", "nosehoover-vxi", "v_eta", "vol0", "nhpres_xi", "nhpres_vxi", "fvir_prev", "fep_state", "MC-rng", "MC-rng-i"
109 eeksEKIN_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksEKINF, eeksEKINO, eeksEKINSCALEF, eeksEKINSCALEH, eeksVSCALE, eeksEKINTOTAL, eeksNR
112 const char *eeks_names[eeksNR] =
114 "Ekin_n", "Ekinh", "dEkindlambda", "mv_cos",
115 "Ekinf", "Ekinh_old", "EkinScaleF_NHC", "EkinScaleH_NHC", "Vscale_NHC", "Ekin_Total"
119 eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
120 eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
121 eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM,
122 eenhENERGY_DELTA_H_NN,
123 eenhENERGY_DELTA_H_LIST,
124 eenhENERGY_DELTA_H_STARTTIME,
125 eenhENERGY_DELTA_H_STARTLAMBDA,
129 const char *eenh_names[eenhNR] =
131 "energy_n", "energy_aver", "energy_sum", "energy_nsum",
132 "energy_sum_sim", "energy_nsum_sim",
133 "energy_nsteps", "energy_nsteps_sim",
135 "energy_delta_h_list",
136 "energy_delta_h_start_time",
137 "energy_delta_h_start_lambda"
140 /* free energy history variables -- need to be preserved over checkpoint */
142 edfhBEQUIL, edfhNATLAMBDA, edfhWLHISTO, edfhWLDELTA, edfhSUMWEIGHTS, edfhSUMDG, edfhSUMMINVAR, edfhSUMVAR,
143 edfhACCUMP, edfhACCUMM, edfhACCUMP2, edfhACCUMM2, edfhTIJ, edfhTIJEMP, edfhNR
145 /* free energy history variable names */
146 const char *edfh_names[edfhNR] =
148 "bEquilibrated", "N_at_state", "Wang-Landau Histogram", "Wang-Landau Delta", "Weights", "Free Energies", "minvar", "variance",
149 "accumulated_plus", "accumulated_minus", "accumulated_plus_2", "accumulated_minus_2", "Tij", "Tij_empirical"
153 ecprREAL, ecprRVEC, ecprMATRIX
157 cptpEST, cptpEEKS, cptpEENH, cptpEDFH
159 /* enums for the different components of checkpoint variables, replacing the hard coded ones.
160 cptpEST - state variables.
161 cptpEEKS - Kinetic energy state variables.
162 cptpEENH - Energy history state variables.
163 cptpEDFH - free energy history variables.
167 static const char *st_names(int cptp, int ecpt)
171 case cptpEST: return est_names [ecpt]; break;
172 case cptpEEKS: return eeks_names[ecpt]; break;
173 case cptpEENH: return eenh_names[ecpt]; break;
174 case cptpEDFH: return edfh_names[ecpt]; break;
180 static void cp_warning(FILE *fp)
182 fprintf(fp, "\nWARNING: Checkpoint file is corrupted or truncated\n\n");
185 static void cp_error()
187 gmx_fatal(FARGS, "Checkpoint file corrupted/truncated, or maybe you are out of disk space?");
190 static void do_cpt_string_err(XDR *xd, gmx_bool bRead, const char *desc, char **s, FILE *list)
198 res = xdr_string(xd, s, CPTSTRLEN);
205 fprintf(list, "%s = %s\n", desc, *s);
210 static int do_cpt_int(XDR *xd, const char *desc, int *i, FILE *list)
214 res = xdr_int(xd, i);
221 fprintf(list, "%s = %d\n", desc, *i);
226 static int do_cpt_u_chars(XDR *xd, const char *desc, int n, unsigned char *i, FILE *list)
232 fprintf(list, "%s = ", desc);
234 for (j = 0; j < n && res; j++)
236 res &= xdr_u_char(xd, &i[j]);
239 fprintf(list, "%02x", i[j]);
254 static void do_cpt_int_err(XDR *xd, const char *desc, int *i, FILE *list)
256 if (do_cpt_int(xd, desc, i, list) < 0)
262 static void do_cpt_step_err(XDR *xd, const char *desc, gmx_int64_t *i, FILE *list)
265 char buf[STEPSTRSIZE];
267 res = xdr_int64(xd, i);
274 fprintf(list, "%s = %s\n", desc, gmx_step_str(*i, buf));
278 static void do_cpt_double_err(XDR *xd, const char *desc, double *f, FILE *list)
282 res = xdr_double(xd, f);
289 fprintf(list, "%s = %f\n", desc, *f);
293 static void do_cpt_real_err(XDR *xd, real *f)
298 res = xdr_double(xd, f);
300 res = xdr_float(xd, f);
308 static void do_cpt_n_rvecs_err(XDR *xd, const char *desc, int n, rvec f[], FILE *list)
312 for (i = 0; i < n; i++)
314 for (j = 0; j < DIM; j++)
316 do_cpt_real_err(xd, &f[i][j]);
322 pr_rvecs(list, 0, desc, f, n);
326 /* If nval >= 0, nval is used; on read this should match the passed value.
327 * If nval n<0, *nptr is used; on read the value is stored in nptr
329 static int do_cpte_reals_low(XDR *xd, int cptp, int ecpt, int sflags,
330 int nval, int *nptr, real **v,
331 FILE *list, int erealtype)
335 int dtc = xdr_datatype_float;
337 int dtc = xdr_datatype_double;
339 real *vp, *va = NULL;
354 gmx_incons("*ntpr=NULL in do_cpte_reals_low");
359 res = xdr_int(xd, &nf);
370 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), nval, nf);
379 res = xdr_int(xd, &dt);
386 fprintf(stderr, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
387 st_names(cptp, ecpt), xdr_datatype_names[dtc],
388 xdr_datatype_names[dt]);
390 if (list || !(sflags & (1<<ecpt)))
403 if (dt == xdr_datatype_float)
405 if (dtc == xdr_datatype_float)
413 res = xdr_vector(xd, (char *)vf, nf,
414 (unsigned int)sizeof(float), (xdrproc_t)xdr_float);
419 if (dtc != xdr_datatype_float)
421 for (i = 0; i < nf; i++)
430 if (dtc == xdr_datatype_double)
432 /* cppcheck-suppress invalidPointerCast
433 * Only executed if real is anyhow double */
440 res = xdr_vector(xd, (char *)vd, nf,
441 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
446 if (dtc != xdr_datatype_double)
448 for (i = 0; i < nf; i++)
461 pr_reals(list, 0, st_names(cptp, ecpt), vp, nf);
464 pr_rvecs(list, 0, st_names(cptp, ecpt), (rvec *)vp, nf/3);
467 gmx_incons("Unknown checkpoint real type");
479 /* This function stores n along with the reals for reading,
480 * but on reading it assumes that n matches the value in the checkpoint file,
481 * a fatal error is generated when this is not the case.
483 static int do_cpte_reals(XDR *xd, int cptp, int ecpt, int sflags,
484 int n, real **v, FILE *list)
486 return do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, v, list, ecprREAL);
489 /* This function does the same as do_cpte_reals,
490 * except that on reading it ignores the passed value of *n
491 * and stored the value read from the checkpoint file in *n.
493 static int do_cpte_n_reals(XDR *xd, int cptp, int ecpt, int sflags,
494 int *n, real **v, FILE *list)
496 return do_cpte_reals_low(xd, cptp, ecpt, sflags, -1, n, v, list, ecprREAL);
499 static int do_cpte_real(XDR *xd, int cptp, int ecpt, int sflags,
502 return do_cpte_reals_low(xd, cptp, ecpt, sflags, 1, NULL, &r, list, ecprREAL);
505 static int do_cpte_ints(XDR *xd, int cptp, int ecpt, int sflags,
506 int n, int **v, FILE *list)
509 int dtc = xdr_datatype_int;
514 res = xdr_int(xd, &nf);
519 if (list == NULL && v != NULL && nf != n)
521 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
524 res = xdr_int(xd, &dt);
531 gmx_fatal(FARGS, "Type mismatch for state entry %s, code type is %s, file type is %s\n",
532 st_names(cptp, ecpt), xdr_datatype_names[dtc],
533 xdr_datatype_names[dt]);
535 if (list || !(sflags & (1<<ecpt)) || v == NULL)
548 res = xdr_vector(xd, (char *)vp, nf,
549 (unsigned int)sizeof(int), (xdrproc_t)xdr_int);
556 pr_ivec(list, 0, st_names(cptp, ecpt), vp, nf, TRUE);
566 static int do_cpte_int(XDR *xd, int cptp, int ecpt, int sflags,
569 return do_cpte_ints(xd, cptp, ecpt, sflags, 1, &i, list);
572 static int do_cpte_doubles(XDR *xd, int cptp, int ecpt, int sflags,
573 int n, double **v, FILE *list)
576 int dtc = xdr_datatype_double;
577 double *vp, *va = NULL;
581 res = xdr_int(xd, &nf);
586 if (list == NULL && nf != n)
588 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
591 res = xdr_int(xd, &dt);
598 gmx_fatal(FARGS, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
599 st_names(cptp, ecpt), xdr_datatype_names[dtc],
600 xdr_datatype_names[dt]);
602 if (list || !(sflags & (1<<ecpt)))
615 res = xdr_vector(xd, (char *)vp, nf,
616 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
623 pr_doubles(list, 0, st_names(cptp, ecpt), vp, nf);
633 static int do_cpte_double(XDR *xd, int cptp, int ecpt, int sflags,
634 double *r, FILE *list)
636 return do_cpte_doubles(xd, cptp, ecpt, sflags, 1, &r, list);
640 static int do_cpte_rvecs(XDR *xd, int cptp, int ecpt, int sflags,
641 int n, rvec **v, FILE *list)
643 return do_cpte_reals_low(xd, cptp, ecpt, sflags,
644 n*DIM, NULL, (real **)v, list, ecprRVEC);
647 static int do_cpte_matrix(XDR *xd, int cptp, int ecpt, int sflags,
648 matrix v, FILE *list)
653 vr = (real *)&(v[0][0]);
654 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
655 DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
657 if (list && ret == 0)
659 pr_rvecs(list, 0, st_names(cptp, ecpt), v, DIM);
666 static int do_cpte_nmatrix(XDR *xd, int cptp, int ecpt, int sflags,
667 int n, real **v, FILE *list)
671 char name[CPTSTRLEN];
678 for (i = 0; i < n; i++)
680 reti = do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, &(v[i]), NULL, ecprREAL);
681 if (list && reti == 0)
683 sprintf(name, "%s[%d]", st_names(cptp, ecpt), i);
684 pr_reals(list, 0, name, v[i], n);
694 static int do_cpte_matrices(XDR *xd, int cptp, int ecpt, int sflags,
695 int n, matrix **v, FILE *list)
698 matrix *vp, *va = NULL;
704 res = xdr_int(xd, &nf);
709 if (list == NULL && nf != n)
711 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
713 if (list || !(sflags & (1<<ecpt)))
726 snew(vr, nf*DIM*DIM);
727 for (i = 0; i < nf; i++)
729 for (j = 0; j < DIM; j++)
731 for (k = 0; k < DIM; k++)
733 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
737 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
738 nf*DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
739 for (i = 0; i < nf; i++)
741 for (j = 0; j < DIM; j++)
743 for (k = 0; k < DIM; k++)
745 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
751 if (list && ret == 0)
753 for (i = 0; i < nf; i++)
755 pr_rvecs(list, 0, st_names(cptp, ecpt), vp[i], DIM);
766 static void do_cpt_header(XDR *xd, gmx_bool bRead, int *file_version,
767 char **version, char **btime, char **buser, char **bhost,
769 char **fprog, char **ftime,
770 int *eIntegrator, int *simulation_part,
771 gmx_int64_t *step, double *t,
772 int *nnodes, int *dd_nc, int *npme,
773 int *natoms, int *ngtc, int *nnhpres, int *nhchainlength,
774 int *nlambda, int *flags_state,
775 int *flags_eks, int *flags_enh, int *flags_dfh,
776 int *nED, int *eSwapCoords,
792 res = xdr_int(xd, &magic);
795 gmx_fatal(FARGS, "The checkpoint file is empty/corrupted, or maybe you are out of disk space?");
797 if (magic != CPT_MAGIC1)
799 gmx_fatal(FARGS, "Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
800 "The checkpoint file is corrupted or not a checkpoint file",
806 gmx_gethostname(fhost, 255);
808 do_cpt_string_err(xd, bRead, "GROMACS version", version, list);
809 do_cpt_string_err(xd, bRead, "GROMACS build time", btime, list);
810 do_cpt_string_err(xd, bRead, "GROMACS build user", buser, list);
811 do_cpt_string_err(xd, bRead, "GROMACS build host", bhost, list);
812 do_cpt_string_err(xd, bRead, "generating program", fprog, list);
813 do_cpt_string_err(xd, bRead, "generation time", ftime, list);
814 *file_version = cpt_version;
815 do_cpt_int_err(xd, "checkpoint file version", file_version, list);
816 if (*file_version > cpt_version)
818 gmx_fatal(FARGS, "Attempting to read a checkpoint file of version %d with code of version %d\n", *file_version, cpt_version);
820 if (*file_version >= 13)
822 do_cpt_int_err(xd, "GROMACS double precision", double_prec, list);
828 if (*file_version >= 12)
830 do_cpt_string_err(xd, bRead, "generating host", &fhost, list);
836 do_cpt_int_err(xd, "#atoms", natoms, list);
837 do_cpt_int_err(xd, "#T-coupling groups", ngtc, list);
838 if (*file_version >= 10)
840 do_cpt_int_err(xd, "#Nose-Hoover T-chains", nhchainlength, list);
846 if (*file_version >= 11)
848 do_cpt_int_err(xd, "#Nose-Hoover T-chains for barostat ", nnhpres, list);
854 if (*file_version >= 14)
856 do_cpt_int_err(xd, "# of total lambda states ", nlambda, list);
862 do_cpt_int_err(xd, "integrator", eIntegrator, list);
863 if (*file_version >= 3)
865 do_cpt_int_err(xd, "simulation part #", simulation_part, list);
869 *simulation_part = 1;
871 if (*file_version >= 5)
873 do_cpt_step_err(xd, "step", step, list);
877 do_cpt_int_err(xd, "step", &idum, list);
880 do_cpt_double_err(xd, "t", t, list);
881 do_cpt_int_err(xd, "#PP-ranks", nnodes, list);
883 do_cpt_int_err(xd, "dd_nc[x]", dd_nc ? &(dd_nc[0]) : &idum, list);
884 do_cpt_int_err(xd, "dd_nc[y]", dd_nc ? &(dd_nc[1]) : &idum, list);
885 do_cpt_int_err(xd, "dd_nc[z]", dd_nc ? &(dd_nc[2]) : &idum, list);
886 do_cpt_int_err(xd, "#PME-only ranks", npme, list);
887 do_cpt_int_err(xd, "state flags", flags_state, list);
888 if (*file_version >= 4)
890 do_cpt_int_err(xd, "ekin data flags", flags_eks, list);
891 do_cpt_int_err(xd, "energy history flags", flags_enh, list);
896 *flags_enh = (*flags_state >> (estORIRE_DTAV+1));
897 *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
898 (1<<(estORIRE_DTAV+2)) |
899 (1<<(estORIRE_DTAV+3))));
901 if (*file_version >= 14)
903 do_cpt_int_err(xd, "df history flags", flags_dfh, list);
910 if (*file_version >= 15)
912 do_cpt_int_err(xd, "ED data sets", nED, list);
918 if (*file_version >= 16)
920 do_cpt_int_err(xd, "swap", eSwapCoords, list);
924 static int do_cpt_footer(XDR *xd, int file_version)
929 if (file_version >= 2)
932 res = xdr_int(xd, &magic);
937 if (magic != CPT_MAGIC2)
946 static int do_cpt_state(XDR *xd, gmx_bool bRead,
947 int fflags, t_state *state,
957 nnht = state->nhchainlength*state->ngtc;
958 nnhtp = state->nhchainlength*state->nnhpres;
960 if (bRead) /* we need to allocate space for dfhist if we are reading */
962 init_df_history(&state->dfhist, state->dfhist.nlambda);
965 sflags = state->flags;
966 for (i = 0; (i < estNR && ret == 0); i++)
972 case estLAMBDA: ret = do_cpte_reals(xd, cptpEST, i, sflags, efptNR, &(state->lambda), list); break;
973 case estFEPSTATE: ret = do_cpte_int (xd, cptpEST, i, sflags, &state->fep_state, list); break;
974 case estBOX: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box, list); break;
975 case estBOX_REL: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box_rel, list); break;
976 case estBOXV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->boxv, list); break;
977 case estPRES_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->pres_prev, list); break;
978 case estSVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->svir_prev, list); break;
979 case estFVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->fvir_prev, list); break;
980 case estNH_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_xi, list); break;
981 case estNH_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_vxi, list); break;
982 case estNHPRES_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_xi, list); break;
983 case estNHPRES_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_vxi, list); break;
984 case estTC_INT: ret = do_cpte_doubles(xd, cptpEST, i, sflags, state->ngtc, &state->therm_integral, list); break;
985 case estVETA: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->veta, list); break;
986 case estVOL0: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->vol0, list); break;
987 case estX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->x, list); break;
988 case estV: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->v, list); break;
989 case estSDX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->sd_X, list); break;
990 /* The RNG entries are no longer written,
991 * the next 4 lines are only for reading old files.
993 case estLD_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
994 case estLD_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
995 case estMC_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
996 case estMC_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
997 case estDISRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.disre_initf, list); break;
998 case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.ndisrepairs, &state->hist.disre_rm3tav, list); break;
999 case estORIRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.orire_initf, list); break;
1000 case estORIRE_DTAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.norire_Dtav, &state->hist.orire_Dtav, list); break;
1002 gmx_fatal(FARGS, "Unknown state entry %d\n"
1003 "You are probably reading a new checkpoint file with old code", i);
1011 static int do_cpt_ekinstate(XDR *xd, int fflags, ekinstate_t *ekins,
1019 for (i = 0; (i < eeksNR && ret == 0); i++)
1021 if (fflags & (1<<i))
1026 case eeksEKIN_N: ret = do_cpte_int(xd, cptpEEKS, i, fflags, &ekins->ekin_n, list); break;
1027 case eeksEKINH: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh, list); break;
1028 case eeksEKINF: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinf, list); break;
1029 case eeksEKINO: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh_old, list); break;
1030 case eeksEKINTOTAL: ret = do_cpte_matrix(xd, cptpEEKS, i, fflags, ekins->ekin_total, list); break;
1031 case eeksEKINSCALEF: ret = do_cpte_doubles(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinscalef_nhc, list); break;
1032 case eeksVSCALE: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->vscale_nhc, list); break;
1033 case eeksEKINSCALEH: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->ekinscaleh_nhc, list); break;
1034 case eeksDEKINDL: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->dekindl, list); break;
1035 case eeksMVCOS: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->mvcos, list); break;
1037 gmx_fatal(FARGS, "Unknown ekin data state entry %d\n"
1038 "You are probably reading a new checkpoint file with old code", i);
1047 static int do_cpt_swapstate(XDR *xd, gmx_bool bRead, swapstate_t *swapstate, FILE *list)
1051 int swap_cpt_version = 1;
1054 if (eswapNO == swapstate->eSwapCoords)
1059 swapstate->bFromCpt = bRead;
1061 do_cpt_int_err(xd, "swap checkpoint version", &swap_cpt_version, list);
1062 do_cpt_int_err(xd, "swap coupling steps", &swapstate->nAverage, list);
1064 /* When reading, init_swapcoords has not been called yet,
1065 * so we have to allocate memory first. */
1067 for (ic = 0; ic < eCompNR; ic++)
1069 for (ii = 0; ii < eIonNR; ii++)
1073 do_cpt_int_err(xd, "swap requested atoms", &swapstate->nat_req[ic][ii], list);
1077 do_cpt_int_err(xd, "swap requested atoms p", swapstate->nat_req_p[ic][ii], list);
1082 do_cpt_int_err(xd, "swap influx netto", &swapstate->inflow_netto[ic][ii], list);
1086 do_cpt_int_err(xd, "swap influx netto p", swapstate->inflow_netto_p[ic][ii], list);
1089 if (bRead && (NULL == swapstate->nat_past[ic][ii]) )
1091 snew(swapstate->nat_past[ic][ii], swapstate->nAverage);
1094 for (j = 0; j < swapstate->nAverage; j++)
1098 do_cpt_int_err(xd, "swap past atom counts", &swapstate->nat_past[ic][ii][j], list);
1102 do_cpt_int_err(xd, "swap past atom counts p", &swapstate->nat_past_p[ic][ii][j], list);
1108 /* Ion flux per channel */
1109 for (ic = 0; ic < eChanNR; ic++)
1111 for (ii = 0; ii < eIonNR; ii++)
1115 do_cpt_int_err(xd, "channel flux", &swapstate->fluxfromAtoB[ic][ii], list);
1119 do_cpt_int_err(xd, "channel flux p", swapstate->fluxfromAtoB_p[ic][ii], list);
1124 /* Ion flux leakage */
1127 snew(swapstate->fluxleak, 1);
1129 do_cpt_int_err(xd, "flux leakage", swapstate->fluxleak, list);
1132 do_cpt_int_err(xd, "number of ions", &swapstate->nions, list);
1136 snew(swapstate->channel_label, swapstate->nions);
1137 snew(swapstate->comp_from, swapstate->nions);
1140 do_cpt_u_chars(xd, "channel history", swapstate->nions, swapstate->channel_label, list);
1141 do_cpt_u_chars(xd, "domain history", swapstate->nions, swapstate->comp_from, list);
1143 /* Save the last known whole positions to checkpoint
1144 * file to be able to also make multimeric channels whole in PBC */
1145 do_cpt_int_err(xd, "Ch0 atoms", &swapstate->nat[eChan0], list);
1146 do_cpt_int_err(xd, "Ch1 atoms", &swapstate->nat[eChan1], list);
1149 snew(swapstate->xc_old_whole[eChan0], swapstate->nat[eChan0]);
1150 snew(swapstate->xc_old_whole[eChan1], swapstate->nat[eChan1]);
1151 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], swapstate->xc_old_whole[eChan0], list);
1152 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], swapstate->xc_old_whole[eChan1], list);
1156 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], *swapstate->xc_old_whole_p[eChan0], list);
1157 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], *swapstate->xc_old_whole_p[eChan1], list);
1164 static int do_cpt_enerhist(XDR *xd, gmx_bool bRead,
1165 int fflags, energyhistory_t *enerhist,
1176 enerhist->nsteps = 0;
1178 enerhist->nsteps_sim = 0;
1179 enerhist->nsum_sim = 0;
1180 enerhist->dht = NULL;
1182 if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
1184 snew(enerhist->dht, 1);
1185 enerhist->dht->ndh = NULL;
1186 enerhist->dht->dh = NULL;
1187 enerhist->dht->start_lambda_set = FALSE;
1191 for (i = 0; (i < eenhNR && ret == 0); i++)
1193 if (fflags & (1<<i))
1197 case eenhENERGY_N: ret = do_cpte_int(xd, cptpEENH, i, fflags, &enerhist->nener, list); break;
1198 case eenhENERGY_AVER: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_ave, list); break;
1199 case eenhENERGY_SUM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum, list); break;
1200 case eenhENERGY_NSUM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum, list); break;
1201 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum_sim, list); break;
1202 case eenhENERGY_NSUM_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum_sim, list); break;
1203 case eenhENERGY_NSTEPS: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps, list); break;
1204 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps_sim, list); break;
1205 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd, eenh_names[i], &(enerhist->dht->nndh), list);
1206 if (bRead) /* now allocate memory for it */
1208 snew(enerhist->dht->dh, enerhist->dht->nndh);
1209 snew(enerhist->dht->ndh, enerhist->dht->nndh);
1210 for (j = 0; j < enerhist->dht->nndh; j++)
1212 enerhist->dht->ndh[j] = 0;
1213 enerhist->dht->dh[j] = NULL;
1217 case eenhENERGY_DELTA_H_LIST:
1218 for (j = 0; j < enerhist->dht->nndh; j++)
1220 ret = do_cpte_n_reals(xd, cptpEENH, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list);
1223 case eenhENERGY_DELTA_H_STARTTIME:
1224 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_time), list); break;
1225 case eenhENERGY_DELTA_H_STARTLAMBDA:
1226 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_lambda), list); break;
1228 gmx_fatal(FARGS, "Unknown energy history entry %d\n"
1229 "You are probably reading a new checkpoint file with old code", i);
1234 if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
1236 /* Assume we have an old file format and copy sum to sum_sim */
1237 srenew(enerhist->ener_sum_sim, enerhist->nener);
1238 for (i = 0; i < enerhist->nener; i++)
1240 enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1244 if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1245 !(fflags & (1<<eenhENERGY_NSTEPS)))
1247 /* Assume we have an old file format and copy nsum to nsteps */
1248 enerhist->nsteps = enerhist->nsum;
1250 if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1251 !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1253 /* Assume we have an old file format and copy nsum to nsteps */
1254 enerhist->nsteps_sim = enerhist->nsum_sim;
1260 static int do_cpt_df_hist(XDR *xd, int fflags, df_history_t *dfhist, FILE *list)
1265 nlambda = dfhist->nlambda;
1268 for (i = 0; (i < edfhNR && ret == 0); i++)
1270 if (fflags & (1<<i))
1274 case edfhBEQUIL: ret = do_cpte_int(xd, cptpEDFH, i, fflags, &dfhist->bEquil, list); break;
1275 case edfhNATLAMBDA: ret = do_cpte_ints(xd, cptpEDFH, i, fflags, nlambda, &dfhist->n_at_lam, list); break;
1276 case edfhWLHISTO: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->wl_histo, list); break;
1277 case edfhWLDELTA: ret = do_cpte_real(xd, cptpEDFH, i, fflags, &dfhist->wl_delta, list); break;
1278 case edfhSUMWEIGHTS: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_weights, list); break;
1279 case edfhSUMDG: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_dg, list); break;
1280 case edfhSUMMINVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_minvar, list); break;
1281 case edfhSUMVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_variance, list); break;
1282 case edfhACCUMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p, list); break;
1283 case edfhACCUMM: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m, list); break;
1284 case edfhACCUMP2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p2, list); break;
1285 case edfhACCUMM2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m2, list); break;
1286 case edfhTIJ: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij, list); break;
1287 case edfhTIJEMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij_empirical, list); break;
1290 gmx_fatal(FARGS, "Unknown df history entry %d\n"
1291 "You are probably reading a new checkpoint file with old code", i);
1300 /* This function stores the last whole configuration of the reference and
1301 * average structure in the .cpt file
1303 static int do_cpt_EDstate(XDR *xd, gmx_bool bRead,
1304 edsamstate_t *EDstate, FILE *list)
1311 EDstate->bFromCpt = bRead;
1313 if (EDstate->nED <= 0)
1318 /* When reading, init_edsam has not been called yet,
1319 * so we have to allocate memory first. */
1322 snew(EDstate->nref, EDstate->nED);
1323 snew(EDstate->old_sref, EDstate->nED);
1324 snew(EDstate->nav, EDstate->nED);
1325 snew(EDstate->old_sav, EDstate->nED);
1328 /* Read/write the last whole conformation of SREF and SAV for each ED dataset (usually only one) */
1329 for (i = 0; i < EDstate->nED; i++)
1331 /* Reference structure SREF */
1332 sprintf(buf, "ED%d # of atoms in reference structure", i+1);
1333 do_cpt_int_err(xd, buf, &EDstate->nref[i], list);
1334 sprintf(buf, "ED%d x_ref", i+1);
1337 snew(EDstate->old_sref[i], EDstate->nref[i]);
1338 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref[i], list);
1342 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref_p[i], list);
1345 /* Average structure SAV */
1346 sprintf(buf, "ED%d # of atoms in average structure", i+1);
1347 do_cpt_int_err(xd, buf, &EDstate->nav[i], list);
1348 sprintf(buf, "ED%d x_av", i+1);
1351 snew(EDstate->old_sav[i], EDstate->nav[i]);
1352 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav[i], list);
1356 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav_p[i], list);
1364 static int do_cpt_files(XDR *xd, gmx_bool bRead,
1365 gmx_file_position_t **p_outputfiles, int *nfiles,
1366 FILE *list, int file_version)
1370 gmx_off_t mask = 0xFFFFFFFFL;
1371 int offset_high, offset_low;
1373 gmx_file_position_t *outputfiles;
1375 if (do_cpt_int(xd, "number of output files", nfiles, list) != 0)
1382 snew(*p_outputfiles, *nfiles);
1385 outputfiles = *p_outputfiles;
1387 for (i = 0; i < *nfiles; i++)
1389 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1392 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1393 strncpy(outputfiles[i].filename, buf, CPTSTRLEN-1);
1399 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1403 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1407 outputfiles[i].offset = ( ((gmx_off_t) offset_high) << 32 ) | ( (gmx_off_t) offset_low & mask );
1411 buf = outputfiles[i].filename;
1412 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1414 offset = outputfiles[i].offset;
1422 offset_low = (int) (offset & mask);
1423 offset_high = (int) ((offset >> 32) & mask);
1425 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1429 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1434 if (file_version >= 8)
1436 if (do_cpt_int(xd, "file_checksum_size", &(outputfiles[i].chksum_size),
1441 if (do_cpt_u_chars(xd, "file_checksum", 16, outputfiles[i].chksum, list) != 0)
1448 outputfiles[i].chksum_size = -1;
1455 void write_checkpoint(const char *fn, gmx_bool bNumberAndKeep,
1456 FILE *fplog, t_commrec *cr,
1457 int eIntegrator, int simulation_part,
1458 gmx_bool bExpanded, int elamstats,
1459 gmx_int64_t step, double t, t_state *state)
1469 char *fntemp; /* the temporary checkpoint file name */
1470 char timebuf[STRLEN];
1471 int nppnodes, npmenodes;
1472 char buf[1024], suffix[5+STEPSTRSIZE], sbuf[STEPSTRSIZE];
1473 gmx_file_position_t *outputfiles;
1476 int flags_eks, flags_enh, flags_dfh;
1479 if (DOMAINDECOMP(cr))
1481 nppnodes = cr->dd->nnodes;
1482 npmenodes = cr->npmenodes;
1491 /* make the new temporary filename */
1492 snew(fntemp, strlen(fn)+5+STEPSTRSIZE);
1494 fntemp[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1495 sprintf(suffix, "_%s%s", "step", gmx_step_str(step, sbuf));
1496 strcat(fntemp, suffix);
1497 strcat(fntemp, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1499 /* if we can't rename, we just overwrite the cpt file.
1500 * dangerous if interrupted.
1502 snew(fntemp, strlen(fn));
1505 gmx_format_current_time(timebuf, STRLEN);
1509 fprintf(fplog, "Writing checkpoint, step %s at %s\n\n",
1510 gmx_step_str(step, buf), timebuf);
1513 /* Get offsets for open files */
1514 gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1516 fp = gmx_fio_open(fntemp, "w");
1518 if (state->ekinstate.bUpToDate)
1521 ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) |
1522 (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) |
1523 (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1531 if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1533 flags_enh |= (1<<eenhENERGY_N) | (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSTEPS_SIM);
1534 if (state->enerhist.nsum > 0)
1536 flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1537 (1<<eenhENERGY_NSUM));
1539 if (state->enerhist.nsum_sim > 0)
1541 flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (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 */
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?");