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)
196 if (xdr_string(xd, s, CPTSTRLEN) == 0)
202 fprintf(list, "%s = %s\n", desc, *s);
207 static int do_cpt_int(XDR *xd, const char *desc, int *i, FILE *list)
209 if (xdr_int(xd, i) == 0)
215 fprintf(list, "%s = %d\n", desc, *i);
220 static int do_cpt_u_chars(XDR *xd, const char *desc, int n, unsigned char *i, FILE *list)
224 fprintf(list, "%s = ", desc);
227 for (int j = 0; j < n && res; j++)
229 res &= xdr_u_char(xd, &i[j]);
232 fprintf(list, "%02x", i[j]);
247 static void do_cpt_int_err(XDR *xd, const char *desc, int *i, FILE *list)
249 if (do_cpt_int(xd, desc, i, list) < 0)
255 static void do_cpt_step_err(XDR *xd, const char *desc, gmx_int64_t *i, FILE *list)
257 char buf[STEPSTRSIZE];
259 if (xdr_int64(xd, i) == 0)
265 fprintf(list, "%s = %s\n", desc, gmx_step_str(*i, buf));
269 static void do_cpt_double_err(XDR *xd, const char *desc, double *f, FILE *list)
271 if (xdr_double(xd, f) == 0)
277 fprintf(list, "%s = %f\n", desc, *f);
281 static void do_cpt_real_err(XDR *xd, real *f)
284 bool_t res = xdr_double(xd, f);
286 bool_t res = xdr_float(xd, f);
294 static void do_cpt_n_rvecs_err(XDR *xd, const char *desc, int n, rvec f[], FILE *list)
296 for (int i = 0; i < n; i++)
298 for (int j = 0; j < DIM; j++)
300 do_cpt_real_err(xd, &f[i][j]);
306 pr_rvecs(list, 0, desc, f, n);
310 /* If nval >= 0, nval is used; on read this should match the passed value.
311 * If nval n<0, *nptr is used; on read the value is stored in nptr
313 static int do_cpte_reals_low(XDR *xd, int cptp, int ecpt, int sflags,
314 int nval, int *nptr, real **v,
315 FILE *list, int erealtype)
319 int dtc = xdr_datatype_float;
321 int dtc = xdr_datatype_double;
323 real *vp, *va = NULL;
338 gmx_incons("*ntpr=NULL in do_cpte_reals_low");
343 res = xdr_int(xd, &nf);
354 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), nval, nf);
363 res = xdr_int(xd, &dt);
370 fprintf(stderr, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
371 st_names(cptp, ecpt), xdr_datatype_names[dtc],
372 xdr_datatype_names[dt]);
374 if (list || !(sflags & (1<<ecpt)))
387 if (dt == xdr_datatype_float)
389 if (dtc == xdr_datatype_float)
391 vf = reinterpret_cast<float *>(vp);
397 res = xdr_vector(xd, reinterpret_cast<char *>(vf), nf,
398 static_cast<unsigned int>(sizeof(float)), (xdrproc_t)xdr_float);
403 if (dtc != xdr_datatype_float)
405 for (i = 0; i < nf; i++)
414 if (dtc == xdr_datatype_double)
416 /* cppcheck-suppress invalidPointerCast
417 * Only executed if real is anyhow double */
424 res = xdr_vector(xd, reinterpret_cast<char *>(vd), nf,
425 static_cast<unsigned int>(sizeof(double)), (xdrproc_t)xdr_double);
430 if (dtc != xdr_datatype_double)
432 for (i = 0; i < nf; i++)
445 pr_reals(list, 0, st_names(cptp, ecpt), vp, nf);
448 pr_rvecs(list, 0, st_names(cptp, ecpt), (rvec *)vp, nf/3);
451 gmx_incons("Unknown checkpoint real type");
463 /* This function stores n along with the reals for reading,
464 * but on reading it assumes that n matches the value in the checkpoint file,
465 * a fatal error is generated when this is not the case.
467 static int do_cpte_reals(XDR *xd, int cptp, int ecpt, int sflags,
468 int n, real **v, FILE *list)
470 return do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, v, list, ecprREAL);
473 /* This function does the same as do_cpte_reals,
474 * except that on reading it ignores the passed value of *n
475 * and stored the value read from the checkpoint file in *n.
477 static int do_cpte_n_reals(XDR *xd, int cptp, int ecpt, int sflags,
478 int *n, real **v, FILE *list)
480 return do_cpte_reals_low(xd, cptp, ecpt, sflags, -1, n, v, list, ecprREAL);
483 static int do_cpte_real(XDR *xd, int cptp, int ecpt, int sflags,
486 return do_cpte_reals_low(xd, cptp, ecpt, sflags, 1, NULL, &r, list, ecprREAL);
489 static int do_cpte_ints(XDR *xd, int cptp, int ecpt, int sflags,
490 int n, int **v, FILE *list)
493 int dtc = xdr_datatype_int;
498 res = xdr_int(xd, &nf);
503 if (list == NULL && v != NULL && nf != n)
505 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
508 res = xdr_int(xd, &dt);
515 gmx_fatal(FARGS, "Type mismatch for state entry %s, code type is %s, file type is %s\n",
516 st_names(cptp, ecpt), xdr_datatype_names[dtc],
517 xdr_datatype_names[dt]);
519 if (list || !(sflags & (1<<ecpt)) || v == NULL)
532 res = xdr_vector(xd, reinterpret_cast<char *>(vp), nf,
533 static_cast<unsigned int>(sizeof(int)), (xdrproc_t)xdr_int);
540 pr_ivec(list, 0, st_names(cptp, ecpt), vp, nf, TRUE);
550 static int do_cpte_int(XDR *xd, int cptp, int ecpt, int sflags,
553 return do_cpte_ints(xd, cptp, ecpt, sflags, 1, &i, list);
556 static int do_cpte_doubles(XDR *xd, int cptp, int ecpt, int sflags,
557 int n, double **v, FILE *list)
560 int dtc = xdr_datatype_double;
561 double *vp, *va = NULL;
565 res = xdr_int(xd, &nf);
570 if (list == NULL && nf != n)
572 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
575 res = xdr_int(xd, &dt);
582 gmx_fatal(FARGS, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
583 st_names(cptp, ecpt), xdr_datatype_names[dtc],
584 xdr_datatype_names[dt]);
586 if (list || !(sflags & (1<<ecpt)))
599 res = xdr_vector(xd, reinterpret_cast<char *>(vp), nf,
600 static_cast<unsigned int>(sizeof(double)), (xdrproc_t)xdr_double);
607 pr_doubles(list, 0, st_names(cptp, ecpt), vp, nf);
617 static int do_cpte_double(XDR *xd, int cptp, int ecpt, int sflags,
618 double *r, FILE *list)
620 return do_cpte_doubles(xd, cptp, ecpt, sflags, 1, &r, list);
624 static int do_cpte_rvecs(XDR *xd, int cptp, int ecpt, int sflags,
625 int n, rvec **v, FILE *list)
627 return do_cpte_reals_low(xd, cptp, ecpt, sflags,
628 n*DIM, NULL, (real **)v, list, ecprRVEC);
631 static int do_cpte_matrix(XDR *xd, int cptp, int ecpt, int sflags,
632 matrix v, FILE *list)
638 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
639 DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
641 if (list && ret == 0)
643 pr_rvecs(list, 0, st_names(cptp, ecpt), v, DIM);
650 static int do_cpte_nmatrix(XDR *xd, int cptp, int ecpt, int sflags,
651 int n, real **v, FILE *list)
655 char name[CPTSTRLEN];
662 for (i = 0; i < n; i++)
664 reti = do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, &(v[i]), NULL, ecprREAL);
665 if (list && reti == 0)
667 sprintf(name, "%s[%d]", st_names(cptp, ecpt), i);
668 pr_reals(list, 0, name, v[i], n);
678 static int do_cpte_matrices(XDR *xd, int cptp, int ecpt, int sflags,
679 int n, matrix **v, FILE *list)
682 matrix *vp, *va = NULL;
688 res = xdr_int(xd, &nf);
693 if (list == NULL && nf != n)
695 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
697 if (list || !(sflags & (1<<ecpt)))
710 snew(vr, nf*DIM*DIM);
711 for (i = 0; i < nf; i++)
713 for (j = 0; j < DIM; j++)
715 for (k = 0; k < DIM; k++)
717 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
721 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
722 nf*DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
723 for (i = 0; i < nf; i++)
725 for (j = 0; j < DIM; j++)
727 for (k = 0; k < DIM; k++)
729 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
735 if (list && ret == 0)
737 for (i = 0; i < nf; i++)
739 pr_rvecs(list, 0, st_names(cptp, ecpt), vp[i], DIM);
750 static void do_cpt_header(XDR *xd, gmx_bool bRead, int *file_version,
751 char **version, char **btime, char **buser, char **bhost,
753 char **fprog, char **ftime,
754 int *eIntegrator, int *simulation_part,
755 gmx_int64_t *step, double *t,
756 int *nnodes, int *dd_nc, int *npme,
757 int *natoms, int *ngtc, int *nnhpres, int *nhchainlength,
758 int *nlambda, int *flags_state,
759 int *flags_eks, int *flags_enh, int *flags_dfh,
760 int *nED, int *eSwapCoords,
776 res = xdr_int(xd, &magic);
779 gmx_fatal(FARGS, "The checkpoint file is empty/corrupted, or maybe you are out of disk space?");
781 if (magic != CPT_MAGIC1)
783 gmx_fatal(FARGS, "Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
784 "The checkpoint file is corrupted or not a checkpoint file",
790 gmx_gethostname(fhost, 255);
792 do_cpt_string_err(xd, bRead, "GROMACS version", version, list);
793 do_cpt_string_err(xd, bRead, "GROMACS build time", btime, list);
794 do_cpt_string_err(xd, bRead, "GROMACS build user", buser, list);
795 do_cpt_string_err(xd, bRead, "GROMACS build host", bhost, list);
796 do_cpt_string_err(xd, bRead, "generating program", fprog, list);
797 do_cpt_string_err(xd, bRead, "generation time", ftime, list);
798 *file_version = cpt_version;
799 do_cpt_int_err(xd, "checkpoint file version", file_version, list);
800 if (*file_version > cpt_version)
802 gmx_fatal(FARGS, "Attempting to read a checkpoint file of version %d with code of version %d\n", *file_version, cpt_version);
804 if (*file_version >= 13)
806 do_cpt_int_err(xd, "GROMACS double precision", double_prec, list);
812 if (*file_version >= 12)
814 do_cpt_string_err(xd, bRead, "generating host", &fhost, list);
820 do_cpt_int_err(xd, "#atoms", natoms, list);
821 do_cpt_int_err(xd, "#T-coupling groups", ngtc, list);
822 if (*file_version >= 10)
824 do_cpt_int_err(xd, "#Nose-Hoover T-chains", nhchainlength, list);
830 if (*file_version >= 11)
832 do_cpt_int_err(xd, "#Nose-Hoover T-chains for barostat ", nnhpres, list);
838 if (*file_version >= 14)
840 do_cpt_int_err(xd, "# of total lambda states ", nlambda, list);
846 do_cpt_int_err(xd, "integrator", eIntegrator, list);
847 if (*file_version >= 3)
849 do_cpt_int_err(xd, "simulation part #", simulation_part, list);
853 *simulation_part = 1;
855 if (*file_version >= 5)
857 do_cpt_step_err(xd, "step", step, list);
861 do_cpt_int_err(xd, "step", &idum, list);
864 do_cpt_double_err(xd, "t", t, list);
865 do_cpt_int_err(xd, "#PP-ranks", nnodes, list);
867 do_cpt_int_err(xd, "dd_nc[x]", dd_nc ? &(dd_nc[0]) : &idum, list);
868 do_cpt_int_err(xd, "dd_nc[y]", dd_nc ? &(dd_nc[1]) : &idum, list);
869 do_cpt_int_err(xd, "dd_nc[z]", dd_nc ? &(dd_nc[2]) : &idum, list);
870 do_cpt_int_err(xd, "#PME-only ranks", npme, list);
871 do_cpt_int_err(xd, "state flags", flags_state, list);
872 if (*file_version >= 4)
874 do_cpt_int_err(xd, "ekin data flags", flags_eks, list);
875 do_cpt_int_err(xd, "energy history flags", flags_enh, list);
880 *flags_enh = (*flags_state >> (estORIRE_DTAV+1));
881 *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
882 (1<<(estORIRE_DTAV+2)) |
883 (1<<(estORIRE_DTAV+3))));
885 if (*file_version >= 14)
887 do_cpt_int_err(xd, "df history flags", flags_dfh, list);
894 if (*file_version >= 15)
896 do_cpt_int_err(xd, "ED data sets", nED, list);
902 if (*file_version >= 16)
904 do_cpt_int_err(xd, "swap", eSwapCoords, list);
908 static int do_cpt_footer(XDR *xd, int file_version)
913 if (file_version >= 2)
916 res = xdr_int(xd, &magic);
921 if (magic != CPT_MAGIC2)
930 static int do_cpt_state(XDR *xd, gmx_bool bRead,
931 int fflags, t_state *state,
941 nnht = state->nhchainlength*state->ngtc;
942 nnhtp = state->nhchainlength*state->nnhpres;
944 if (bRead) /* we need to allocate space for dfhist if we are reading */
946 init_df_history(&state->dfhist, state->dfhist.nlambda);
949 sflags = state->flags;
950 for (i = 0; (i < estNR && ret == 0); i++)
956 case estLAMBDA: ret = do_cpte_reals(xd, cptpEST, i, sflags, efptNR, &(state->lambda), list); break;
957 case estFEPSTATE: ret = do_cpte_int (xd, cptpEST, i, sflags, &state->fep_state, list); break;
958 case estBOX: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box, list); break;
959 case estBOX_REL: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box_rel, list); break;
960 case estBOXV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->boxv, list); break;
961 case estPRES_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->pres_prev, list); break;
962 case estSVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->svir_prev, list); break;
963 case estFVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->fvir_prev, list); break;
964 case estNH_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_xi, list); break;
965 case estNH_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_vxi, list); break;
966 case estNHPRES_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_xi, list); break;
967 case estNHPRES_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_vxi, list); break;
968 case estTC_INT: ret = do_cpte_doubles(xd, cptpEST, i, sflags, state->ngtc, &state->therm_integral, list); break;
969 case estVETA: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->veta, list); break;
970 case estVOL0: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->vol0, list); break;
971 case estX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->x, list); break;
972 case estV: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->v, list); break;
973 case estSDX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->sd_X, list); break;
974 /* The RNG entries are no longer written,
975 * the next 4 lines are only for reading old files.
977 case estLD_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
978 case estLD_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
979 case estMC_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
980 case estMC_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
981 case estDISRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.disre_initf, list); break;
982 case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.ndisrepairs, &state->hist.disre_rm3tav, list); break;
983 case estORIRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.orire_initf, list); break;
984 case estORIRE_DTAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.norire_Dtav, &state->hist.orire_Dtav, list); break;
986 gmx_fatal(FARGS, "Unknown state entry %d\n"
987 "You are probably reading a new checkpoint file with old code", i);
995 static int do_cpt_ekinstate(XDR *xd, int fflags, ekinstate_t *ekins,
1003 for (i = 0; (i < eeksNR && ret == 0); i++)
1005 if (fflags & (1<<i))
1010 case eeksEKIN_N: ret = do_cpte_int(xd, cptpEEKS, i, fflags, &ekins->ekin_n, list); break;
1011 case eeksEKINH: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh, list); break;
1012 case eeksEKINF: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinf, list); break;
1013 case eeksEKINO: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh_old, list); break;
1014 case eeksEKINTOTAL: ret = do_cpte_matrix(xd, cptpEEKS, i, fflags, ekins->ekin_total, list); break;
1015 case eeksEKINSCALEF: ret = do_cpte_doubles(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinscalef_nhc, list); break;
1016 case eeksVSCALE: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->vscale_nhc, list); break;
1017 case eeksEKINSCALEH: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->ekinscaleh_nhc, list); break;
1018 case eeksDEKINDL: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->dekindl, list); break;
1019 case eeksMVCOS: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->mvcos, list); break;
1021 gmx_fatal(FARGS, "Unknown ekin data state entry %d\n"
1022 "You are probably reading a new checkpoint file with old code", i);
1031 static int do_cpt_swapstate(XDR *xd, gmx_bool bRead, swapstate_t *swapstate, FILE *list)
1035 int swap_cpt_version = 1;
1038 if (eswapNO == swapstate->eSwapCoords)
1043 swapstate->bFromCpt = bRead;
1045 do_cpt_int_err(xd, "swap checkpoint version", &swap_cpt_version, list);
1046 do_cpt_int_err(xd, "swap coupling steps", &swapstate->nAverage, list);
1048 /* When reading, init_swapcoords has not been called yet,
1049 * so we have to allocate memory first. */
1051 for (ic = 0; ic < eCompNR; ic++)
1053 for (ii = 0; ii < eIonNR; ii++)
1057 do_cpt_int_err(xd, "swap requested atoms", &swapstate->nat_req[ic][ii], list);
1061 do_cpt_int_err(xd, "swap requested atoms p", swapstate->nat_req_p[ic][ii], list);
1066 do_cpt_int_err(xd, "swap influx netto", &swapstate->inflow_netto[ic][ii], list);
1070 do_cpt_int_err(xd, "swap influx netto p", swapstate->inflow_netto_p[ic][ii], list);
1073 if (bRead && (NULL == swapstate->nat_past[ic][ii]) )
1075 snew(swapstate->nat_past[ic][ii], swapstate->nAverage);
1078 for (j = 0; j < swapstate->nAverage; j++)
1082 do_cpt_int_err(xd, "swap past atom counts", &swapstate->nat_past[ic][ii][j], list);
1086 do_cpt_int_err(xd, "swap past atom counts p", &swapstate->nat_past_p[ic][ii][j], list);
1092 /* Ion flux per channel */
1093 for (ic = 0; ic < eChanNR; ic++)
1095 for (ii = 0; ii < eIonNR; ii++)
1099 do_cpt_int_err(xd, "channel flux", &swapstate->fluxfromAtoB[ic][ii], list);
1103 do_cpt_int_err(xd, "channel flux p", swapstate->fluxfromAtoB_p[ic][ii], list);
1108 /* Ion flux leakage */
1111 snew(swapstate->fluxleak, 1);
1113 do_cpt_int_err(xd, "flux leakage", swapstate->fluxleak, list);
1116 do_cpt_int_err(xd, "number of ions", &swapstate->nions, list);
1120 snew(swapstate->channel_label, swapstate->nions);
1121 snew(swapstate->comp_from, swapstate->nions);
1124 do_cpt_u_chars(xd, "channel history", swapstate->nions, swapstate->channel_label, list);
1125 do_cpt_u_chars(xd, "domain history", swapstate->nions, swapstate->comp_from, list);
1127 /* Save the last known whole positions to checkpoint
1128 * file to be able to also make multimeric channels whole in PBC */
1129 do_cpt_int_err(xd, "Ch0 atoms", &swapstate->nat[eChan0], list);
1130 do_cpt_int_err(xd, "Ch1 atoms", &swapstate->nat[eChan1], list);
1133 snew(swapstate->xc_old_whole[eChan0], swapstate->nat[eChan0]);
1134 snew(swapstate->xc_old_whole[eChan1], swapstate->nat[eChan1]);
1135 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], swapstate->xc_old_whole[eChan0], list);
1136 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], swapstate->xc_old_whole[eChan1], list);
1140 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], *swapstate->xc_old_whole_p[eChan0], list);
1141 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], *swapstate->xc_old_whole_p[eChan1], list);
1148 static int do_cpt_enerhist(XDR *xd, gmx_bool bRead,
1149 int fflags, energyhistory_t *enerhist,
1160 enerhist->nsteps = 0;
1162 enerhist->nsteps_sim = 0;
1163 enerhist->nsum_sim = 0;
1164 enerhist->dht = NULL;
1166 if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
1168 snew(enerhist->dht, 1);
1169 enerhist->dht->ndh = NULL;
1170 enerhist->dht->dh = NULL;
1171 enerhist->dht->start_lambda_set = FALSE;
1175 for (i = 0; (i < eenhNR && ret == 0); i++)
1177 if (fflags & (1<<i))
1181 case eenhENERGY_N: ret = do_cpte_int(xd, cptpEENH, i, fflags, &enerhist->nener, list); break;
1182 case eenhENERGY_AVER: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_ave, list); break;
1183 case eenhENERGY_SUM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum, list); break;
1184 case eenhENERGY_NSUM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum, list); break;
1185 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum_sim, list); break;
1186 case eenhENERGY_NSUM_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum_sim, list); break;
1187 case eenhENERGY_NSTEPS: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps, list); break;
1188 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps_sim, list); break;
1189 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd, eenh_names[i], &(enerhist->dht->nndh), list);
1190 if (bRead) /* now allocate memory for it */
1192 snew(enerhist->dht->dh, enerhist->dht->nndh);
1193 snew(enerhist->dht->ndh, enerhist->dht->nndh);
1194 for (j = 0; j < enerhist->dht->nndh; j++)
1196 enerhist->dht->ndh[j] = 0;
1197 enerhist->dht->dh[j] = NULL;
1201 case eenhENERGY_DELTA_H_LIST:
1202 for (j = 0; j < enerhist->dht->nndh; j++)
1204 ret = do_cpte_n_reals(xd, cptpEENH, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list);
1207 case eenhENERGY_DELTA_H_STARTTIME:
1208 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_time), list); break;
1209 case eenhENERGY_DELTA_H_STARTLAMBDA:
1210 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_lambda), list); break;
1212 gmx_fatal(FARGS, "Unknown energy history entry %d\n"
1213 "You are probably reading a new checkpoint file with old code", i);
1218 if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
1220 /* Assume we have an old file format and copy sum to sum_sim */
1221 srenew(enerhist->ener_sum_sim, enerhist->nener);
1222 for (i = 0; i < enerhist->nener; i++)
1224 enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1228 if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1229 !(fflags & (1<<eenhENERGY_NSTEPS)))
1231 /* Assume we have an old file format and copy nsum to nsteps */
1232 enerhist->nsteps = enerhist->nsum;
1234 if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1235 !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1237 /* Assume we have an old file format and copy nsum to nsteps */
1238 enerhist->nsteps_sim = enerhist->nsum_sim;
1244 static int do_cpt_df_hist(XDR *xd, int fflags, df_history_t *dfhist, FILE *list)
1249 nlambda = dfhist->nlambda;
1252 for (i = 0; (i < edfhNR && ret == 0); i++)
1254 if (fflags & (1<<i))
1258 case edfhBEQUIL: ret = do_cpte_int(xd, cptpEDFH, i, fflags, &dfhist->bEquil, list); break;
1259 case edfhNATLAMBDA: ret = do_cpte_ints(xd, cptpEDFH, i, fflags, nlambda, &dfhist->n_at_lam, list); break;
1260 case edfhWLHISTO: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->wl_histo, list); break;
1261 case edfhWLDELTA: ret = do_cpte_real(xd, cptpEDFH, i, fflags, &dfhist->wl_delta, list); break;
1262 case edfhSUMWEIGHTS: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_weights, list); break;
1263 case edfhSUMDG: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_dg, list); break;
1264 case edfhSUMMINVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_minvar, list); break;
1265 case edfhSUMVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_variance, list); break;
1266 case edfhACCUMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p, list); break;
1267 case edfhACCUMM: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m, list); break;
1268 case edfhACCUMP2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p2, list); break;
1269 case edfhACCUMM2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m2, list); break;
1270 case edfhTIJ: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij, list); break;
1271 case edfhTIJEMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij_empirical, list); break;
1274 gmx_fatal(FARGS, "Unknown df history entry %d\n"
1275 "You are probably reading a new checkpoint file with old code", i);
1284 /* This function stores the last whole configuration of the reference and
1285 * average structure in the .cpt file
1287 static int do_cpt_EDstate(XDR *xd, gmx_bool bRead,
1288 edsamstate_t *EDstate, FILE *list)
1295 EDstate->bFromCpt = bRead;
1297 if (EDstate->nED <= 0)
1302 /* When reading, init_edsam has not been called yet,
1303 * so we have to allocate memory first. */
1306 snew(EDstate->nref, EDstate->nED);
1307 snew(EDstate->old_sref, EDstate->nED);
1308 snew(EDstate->nav, EDstate->nED);
1309 snew(EDstate->old_sav, EDstate->nED);
1312 /* Read/write the last whole conformation of SREF and SAV for each ED dataset (usually only one) */
1313 for (i = 0; i < EDstate->nED; i++)
1315 /* Reference structure SREF */
1316 sprintf(buf, "ED%d # of atoms in reference structure", i+1);
1317 do_cpt_int_err(xd, buf, &EDstate->nref[i], list);
1318 sprintf(buf, "ED%d x_ref", i+1);
1321 snew(EDstate->old_sref[i], EDstate->nref[i]);
1322 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref[i], list);
1326 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref_p[i], list);
1329 /* Average structure SAV */
1330 sprintf(buf, "ED%d # of atoms in average structure", i+1);
1331 do_cpt_int_err(xd, buf, &EDstate->nav[i], list);
1332 sprintf(buf, "ED%d x_av", i+1);
1335 snew(EDstate->old_sav[i], EDstate->nav[i]);
1336 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav[i], list);
1340 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav_p[i], list);
1348 static int do_cpt_files(XDR *xd, gmx_bool bRead,
1349 gmx_file_position_t **p_outputfiles, int *nfiles,
1350 FILE *list, int file_version)
1354 gmx_off_t mask = 0xFFFFFFFFL;
1355 int offset_high, offset_low;
1357 gmx_file_position_t *outputfiles;
1359 if (do_cpt_int(xd, "number of output files", nfiles, list) != 0)
1366 snew(*p_outputfiles, *nfiles);
1369 outputfiles = *p_outputfiles;
1371 for (i = 0; i < *nfiles; i++)
1373 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1376 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1377 std::strncpy(outputfiles[i].filename, buf, CPTSTRLEN-1);
1383 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1387 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1391 outputfiles[i].offset = (static_cast<gmx_off_t>(offset_high) << 32 ) | ( static_cast<gmx_off_t>(offset_low) & mask );
1395 buf = outputfiles[i].filename;
1396 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1398 offset = outputfiles[i].offset;
1406 offset_low = static_cast<int>(offset & mask);
1407 offset_high = static_cast<int>((offset >> 32) & mask);
1409 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1413 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1418 if (file_version >= 8)
1420 if (do_cpt_int(xd, "file_checksum_size", &(outputfiles[i].chksum_size),
1425 if (do_cpt_u_chars(xd, "file_checksum", 16, outputfiles[i].chksum, list) != 0)
1432 outputfiles[i].chksum_size = -1;
1439 void write_checkpoint(const char *fn, gmx_bool bNumberAndKeep,
1440 FILE *fplog, t_commrec *cr,
1441 int eIntegrator, int simulation_part,
1442 gmx_bool bExpanded, int elamstats,
1443 gmx_int64_t step, double t, t_state *state)
1453 char *fntemp; /* the temporary checkpoint file name */
1454 char timebuf[STRLEN];
1455 int nppnodes, npmenodes;
1456 char buf[1024], suffix[5+STEPSTRSIZE], sbuf[STEPSTRSIZE];
1457 gmx_file_position_t *outputfiles;
1460 int flags_eks, flags_enh, flags_dfh;
1463 if (DOMAINDECOMP(cr))
1465 nppnodes = cr->dd->nnodes;
1466 npmenodes = cr->npmenodes;
1475 /* make the new temporary filename */
1476 snew(fntemp, std::strlen(fn)+5+STEPSTRSIZE);
1477 std::strcpy(fntemp, fn);
1478 fntemp[std::strlen(fn) - std::strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1479 sprintf(suffix, "_%s%s", "step", gmx_step_str(step, sbuf));
1480 std::strcat(fntemp, suffix);
1481 std::strcat(fntemp, fn+std::strlen(fn) - std::strlen(ftp2ext(fn2ftp(fn))) - 1);
1483 /* if we can't rename, we just overwrite the cpt file.
1484 * dangerous if interrupted.
1486 snew(fntemp, std::strlen(fn));
1487 std::strcpy(fntemp, fn);
1489 gmx_format_current_time(timebuf, STRLEN);
1493 fprintf(fplog, "Writing checkpoint, step %s at %s\n\n",
1494 gmx_step_str(step, buf), timebuf);
1497 /* Get offsets for open files */
1498 gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1500 fp = gmx_fio_open(fntemp, "w");
1502 if (state->ekinstate.bUpToDate)
1505 ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) |
1506 (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) |
1507 (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1515 if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1517 flags_enh |= (1<<eenhENERGY_N) | (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSTEPS_SIM);
1518 if (state->enerhist.nsum > 0)
1520 flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1521 (1<<eenhENERGY_NSUM));
1523 if (state->enerhist.nsum_sim > 0)
1525 flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSUM_SIM));
1527 if (state->enerhist.dht)
1529 flags_enh |= ( (1<< eenhENERGY_DELTA_H_NN) |
1530 (1<< eenhENERGY_DELTA_H_LIST) |
1531 (1<< eenhENERGY_DELTA_H_STARTTIME) |
1532 (1<< eenhENERGY_DELTA_H_STARTLAMBDA) );
1538 flags_dfh = ((1<<edfhBEQUIL) | (1<<edfhNATLAMBDA) | (1<<edfhSUMWEIGHTS) | (1<<edfhSUMDG) |
1539 (1<<edfhTIJ) | (1<<edfhTIJEMP));
1542 flags_dfh |= ((1<<edfhWLDELTA) | (1<<edfhWLHISTO));
1544 if ((elamstats == elamstatsMINVAR) || (elamstats == elamstatsBARKER) || (elamstats == elamstatsMETROPOLIS))
1546 flags_dfh |= ((1<<edfhACCUMP) | (1<<edfhACCUMM) | (1<<edfhACCUMP2) | (1<<edfhACCUMM2)
1547 | (1<<edfhSUMMINVAR) | (1<<edfhSUMVAR));
1555 /* We can check many more things now (CPU, acceleration, etc), but
1556 * it is highly unlikely to have two separate builds with exactly
1557 * the same version, user, time, and build host!
1560 version = gmx_strdup(gmx_version());
1561 btime = gmx_strdup(BUILD_TIME);
1562 buser = gmx_strdup(BUILD_USER);
1563 bhost = gmx_strdup(BUILD_HOST);
1565 double_prec = GMX_CPT_BUILD_DP;
1566 fprog = gmx_strdup(Program());
1568 ftime = &(timebuf[0]);
1570 do_cpt_header(gmx_fio_getxdr(fp), FALSE, &file_version,
1571 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1572 &eIntegrator, &simulation_part, &step, &t, &nppnodes,
1573 DOMAINDECOMP(cr) ? cr->dd->nc : NULL, &npmenodes,
1574 &state->natoms, &state->ngtc, &state->nnhpres,
1575 &state->nhchainlength, &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
1576 &state->edsamstate.nED, &state->swapstate.eSwapCoords,
1585 if ((do_cpt_state(gmx_fio_getxdr(fp), FALSE, state->flags, state, NULL) < 0) ||
1586 (do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL) < 0) ||
1587 (do_cpt_enerhist(gmx_fio_getxdr(fp), FALSE, flags_enh, &state->enerhist, NULL) < 0) ||
1588 (do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL) < 0) ||
1589 (do_cpt_EDstate(gmx_fio_getxdr(fp), FALSE, &state->edsamstate, NULL) < 0) ||
1590 (do_cpt_swapstate(gmx_fio_getxdr(fp), FALSE, &state->swapstate, NULL) < 0) ||
1591 (do_cpt_files(gmx_fio_getxdr(fp), FALSE, &outputfiles, &noutputfiles, NULL,
1594 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1597 do_cpt_footer(gmx_fio_getxdr(fp), file_version);
1599 /* we really, REALLY, want to make sure to physically write the checkpoint,
1600 and all the files it depends on, out to disk. Because we've
1601 opened the checkpoint with gmx_fio_open(), it's in our list
1603 ret = gmx_fio_all_output_fsync();
1609 "Cannot fsync '%s'; maybe you are out of disk space?",
1610 gmx_fio_getname(ret));
1612 if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV) == NULL)
1622 if (gmx_fio_close(fp) != 0)
1624 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1627 /* we don't move the checkpoint if the user specified they didn't want it,
1628 or if the fsyncs failed */
1630 if (!bNumberAndKeep && !ret)
1634 /* Rename the previous checkpoint file */
1635 std::strcpy(buf, fn);
1636 buf[std::strlen(fn) - std::strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1637 std::strcat(buf, "_prev");
1638 std::strcat(buf, fn+std::strlen(fn) - std::strlen(ftp2ext(fn2ftp(fn))) - 1);
1640 /* we copy here so that if something goes wrong between now and
1641 * the rename below, there's always a state.cpt.
1642 * If renames are atomic (such as in POSIX systems),
1643 * this copying should be unneccesary.
1645 gmx_file_copy(fn, buf, FALSE);
1646 /* We don't really care if this fails:
1647 * there's already a new checkpoint.
1650 gmx_file_rename(fn, buf);
1653 if (gmx_file_rename(fntemp, fn) != 0)
1655 gmx_file("Cannot rename checkpoint file; maybe you are out of disk space?");
1658 #endif /* GMX_NO_RENAME */
1664 /*code for alternate checkpointing scheme. moved from top of loop over
1666 fcRequestCheckPoint();
1667 if (fcCheckPointParallel( cr->nodeid, NULL, 0) == 0)
1669 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", step );
1671 #endif /* end GMX_FAHCORE block */
1674 static void print_flag_mismatch(FILE *fplog, int sflags, int fflags)
1678 fprintf(fplog, "\nState entry mismatch between the simulation and the checkpoint file\n");
1679 fprintf(fplog, "Entries which are not present in the checkpoint file will not be updated\n");
1680 fprintf(fplog, " %24s %11s %11s\n", "", "simulation", "checkpoint");
1681 for (i = 0; i < estNR; i++)
1683 if ((sflags & (1<<i)) || (fflags & (1<<i)))
1685 fprintf(fplog, " %24s %11s %11s\n",
1687 (sflags & (1<<i)) ? " present " : "not present",
1688 (fflags & (1<<i)) ? " present " : "not present");
1693 static void check_int(FILE *fplog, const char *type, int p, int f, gmx_bool *mm)
1695 FILE *fp = fplog ? fplog : stderr;
1699 fprintf(fp, " %s mismatch,\n", type);
1700 fprintf(fp, " current program: %d\n", p);
1701 fprintf(fp, " checkpoint file: %d\n", f);
1707 static void check_string(FILE *fplog, const char *type, const char *p,
1708 const char *f, gmx_bool *mm)
1710 FILE *fp = fplog ? fplog : stderr;
1712 if (std::strcmp(p, f) != 0)
1714 fprintf(fp, " %s mismatch,\n", type);
1715 fprintf(fp, " current program: %s\n", p);
1716 fprintf(fp, " checkpoint file: %s\n", f);
1722 static void check_match(FILE *fplog,
1724 char *btime, char *buser, char *bhost, int double_prec,
1726 t_commrec *cr, int npp_f, int npme_f,
1727 ivec dd_nc, ivec dd_nc_f)
1730 gmx_bool mm = FALSE;
1731 gmx_bool patchlevel_differs = FALSE;
1732 gmx_bool version_differs = FALSE;
1734 check_string(fplog, "Version", gmx_version(), version, &mm);
1735 patchlevel_differs = mm;
1737 if (patchlevel_differs)
1739 /* Gromacs should be able to continue from checkpoints between
1740 * different patch level versions, but we do not guarantee
1741 * compatibility between different major/minor versions - check this.
1743 int gmx_major, gmx_minor;
1744 int cpt_major, cpt_minor;
1745 sscanf(gmx_version(), "VERSION %5d.%5d", &gmx_major, &gmx_minor);
1746 int ret = sscanf(version, "VERSION %5d.%5d", &cpt_major, &cpt_minor);
1747 version_differs = (ret < 2 || gmx_major != cpt_major ||
1748 gmx_minor != cpt_minor);
1751 check_string(fplog, "Build time", BUILD_TIME, btime, &mm);
1752 check_string(fplog, "Build user", BUILD_USER, buser, &mm);
1753 check_string(fplog, "Build host", BUILD_HOST, bhost, &mm);
1754 check_int (fplog, "Double prec.", GMX_CPT_BUILD_DP, double_prec, &mm);
1755 check_string(fplog, "Program name", Program(), fprog, &mm);
1757 check_int (fplog, "#ranks", cr->nnodes, npp_f+npme_f, &mm);
1760 check_int (fplog, "#PME-ranks", cr->npmenodes, npme_f, &mm);
1763 if (cr->npmenodes >= 0)
1765 npp -= cr->npmenodes;
1769 check_int (fplog, "#DD-cells[x]", dd_nc[XX], dd_nc_f[XX], &mm);
1770 check_int (fplog, "#DD-cells[y]", dd_nc[YY], dd_nc_f[YY], &mm);
1771 check_int (fplog, "#DD-cells[z]", dd_nc[ZZ], dd_nc_f[ZZ], &mm);
1777 const char msg_version_difference[] =
1778 "The current GROMACS major & minor version are not identical to those that\n"
1779 "generated the checkpoint file. In principle GROMACS does not support\n"
1780 "continuation from checkpoints between different versions, so we advise\n"
1781 "against this. If you still want to try your luck we recommend that you use\n"
1782 "the -noappend flag to keep your output files from the two versions separate.\n"
1783 "This might also work around errors where the output fields in the energy\n"
1784 "file have changed between the different major & minor versions.\n";
1786 const char msg_mismatch_notice[] =
1787 "GROMACS patchlevel, binary or parallel settings differ from previous run.\n"
1788 "Continuation is exact, but not guaranteed to be binary identical.\n";
1790 const char msg_logdetails[] =
1791 "See the log file for details.\n";
1793 if (version_differs)
1795 fprintf(stderr, "%s%s\n", msg_version_difference, fplog ? msg_logdetails : "");
1799 fprintf(fplog, "%s\n", msg_version_difference);
1804 /* Major & minor versions match at least, but something is different. */
1805 fprintf(stderr, "%s%s\n", msg_mismatch_notice, fplog ? msg_logdetails : "");
1808 fprintf(fplog, "%s\n", msg_mismatch_notice);
1814 static void read_checkpoint(const char *fn, FILE **pfplog,
1815 t_commrec *cr, ivec dd_nc,
1816 int eIntegrator, int *init_fep_state, gmx_int64_t *step, double *t,
1817 t_state *state, gmx_bool *bReadEkin,
1818 int *simulation_part,
1819 gmx_bool bAppendOutputFiles, gmx_bool bForceAppend)
1824 char *version, *btime, *buser, *bhost, *fprog, *ftime;
1826 char buf[STEPSTRSIZE];
1827 int eIntegrator_f, nppnodes_f, npmenodes_f;
1829 int natoms, ngtc, nnhpres, nhchainlength, nlambda, fflags, flags_eks, flags_enh, flags_dfh;
1832 gmx_file_position_t *outputfiles;
1834 t_fileio *chksum_file;
1835 FILE * fplog = *pfplog;
1836 unsigned char digest[16];
1837 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1838 struct flock fl; /* don't initialize here: the struct order is OS
1842 const char *int_warn =
1843 "WARNING: The checkpoint file was generated with integrator %s,\n"
1844 " while the simulation uses integrator %s\n\n";
1846 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1847 fl.l_type = F_WRLCK;
1848 fl.l_whence = SEEK_SET;
1854 fp = gmx_fio_open(fn, "r");
1855 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
1856 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1857 &eIntegrator_f, simulation_part, step, t,
1858 &nppnodes_f, dd_nc_f, &npmenodes_f,
1859 &natoms, &ngtc, &nnhpres, &nhchainlength, &nlambda,
1860 &fflags, &flags_eks, &flags_enh, &flags_dfh,
1861 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
1863 if (bAppendOutputFiles &&
1864 file_version >= 13 && double_prec != GMX_CPT_BUILD_DP)
1866 gmx_fatal(FARGS, "Output file appending requested, but the code and checkpoint file precision (single/double) don't match");
1869 if (cr == NULL || MASTER(cr))
1871 fprintf(stderr, "\nReading checkpoint file %s generated: %s\n\n",
1875 /* This will not be written if we do appending, since fplog is still NULL then */
1878 fprintf(fplog, "\n");
1879 fprintf(fplog, "Reading checkpoint file %s\n", fn);
1880 fprintf(fplog, " file generated by: %s\n", fprog);
1881 fprintf(fplog, " file generated at: %s\n", ftime);
1882 fprintf(fplog, " GROMACS build time: %s\n", btime);
1883 fprintf(fplog, " GROMACS build user: %s\n", buser);
1884 fprintf(fplog, " GROMACS build host: %s\n", bhost);
1885 fprintf(fplog, " GROMACS double prec.: %d\n", double_prec);
1886 fprintf(fplog, " simulation part #: %d\n", *simulation_part);
1887 fprintf(fplog, " step: %s\n", gmx_step_str(*step, buf));
1888 fprintf(fplog, " time: %f\n", *t);
1889 fprintf(fplog, "\n");
1892 if (natoms != state->natoms)
1894 gmx_fatal(FARGS, "Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms", natoms, state->natoms);
1896 if (ngtc != state->ngtc)
1898 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);
1900 if (nnhpres != state->nnhpres)
1902 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);
1905 if (nlambda != state->dfhist.nlambda)
1907 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);
1910 init_gtc_state(state, state->ngtc, state->nnhpres, nhchainlength); /* need to keep this here to keep the tpr format working */
1911 /* write over whatever was read; we use the number of Nose-Hoover chains from the checkpoint */
1913 if (eIntegrator_f != eIntegrator)
1917 fprintf(stderr, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1919 if (bAppendOutputFiles)
1922 "Output file appending requested, but input/checkpoint integrators do not match.\n"
1923 "Stopping the run to prevent you from ruining all your data...\n"
1924 "If you _really_ know what you are doing, try with the -noappend option.\n");
1928 fprintf(fplog, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1936 else if (cr->nnodes == nppnodes_f + npmenodes_f)
1938 if (cr->npmenodes < 0)
1940 cr->npmenodes = npmenodes_f;
1942 int nppnodes = cr->nnodes - cr->npmenodes;
1943 if (nppnodes == nppnodes_f)
1945 for (d = 0; d < DIM; d++)
1949 dd_nc[d] = dd_nc_f[d];
1955 if (fflags != state->flags)
1960 if (bAppendOutputFiles)
1963 "Output file appending requested, but input and checkpoint states are not identical.\n"
1964 "Stopping the run to prevent you from ruining all your data...\n"
1965 "You can try with the -noappend option, and get more info in the log file.\n");
1968 if (getenv("GMX_ALLOW_CPT_MISMATCH") == NULL)
1970 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");
1975 "WARNING: The checkpoint state entries do not match the simulation,\n"
1976 " see the log file for details\n\n");
1982 print_flag_mismatch(fplog, state->flags, fflags);
1989 check_match(fplog, version, btime, buser, bhost, double_prec, fprog,
1990 cr, nppnodes_f, npmenodes_f, dd_nc, dd_nc_f);
1993 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, fflags, state, NULL);
1994 *init_fep_state = state->fep_state; /* there should be a better way to do this than setting it here.
1995 Investigate for 5.0. */
2000 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2005 *bReadEkin = ((flags_eks & (1<<eeksEKINH)) || (flags_eks & (1<<eeksEKINF)) || (flags_eks & (1<<eeksEKINO)) ||
2006 ((flags_eks & (1<<eeksEKINSCALEF)) | (flags_eks & (1<<eeksEKINSCALEH)) | (flags_eks & (1<<eeksVSCALE))));
2008 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2009 flags_enh, &state->enerhist, NULL);
2015 if (file_version < 6)
2017 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.";
2019 fprintf(stderr, "\nWARNING: %s\n\n", warn);
2022 fprintf(fplog, "\nWARNING: %s\n\n", warn);
2024 state->enerhist.nsum = *step;
2025 state->enerhist.nsum_sim = *step;
2028 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2034 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2040 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2046 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, NULL, file_version);
2052 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2057 if (gmx_fio_close(fp) != 0)
2059 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2068 /* If the user wants to append to output files,
2069 * we use the file pointer positions of the output files stored
2070 * in the checkpoint file and truncate the files such that any frames
2071 * written after the checkpoint time are removed.
2072 * All files are md5sum checked such that we can be sure that
2073 * we do not truncate other (maybe imprortant) files.
2075 if (bAppendOutputFiles)
2077 if (fn2ftp(outputfiles[0].filename) != efLOG)
2079 /* make sure first file is log file so that it is OK to use it for
2082 gmx_fatal(FARGS, "The first output file should always be the log "
2083 "file but instead is: %s. Cannot do appending because of this condition.", outputfiles[0].filename);
2085 for (i = 0; i < nfiles; i++)
2087 if (outputfiles[i].offset < 0)
2089 gmx_fatal(FARGS, "The original run wrote a file called '%s' which "
2090 "is larger than 2 GB, but mdrun did not support large file"
2091 " offsets. Can not append. Run mdrun with -noappend",
2092 outputfiles[i].filename);
2095 chksum_file = gmx_fio_open(outputfiles[i].filename, "a");
2098 chksum_file = gmx_fio_open(outputfiles[i].filename, "r+");
2103 /* Note that there are systems where the lock operation
2104 * will succeed, but a second process can also lock the file.
2105 * We should probably try to detect this.
2107 #if defined __native_client__
2111 #elif defined GMX_NATIVE_WINDOWS
2112 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX) == -1)
2114 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl) == -1)
2117 if (errno == ENOSYS)
2121 gmx_fatal(FARGS, "File locking is not supported on this system. Use -noappend or specify -append explicitly to append anyhow.");
2125 fprintf(stderr, "\nNOTE: File locking is not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2128 fprintf(fplog, "\nNOTE: File locking not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2132 else if (errno == EACCES || errno == EAGAIN)
2134 gmx_fatal(FARGS, "Failed to lock: %s. Already running "
2135 "simulation?", outputfiles[i].filename);
2139 gmx_fatal(FARGS, "Failed to lock: %s. %s.",
2140 outputfiles[i].filename, std::strerror(errno));
2145 /* compute md5 chksum */
2146 if (outputfiles[i].chksum_size != -1)
2148 if (gmx_fio_get_file_md5(chksum_file, outputfiles[i].offset,
2149 digest) != outputfiles[i].chksum_size) /*at the end of the call the file position is at the end of the file*/
2151 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.",
2152 outputfiles[i].chksum_size,
2153 outputfiles[i].filename);
2156 if (i == 0) /*log file needs to be seeked in case we need to truncate (other files are truncated below)*/
2158 if (gmx_fio_seek(chksum_file, outputfiles[i].offset))
2160 gmx_fatal(FARGS, "Seek error! Failed to truncate log-file: %s.", std::strerror(errno));
2165 if (i == 0) /*open log file here - so that lock is never lifted
2166 after chksum is calculated */
2168 *pfplog = gmx_fio_getfp(chksum_file);
2172 gmx_fio_close(chksum_file);
2175 /* compare md5 chksum */
2176 if (outputfiles[i].chksum_size != -1 &&
2177 memcmp(digest, outputfiles[i].chksum, 16) != 0)
2181 fprintf(debug, "chksum for %s: ", outputfiles[i].filename);
2182 for (j = 0; j < 16; j++)
2184 fprintf(debug, "%02x", digest[j]);
2186 fprintf(debug, "\n");
2188 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.",
2189 outputfiles[i].filename);
2194 if (i != 0) /*log file is already seeked to correct position */
2196 #if !defined(GMX_NATIVE_WINDOWS) || !defined(GMX_FAHCORE)
2197 /* For FAHCORE, we do this elsewhere*/
2198 rc = gmx_truncate(outputfiles[i].filename, outputfiles[i].offset);
2201 gmx_fatal(FARGS, "Truncation of file %s failed. Cannot do appending because of this failure.", outputfiles[i].filename);
2212 void load_checkpoint(const char *fn, FILE **fplog,
2213 t_commrec *cr, ivec dd_nc,
2214 t_inputrec *ir, t_state *state,
2215 gmx_bool *bReadEkin,
2216 gmx_bool bAppend, gmx_bool bForceAppend)
2223 /* Read the state from the checkpoint file */
2224 read_checkpoint(fn, fplog,
2226 ir->eI, &(ir->fepvals->init_fep_state), &step, &t, state, bReadEkin,
2227 &ir->simulation_part, bAppend, bForceAppend);
2231 gmx_bcast(sizeof(cr->npmenodes), &cr->npmenodes, cr);
2232 gmx_bcast(DIM*sizeof(dd_nc[0]), dd_nc, cr);
2233 gmx_bcast(sizeof(step), &step, cr);
2234 gmx_bcast(sizeof(*bReadEkin), bReadEkin, cr);
2236 ir->bContinuation = TRUE;
2237 if (ir->nsteps >= 0)
2239 ir->nsteps += ir->init_step - step;
2241 ir->init_step = step;
2242 ir->simulation_part += 1;
2245 void read_checkpoint_part_and_step(const char *filename,
2246 int *simulation_part,
2250 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2255 int flags_eks, flags_enh, flags_dfh;
2260 if (filename == NULL ||
2261 !gmx_fexist(filename) ||
2262 (!(fp = gmx_fio_open(filename, "r"))))
2264 *simulation_part = 0;
2269 /* Not calling initializing state before use is nasty, but all we
2270 do is read into its member variables and throw the struct away
2271 again immediately. */
2273 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2274 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2275 &eIntegrator, simulation_part, step, &t, &nppnodes, dd_nc, &npme,
2276 &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
2277 &(state.dfhist.nlambda), &state.flags, &flags_eks, &flags_enh, &flags_dfh,
2278 &state.edsamstate.nED, &state.swapstate.eSwapCoords, NULL);
2283 static void read_checkpoint_data(t_fileio *fp, int *simulation_part,
2284 gmx_int64_t *step, double *t, t_state *state,
2285 int *nfiles, gmx_file_position_t **outputfiles)
2288 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2293 int flags_eks, flags_enh, flags_dfh;
2295 gmx_file_position_t *files_loc = NULL;
2298 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2299 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2300 &eIntegrator, simulation_part, step, t, &nppnodes, dd_nc, &npme,
2301 &state->natoms, &state->ngtc, &state->nnhpres, &state->nhchainlength,
2302 &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
2303 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
2305 do_cpt_state(gmx_fio_getxdr(fp), TRUE, state->flags, state, NULL);
2310 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2315 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2316 flags_enh, &state->enerhist, NULL);
2321 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2327 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2333 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2339 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE,
2340 outputfiles != NULL ? outputfiles : &files_loc,
2341 outputfiles != NULL ? nfiles : &nfiles_loc,
2342 NULL, file_version);
2343 if (files_loc != NULL)
2353 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2367 read_checkpoint_state(const char *fn, int *simulation_part,
2368 gmx_int64_t *step, double *t, t_state *state)
2372 fp = gmx_fio_open(fn, "r");
2373 read_checkpoint_data(fp, simulation_part, step, t, state, NULL, NULL);
2374 if (gmx_fio_close(fp) != 0)
2376 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2380 void read_checkpoint_trxframe(t_fileio *fp, t_trxframe *fr)
2382 /* This next line is nasty because the sub-structures of t_state
2383 * cannot be assumed to be zeroed (or even initialized in ways the
2384 * rest of the code might assume). Using snew would be better, but
2385 * this will all go away for 5.0. */
2387 int simulation_part;
2391 init_state(&state, 0, 0, 0, 0, 0);
2393 read_checkpoint_data(fp, &simulation_part, &step, &t, &state, NULL, NULL);
2395 fr->natoms = state.natoms;
2398 fr->step = gmx_int64_to_int(step,
2399 "conversion of checkpoint to trajectory");
2403 fr->lambda = state.lambda[efptFEP];
2404 fr->fep_state = state.fep_state;
2406 fr->bX = (state.flags & (1<<estX));
2412 fr->bV = (state.flags & (1<<estV));
2419 fr->bBox = (state.flags & (1<<estBOX));
2422 copy_mat(state.box, fr->box);
2427 void list_checkpoint(const char *fn, FILE *out)
2431 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2433 int eIntegrator, simulation_part, nppnodes, npme;
2438 int flags_eks, flags_enh, flags_dfh;
2440 gmx_file_position_t *outputfiles;
2443 init_state(&state, -1, -1, -1, -1, 0);
2445 fp = gmx_fio_open(fn, "r");
2446 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2447 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2448 &eIntegrator, &simulation_part, &step, &t, &nppnodes, dd_nc, &npme,
2449 &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
2450 &(state.dfhist.nlambda), &state.flags,
2451 &flags_eks, &flags_enh, &flags_dfh, &state.edsamstate.nED,
2452 &state.swapstate.eSwapCoords, out);
2453 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, state.flags, &state, out);
2458 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state.ekinstate, out);
2463 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2464 flags_enh, &state.enerhist, out);
2468 ret = do_cpt_df_hist(gmx_fio_getxdr(fp),
2469 flags_dfh, &state.dfhist, out);
2474 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state.edsamstate, out);
2479 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state.swapstate, out);
2484 do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, out, file_version);
2489 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2496 if (gmx_fio_close(fp) != 0)
2498 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2504 /* This routine cannot print tons of data, since it is called before the log file is opened. */
2506 read_checkpoint_simulation_part_and_filenames(t_fileio *fp,
2507 int *simulation_part,
2509 gmx_file_position_t **outputfiles)
2511 gmx_int64_t step = 0;
2515 init_state(&state, 0, 0, 0, 0, 0);
2517 read_checkpoint_data(fp, simulation_part, &step, &t, &state,
2518 nfiles, outputfiles);
2519 if (gmx_fio_close(fp) != 0)
2521 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");