2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2008,2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 /* The source code in this file should be thread-safe.
37 Please keep it that way. */
38 #include "checkpoint.h"
50 #ifdef HAVE_SYS_TIME_H
58 #ifdef GMX_NATIVE_WINDOWS
61 #include <sys/locking.h>
67 #include "types/commrec.h"
69 #include "gromacs/math/vec.h"
72 #include "gromacs/fileio/filenm.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/fileio/gmxfio.h"
75 #include "gromacs/fileio/xdrf.h"
76 #include "gromacs/fileio/xdr_datatype.h"
77 #include "gromacs/utility/basenetwork.h"
78 #include "gromacs/utility/baseversion.h"
79 #include "gromacs/utility/cstringutil.h"
80 #include "gromacs/utility/fatalerror.h"
81 #include "gromacs/utility/smalloc.h"
83 #include "buildinfo.h"
89 #define CPT_MAGIC1 171817
90 #define CPT_MAGIC2 171819
91 #define CPTSTRLEN 1024
94 #define GMX_CPT_BUILD_DP 1
96 #define GMX_CPT_BUILD_DP 0
99 /* cpt_version should normally only be changed
100 * when the header of footer format changes.
101 * The state data format itself is backward and forward compatible.
102 * But old code can not read a new entry that is present in the file
103 * (but can read a new format when new entries are not present).
105 static const int cpt_version = 16;
108 const char *est_names[estNR] =
111 "box", "box-rel", "box-v", "pres_prev",
112 "nosehoover-xi", "thermostat-integral",
113 "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
114 "disre_initf", "disre_rm3tav",
115 "orire_initf", "orire_Dtav",
116 "svir_prev", "nosehoover-vxi", "v_eta", "vol0", "nhpres_xi", "nhpres_vxi", "fvir_prev", "fep_state", "MC-rng", "MC-rng-i"
120 eeksEKIN_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksEKINF, eeksEKINO, eeksEKINSCALEF, eeksEKINSCALEH, eeksVSCALE, eeksEKINTOTAL, eeksNR
123 const char *eeks_names[eeksNR] =
125 "Ekin_n", "Ekinh", "dEkindlambda", "mv_cos",
126 "Ekinf", "Ekinh_old", "EkinScaleF_NHC", "EkinScaleH_NHC", "Vscale_NHC", "Ekin_Total"
130 eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
131 eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
132 eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM,
133 eenhENERGY_DELTA_H_NN,
134 eenhENERGY_DELTA_H_LIST,
135 eenhENERGY_DELTA_H_STARTTIME,
136 eenhENERGY_DELTA_H_STARTLAMBDA,
140 const char *eenh_names[eenhNR] =
142 "energy_n", "energy_aver", "energy_sum", "energy_nsum",
143 "energy_sum_sim", "energy_nsum_sim",
144 "energy_nsteps", "energy_nsteps_sim",
146 "energy_delta_h_list",
147 "energy_delta_h_start_time",
148 "energy_delta_h_start_lambda"
151 /* free energy history variables -- need to be preserved over checkpoint */
153 edfhBEQUIL, edfhNATLAMBDA, edfhWLHISTO, edfhWLDELTA, edfhSUMWEIGHTS, edfhSUMDG, edfhSUMMINVAR, edfhSUMVAR,
154 edfhACCUMP, edfhACCUMM, edfhACCUMP2, edfhACCUMM2, edfhTIJ, edfhTIJEMP, edfhNR
156 /* free energy history variable names */
157 const char *edfh_names[edfhNR] =
159 "bEquilibrated", "N_at_state", "Wang-Landau Histogram", "Wang-Landau Delta", "Weights", "Free Energies", "minvar", "variance",
160 "accumulated_plus", "accumulated_minus", "accumulated_plus_2", "accumulated_minus_2", "Tij", "Tij_empirical"
163 #ifdef GMX_NATIVE_WINDOWS
165 gmx_wintruncate(const char *filename, __int64 size)
168 /*we do this elsewhere*/
173 fp = fopen(filename, "rb+");
180 return _chsize_s( fileno(fp), size);
187 ecprREAL, ecprRVEC, ecprMATRIX
191 cptpEST, cptpEEKS, cptpEENH, cptpEDFH
193 /* enums for the different components of checkpoint variables, replacing the hard coded ones.
194 cptpEST - state variables.
195 cptpEEKS - Kinetic energy state variables.
196 cptpEENH - Energy history state variables.
197 cptpEDFH - free energy history variables.
201 static const char *st_names(int cptp, int ecpt)
205 case cptpEST: return est_names [ecpt]; break;
206 case cptpEEKS: return eeks_names[ecpt]; break;
207 case cptpEENH: return eenh_names[ecpt]; break;
208 case cptpEDFH: return edfh_names[ecpt]; break;
214 static void cp_warning(FILE *fp)
216 fprintf(fp, "\nWARNING: Checkpoint file is corrupted or truncated\n\n");
219 static void cp_error()
221 gmx_fatal(FARGS, "Checkpoint file corrupted/truncated, or maybe you are out of disk space?");
224 static void do_cpt_string_err(XDR *xd, gmx_bool bRead, const char *desc, char **s, FILE *list)
232 res = xdr_string(xd, s, CPTSTRLEN);
239 fprintf(list, "%s = %s\n", desc, *s);
244 static int do_cpt_int(XDR *xd, const char *desc, int *i, FILE *list)
248 res = xdr_int(xd, i);
255 fprintf(list, "%s = %d\n", desc, *i);
260 static int do_cpt_u_chars(XDR *xd, const char *desc, int n, unsigned char *i, FILE *list)
266 fprintf(list, "%s = ", desc);
268 for (j = 0; j < n && res; j++)
270 res &= xdr_u_char(xd, &i[j]);
273 fprintf(list, "%02x", i[j]);
288 static void do_cpt_int_err(XDR *xd, const char *desc, int *i, FILE *list)
290 if (do_cpt_int(xd, desc, i, list) < 0)
296 static void do_cpt_step_err(XDR *xd, const char *desc, gmx_int64_t *i, FILE *list)
299 char buf[STEPSTRSIZE];
301 res = xdr_int64(xd, i);
308 fprintf(list, "%s = %s\n", desc, gmx_step_str(*i, buf));
312 static void do_cpt_double_err(XDR *xd, const char *desc, double *f, FILE *list)
316 res = xdr_double(xd, f);
323 fprintf(list, "%s = %f\n", desc, *f);
327 static void do_cpt_real_err(XDR *xd, real *f)
332 res = xdr_double(xd, f);
334 res = xdr_float(xd, f);
342 static void do_cpt_n_rvecs_err(XDR *xd, const char *desc, int n, rvec f[], FILE *list)
346 for (i = 0; i < n; i++)
348 for (j = 0; j < DIM; j++)
350 do_cpt_real_err(xd, &f[i][j]);
356 pr_rvecs(list, 0, desc, f, n);
360 /* If nval >= 0, nval is used; on read this should match the passed value.
361 * If nval n<0, *nptr is used; on read the value is stored in nptr
363 static int do_cpte_reals_low(XDR *xd, int cptp, int ecpt, int sflags,
364 int nval, int *nptr, real **v,
365 FILE *list, int erealtype)
369 int dtc = xdr_datatype_float;
371 int dtc = xdr_datatype_double;
373 real *vp, *va = NULL;
388 gmx_incons("*ntpr=NULL in do_cpte_reals_low");
393 res = xdr_int(xd, &nf);
404 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), nval, nf);
413 res = xdr_int(xd, &dt);
420 fprintf(stderr, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
421 st_names(cptp, ecpt), xdr_datatype_names[dtc],
422 xdr_datatype_names[dt]);
424 if (list || !(sflags & (1<<ecpt)))
437 if (dt == xdr_datatype_float)
439 if (dtc == xdr_datatype_float)
447 res = xdr_vector(xd, (char *)vf, nf,
448 (unsigned int)sizeof(float), (xdrproc_t)xdr_float);
453 if (dtc != xdr_datatype_float)
455 for (i = 0; i < nf; i++)
464 if (dtc == xdr_datatype_double)
466 /* cppcheck-suppress invalidPointerCast
467 * Only executed if real is anyhow double */
474 res = xdr_vector(xd, (char *)vd, nf,
475 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
480 if (dtc != xdr_datatype_double)
482 for (i = 0; i < nf; i++)
495 pr_reals(list, 0, st_names(cptp, ecpt), vp, nf);
498 pr_rvecs(list, 0, st_names(cptp, ecpt), (rvec *)vp, nf/3);
501 gmx_incons("Unknown checkpoint real type");
513 /* This function stores n along with the reals for reading,
514 * but on reading it assumes that n matches the value in the checkpoint file,
515 * a fatal error is generated when this is not the case.
517 static int do_cpte_reals(XDR *xd, int cptp, int ecpt, int sflags,
518 int n, real **v, FILE *list)
520 return do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, v, list, ecprREAL);
523 /* This function does the same as do_cpte_reals,
524 * except that on reading it ignores the passed value of *n
525 * and stored the value read from the checkpoint file in *n.
527 static int do_cpte_n_reals(XDR *xd, int cptp, int ecpt, int sflags,
528 int *n, real **v, FILE *list)
530 return do_cpte_reals_low(xd, cptp, ecpt, sflags, -1, n, v, list, ecprREAL);
533 static int do_cpte_real(XDR *xd, int cptp, int ecpt, int sflags,
536 return do_cpte_reals_low(xd, cptp, ecpt, sflags, 1, NULL, &r, list, ecprREAL);
539 static int do_cpte_ints(XDR *xd, int cptp, int ecpt, int sflags,
540 int n, int **v, FILE *list)
543 int dtc = xdr_datatype_int;
548 res = xdr_int(xd, &nf);
553 if (list == NULL && v != NULL && nf != n)
555 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
558 res = xdr_int(xd, &dt);
565 gmx_fatal(FARGS, "Type mismatch for state entry %s, code type is %s, file type is %s\n",
566 st_names(cptp, ecpt), xdr_datatype_names[dtc],
567 xdr_datatype_names[dt]);
569 if (list || !(sflags & (1<<ecpt)) || v == NULL)
582 res = xdr_vector(xd, (char *)vp, nf,
583 (unsigned int)sizeof(int), (xdrproc_t)xdr_int);
590 pr_ivec(list, 0, st_names(cptp, ecpt), vp, nf, TRUE);
600 static int do_cpte_int(XDR *xd, int cptp, int ecpt, int sflags,
603 return do_cpte_ints(xd, cptp, ecpt, sflags, 1, &i, list);
606 static int do_cpte_doubles(XDR *xd, int cptp, int ecpt, int sflags,
607 int n, double **v, FILE *list)
610 int dtc = xdr_datatype_double;
611 double *vp, *va = NULL;
615 res = xdr_int(xd, &nf);
620 if (list == NULL && nf != n)
622 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
625 res = xdr_int(xd, &dt);
632 gmx_fatal(FARGS, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
633 st_names(cptp, ecpt), xdr_datatype_names[dtc],
634 xdr_datatype_names[dt]);
636 if (list || !(sflags & (1<<ecpt)))
649 res = xdr_vector(xd, (char *)vp, nf,
650 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
657 pr_doubles(list, 0, st_names(cptp, ecpt), vp, nf);
667 static int do_cpte_double(XDR *xd, int cptp, int ecpt, int sflags,
668 double *r, FILE *list)
670 return do_cpte_doubles(xd, cptp, ecpt, sflags, 1, &r, list);
674 static int do_cpte_rvecs(XDR *xd, int cptp, int ecpt, int sflags,
675 int n, rvec **v, FILE *list)
677 return do_cpte_reals_low(xd, cptp, ecpt, sflags,
678 n*DIM, NULL, (real **)v, list, ecprRVEC);
681 static int do_cpte_matrix(XDR *xd, int cptp, int ecpt, int sflags,
682 matrix v, FILE *list)
687 vr = (real *)&(v[0][0]);
688 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
689 DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
691 if (list && ret == 0)
693 pr_rvecs(list, 0, st_names(cptp, ecpt), v, DIM);
700 static int do_cpte_nmatrix(XDR *xd, int cptp, int ecpt, int sflags,
701 int n, real **v, FILE *list)
705 char name[CPTSTRLEN];
712 for (i = 0; i < n; i++)
714 reti = do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, &(v[i]), NULL, ecprREAL);
715 if (list && reti == 0)
717 sprintf(name, "%s[%d]", st_names(cptp, ecpt), i);
718 pr_reals(list, 0, name, v[i], n);
728 static int do_cpte_matrices(XDR *xd, int cptp, int ecpt, int sflags,
729 int n, matrix **v, FILE *list)
732 matrix *vp, *va = NULL;
738 res = xdr_int(xd, &nf);
743 if (list == NULL && nf != n)
745 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
747 if (list || !(sflags & (1<<ecpt)))
760 snew(vr, nf*DIM*DIM);
761 for (i = 0; i < nf; i++)
763 for (j = 0; j < DIM; j++)
765 for (k = 0; k < DIM; k++)
767 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
771 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
772 nf*DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
773 for (i = 0; i < nf; i++)
775 for (j = 0; j < DIM; j++)
777 for (k = 0; k < DIM; k++)
779 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
785 if (list && ret == 0)
787 for (i = 0; i < nf; i++)
789 pr_rvecs(list, 0, st_names(cptp, ecpt), vp[i], DIM);
800 static void do_cpt_header(XDR *xd, gmx_bool bRead, int *file_version,
801 char **version, char **btime, char **buser, char **bhost,
803 char **fprog, char **ftime,
804 int *eIntegrator, int *simulation_part,
805 gmx_int64_t *step, double *t,
806 int *nnodes, int *dd_nc, int *npme,
807 int *natoms, int *ngtc, int *nnhpres, int *nhchainlength,
808 int *nlambda, int *flags_state,
809 int *flags_eks, int *flags_enh, int *flags_dfh,
810 int *nED, int *eSwapCoords,
826 res = xdr_int(xd, &magic);
829 gmx_fatal(FARGS, "The checkpoint file is empty/corrupted, or maybe you are out of disk space?");
831 if (magic != CPT_MAGIC1)
833 gmx_fatal(FARGS, "Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
834 "The checkpoint file is corrupted or not a checkpoint file",
840 gmx_gethostname(fhost, 255);
842 do_cpt_string_err(xd, bRead, "GROMACS version", version, list);
843 do_cpt_string_err(xd, bRead, "GROMACS build time", btime, list);
844 do_cpt_string_err(xd, bRead, "GROMACS build user", buser, list);
845 do_cpt_string_err(xd, bRead, "GROMACS build host", bhost, list);
846 do_cpt_string_err(xd, bRead, "generating program", fprog, list);
847 do_cpt_string_err(xd, bRead, "generation time", ftime, list);
848 *file_version = cpt_version;
849 do_cpt_int_err(xd, "checkpoint file version", file_version, list);
850 if (*file_version > cpt_version)
852 gmx_fatal(FARGS, "Attempting to read a checkpoint file of version %d with code of version %d\n", *file_version, cpt_version);
854 if (*file_version >= 13)
856 do_cpt_int_err(xd, "GROMACS double precision", double_prec, list);
862 if (*file_version >= 12)
864 do_cpt_string_err(xd, bRead, "generating host", &fhost, list);
870 do_cpt_int_err(xd, "#atoms", natoms, list);
871 do_cpt_int_err(xd, "#T-coupling groups", ngtc, list);
872 if (*file_version >= 10)
874 do_cpt_int_err(xd, "#Nose-Hoover T-chains", nhchainlength, list);
880 if (*file_version >= 11)
882 do_cpt_int_err(xd, "#Nose-Hoover T-chains for barostat ", nnhpres, list);
888 if (*file_version >= 14)
890 do_cpt_int_err(xd, "# of total lambda states ", nlambda, list);
896 do_cpt_int_err(xd, "integrator", eIntegrator, list);
897 if (*file_version >= 3)
899 do_cpt_int_err(xd, "simulation part #", simulation_part, list);
903 *simulation_part = 1;
905 if (*file_version >= 5)
907 do_cpt_step_err(xd, "step", step, list);
911 do_cpt_int_err(xd, "step", &idum, list);
914 do_cpt_double_err(xd, "t", t, list);
915 do_cpt_int_err(xd, "#PP-nodes", nnodes, list);
917 do_cpt_int_err(xd, "dd_nc[x]", dd_nc ? &(dd_nc[0]) : &idum, list);
918 do_cpt_int_err(xd, "dd_nc[y]", dd_nc ? &(dd_nc[1]) : &idum, list);
919 do_cpt_int_err(xd, "dd_nc[z]", dd_nc ? &(dd_nc[2]) : &idum, list);
920 do_cpt_int_err(xd, "#PME-only nodes", npme, list);
921 do_cpt_int_err(xd, "state flags", flags_state, list);
922 if (*file_version >= 4)
924 do_cpt_int_err(xd, "ekin data flags", flags_eks, list);
925 do_cpt_int_err(xd, "energy history flags", flags_enh, list);
930 *flags_enh = (*flags_state >> (estORIRE_DTAV+1));
931 *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
932 (1<<(estORIRE_DTAV+2)) |
933 (1<<(estORIRE_DTAV+3))));
935 if (*file_version >= 14)
937 do_cpt_int_err(xd, "df history flags", flags_dfh, list);
944 if (*file_version >= 15)
946 do_cpt_int_err(xd, "ED data sets", nED, list);
952 if (*file_version >= 16)
954 do_cpt_int_err(xd, "swap", eSwapCoords, list);
958 static int do_cpt_footer(XDR *xd, int file_version)
963 if (file_version >= 2)
966 res = xdr_int(xd, &magic);
971 if (magic != CPT_MAGIC2)
980 static int do_cpt_state(XDR *xd, gmx_bool bRead,
981 int fflags, t_state *state,
991 nnht = state->nhchainlength*state->ngtc;
992 nnhtp = state->nhchainlength*state->nnhpres;
994 if (bRead) /* we need to allocate space for dfhist if we are reading */
996 init_df_history(&state->dfhist, state->dfhist.nlambda);
999 sflags = state->flags;
1000 for (i = 0; (i < estNR && ret == 0); i++)
1002 if (fflags & (1<<i))
1006 case estLAMBDA: ret = do_cpte_reals(xd, cptpEST, i, sflags, efptNR, &(state->lambda), list); break;
1007 case estFEPSTATE: ret = do_cpte_int (xd, cptpEST, i, sflags, &state->fep_state, list); break;
1008 case estBOX: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box, list); break;
1009 case estBOX_REL: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box_rel, list); break;
1010 case estBOXV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->boxv, list); break;
1011 case estPRES_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->pres_prev, list); break;
1012 case estSVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->svir_prev, list); break;
1013 case estFVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->fvir_prev, list); break;
1014 case estNH_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_xi, list); break;
1015 case estNH_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_vxi, list); break;
1016 case estNHPRES_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_xi, list); break;
1017 case estNHPRES_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_vxi, list); break;
1018 case estTC_INT: ret = do_cpte_doubles(xd, cptpEST, i, sflags, state->ngtc, &state->therm_integral, list); break;
1019 case estVETA: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->veta, list); break;
1020 case estVOL0: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->vol0, list); break;
1021 case estX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->x, list); break;
1022 case estV: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->v, list); break;
1023 case estSDX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->sd_X, list); break;
1024 /* The RNG entries are no longer written,
1025 * the next 4 lines are only for reading old files.
1027 case estLD_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1028 case estLD_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1029 case estMC_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1030 case estMC_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1031 case estDISRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.disre_initf, list); break;
1032 case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.ndisrepairs, &state->hist.disre_rm3tav, list); break;
1033 case estORIRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.orire_initf, list); break;
1034 case estORIRE_DTAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.norire_Dtav, &state->hist.orire_Dtav, list); break;
1036 gmx_fatal(FARGS, "Unknown state entry %d\n"
1037 "You are probably reading a new checkpoint file with old code", i);
1045 static int do_cpt_ekinstate(XDR *xd, int fflags, ekinstate_t *ekins,
1053 for (i = 0; (i < eeksNR && ret == 0); i++)
1055 if (fflags & (1<<i))
1060 case eeksEKIN_N: ret = do_cpte_int(xd, cptpEEKS, i, fflags, &ekins->ekin_n, list); break;
1061 case eeksEKINH: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh, list); break;
1062 case eeksEKINF: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinf, list); break;
1063 case eeksEKINO: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh_old, list); break;
1064 case eeksEKINTOTAL: ret = do_cpte_matrix(xd, cptpEEKS, i, fflags, ekins->ekin_total, list); break;
1065 case eeksEKINSCALEF: ret = do_cpte_doubles(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinscalef_nhc, list); break;
1066 case eeksVSCALE: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->vscale_nhc, list); break;
1067 case eeksEKINSCALEH: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->ekinscaleh_nhc, list); break;
1068 case eeksDEKINDL: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->dekindl, list); break;
1069 case eeksMVCOS: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->mvcos, list); break;
1071 gmx_fatal(FARGS, "Unknown ekin data state entry %d\n"
1072 "You are probably reading a new checkpoint file with old code", i);
1081 static int do_cpt_swapstate(XDR *xd, gmx_bool bRead, swapstate_t *swapstate, FILE *list)
1085 int swap_cpt_version = 1;
1088 if (eswapNO == swapstate->eSwapCoords)
1093 swapstate->bFromCpt = bRead;
1095 do_cpt_int_err(xd, "swap checkpoint version", &swap_cpt_version, list);
1096 do_cpt_int_err(xd, "swap coupling steps", &swapstate->nAverage, list);
1098 /* When reading, init_swapcoords has not been called yet,
1099 * so we have to allocate memory first. */
1101 for (ic = 0; ic < eCompNR; ic++)
1103 for (ii = 0; ii < eIonNR; ii++)
1107 do_cpt_int_err(xd, "swap requested atoms", &swapstate->nat_req[ic][ii], list);
1111 do_cpt_int_err(xd, "swap requested atoms p", swapstate->nat_req_p[ic][ii], list);
1116 do_cpt_int_err(xd, "swap influx netto", &swapstate->inflow_netto[ic][ii], list);
1120 do_cpt_int_err(xd, "swap influx netto p", swapstate->inflow_netto_p[ic][ii], list);
1123 if (bRead && (NULL == swapstate->nat_past[ic][ii]) )
1125 snew(swapstate->nat_past[ic][ii], swapstate->nAverage);
1128 for (j = 0; j < swapstate->nAverage; j++)
1132 do_cpt_int_err(xd, "swap past atom counts", &swapstate->nat_past[ic][ii][j], list);
1136 do_cpt_int_err(xd, "swap past atom counts p", &swapstate->nat_past_p[ic][ii][j], list);
1142 /* Ion flux per channel */
1143 for (ic = 0; ic < eChanNR; ic++)
1145 for (ii = 0; ii < eIonNR; ii++)
1149 do_cpt_int_err(xd, "channel flux", &swapstate->fluxfromAtoB[ic][ii], list);
1153 do_cpt_int_err(xd, "channel flux p", swapstate->fluxfromAtoB_p[ic][ii], list);
1158 /* Ion flux leakage */
1161 snew(swapstate->fluxleak, 1);
1163 do_cpt_int_err(xd, "flux leakage", swapstate->fluxleak, list);
1166 do_cpt_int_err(xd, "number of ions", &swapstate->nions, list);
1170 snew(swapstate->channel_label, swapstate->nions);
1171 snew(swapstate->comp_from, swapstate->nions);
1174 do_cpt_u_chars(xd, "channel history", swapstate->nions, swapstate->channel_label, list);
1175 do_cpt_u_chars(xd, "domain history", swapstate->nions, swapstate->comp_from, list);
1177 /* Save the last known whole positions to checkpoint
1178 * file to be able to also make multimeric channels whole in PBC */
1179 do_cpt_int_err(xd, "Ch0 atoms", &swapstate->nat[eChan0], list);
1180 do_cpt_int_err(xd, "Ch1 atoms", &swapstate->nat[eChan1], list);
1183 snew(swapstate->xc_old_whole[eChan0], swapstate->nat[eChan0]);
1184 snew(swapstate->xc_old_whole[eChan1], swapstate->nat[eChan1]);
1185 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], swapstate->xc_old_whole[eChan0], list);
1186 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], swapstate->xc_old_whole[eChan1], list);
1190 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], *swapstate->xc_old_whole_p[eChan0], list);
1191 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], *swapstate->xc_old_whole_p[eChan1], list);
1198 static int do_cpt_enerhist(XDR *xd, gmx_bool bRead,
1199 int fflags, energyhistory_t *enerhist,
1210 enerhist->nsteps = 0;
1212 enerhist->nsteps_sim = 0;
1213 enerhist->nsum_sim = 0;
1214 enerhist->dht = NULL;
1216 if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
1218 snew(enerhist->dht, 1);
1219 enerhist->dht->ndh = NULL;
1220 enerhist->dht->dh = NULL;
1221 enerhist->dht->start_lambda_set = FALSE;
1225 for (i = 0; (i < eenhNR && ret == 0); i++)
1227 if (fflags & (1<<i))
1231 case eenhENERGY_N: ret = do_cpte_int(xd, cptpEENH, i, fflags, &enerhist->nener, list); break;
1232 case eenhENERGY_AVER: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_ave, list); break;
1233 case eenhENERGY_SUM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum, list); break;
1234 case eenhENERGY_NSUM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum, list); break;
1235 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum_sim, list); break;
1236 case eenhENERGY_NSUM_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum_sim, list); break;
1237 case eenhENERGY_NSTEPS: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps, list); break;
1238 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps_sim, list); break;
1239 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd, eenh_names[i], &(enerhist->dht->nndh), list);
1240 if (bRead) /* now allocate memory for it */
1242 snew(enerhist->dht->dh, enerhist->dht->nndh);
1243 snew(enerhist->dht->ndh, enerhist->dht->nndh);
1244 for (j = 0; j < enerhist->dht->nndh; j++)
1246 enerhist->dht->ndh[j] = 0;
1247 enerhist->dht->dh[j] = NULL;
1251 case eenhENERGY_DELTA_H_LIST:
1252 for (j = 0; j < enerhist->dht->nndh; j++)
1254 ret = do_cpte_n_reals(xd, cptpEENH, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list);
1257 case eenhENERGY_DELTA_H_STARTTIME:
1258 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_time), list); break;
1259 case eenhENERGY_DELTA_H_STARTLAMBDA:
1260 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_lambda), list); break;
1262 gmx_fatal(FARGS, "Unknown energy history entry %d\n"
1263 "You are probably reading a new checkpoint file with old code", i);
1268 if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
1270 /* Assume we have an old file format and copy sum to sum_sim */
1271 srenew(enerhist->ener_sum_sim, enerhist->nener);
1272 for (i = 0; i < enerhist->nener; i++)
1274 enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1278 if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1279 !(fflags & (1<<eenhENERGY_NSTEPS)))
1281 /* Assume we have an old file format and copy nsum to nsteps */
1282 enerhist->nsteps = enerhist->nsum;
1284 if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1285 !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1287 /* Assume we have an old file format and copy nsum to nsteps */
1288 enerhist->nsteps_sim = enerhist->nsum_sim;
1294 static int do_cpt_df_hist(XDR *xd, int fflags, df_history_t *dfhist, FILE *list)
1299 nlambda = dfhist->nlambda;
1302 for (i = 0; (i < edfhNR && ret == 0); i++)
1304 if (fflags & (1<<i))
1308 case edfhBEQUIL: ret = do_cpte_int(xd, cptpEDFH, i, fflags, &dfhist->bEquil, list); break;
1309 case edfhNATLAMBDA: ret = do_cpte_ints(xd, cptpEDFH, i, fflags, nlambda, &dfhist->n_at_lam, list); break;
1310 case edfhWLHISTO: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->wl_histo, list); break;
1311 case edfhWLDELTA: ret = do_cpte_real(xd, cptpEDFH, i, fflags, &dfhist->wl_delta, list); break;
1312 case edfhSUMWEIGHTS: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_weights, list); break;
1313 case edfhSUMDG: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_dg, list); break;
1314 case edfhSUMMINVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_minvar, list); break;
1315 case edfhSUMVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_variance, list); break;
1316 case edfhACCUMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p, list); break;
1317 case edfhACCUMM: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m, list); break;
1318 case edfhACCUMP2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p2, list); break;
1319 case edfhACCUMM2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m2, list); break;
1320 case edfhTIJ: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij, list); break;
1321 case edfhTIJEMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij_empirical, list); break;
1324 gmx_fatal(FARGS, "Unknown df history entry %d\n"
1325 "You are probably reading a new checkpoint file with old code", i);
1334 /* This function stores the last whole configuration of the reference and
1335 * average structure in the .cpt file
1337 static int do_cpt_EDstate(XDR *xd, gmx_bool bRead,
1338 edsamstate_t *EDstate, FILE *list)
1345 EDstate->bFromCpt = bRead;
1347 if (EDstate->nED <= 0)
1352 /* When reading, init_edsam has not been called yet,
1353 * so we have to allocate memory first. */
1356 snew(EDstate->nref, EDstate->nED);
1357 snew(EDstate->old_sref, EDstate->nED);
1358 snew(EDstate->nav, EDstate->nED);
1359 snew(EDstate->old_sav, EDstate->nED);
1362 /* Read/write the last whole conformation of SREF and SAV for each ED dataset (usually only one) */
1363 for (i = 0; i < EDstate->nED; i++)
1365 /* Reference structure SREF */
1366 sprintf(buf, "ED%d # of atoms in reference structure", i+1);
1367 do_cpt_int_err(xd, buf, &EDstate->nref[i], list);
1368 sprintf(buf, "ED%d x_ref", i+1);
1371 snew(EDstate->old_sref[i], EDstate->nref[i]);
1372 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref[i], list);
1376 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref_p[i], list);
1379 /* Average structure SAV */
1380 sprintf(buf, "ED%d # of atoms in average structure", i+1);
1381 do_cpt_int_err(xd, buf, &EDstate->nav[i], list);
1382 sprintf(buf, "ED%d x_av", i+1);
1385 snew(EDstate->old_sav[i], EDstate->nav[i]);
1386 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav[i], list);
1390 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav_p[i], list);
1398 static int do_cpt_files(XDR *xd, gmx_bool bRead,
1399 gmx_file_position_t **p_outputfiles, int *nfiles,
1400 FILE *list, int file_version)
1404 gmx_off_t mask = 0xFFFFFFFFL;
1405 int offset_high, offset_low;
1407 gmx_file_position_t *outputfiles;
1409 if (do_cpt_int(xd, "number of output files", nfiles, list) != 0)
1416 snew(*p_outputfiles, *nfiles);
1419 outputfiles = *p_outputfiles;
1421 for (i = 0; i < *nfiles; i++)
1423 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1426 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1427 strncpy(outputfiles[i].filename, buf, CPTSTRLEN-1);
1433 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1437 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1441 outputfiles[i].offset = ( ((gmx_off_t) offset_high) << 32 ) | ( (gmx_off_t) offset_low & mask );
1445 buf = outputfiles[i].filename;
1446 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1448 offset = outputfiles[i].offset;
1456 offset_low = (int) (offset & mask);
1457 offset_high = (int) ((offset >> 32) & mask);
1459 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1463 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1468 if (file_version >= 8)
1470 if (do_cpt_int(xd, "file_checksum_size", &(outputfiles[i].chksum_size),
1475 if (do_cpt_u_chars(xd, "file_checksum", 16, outputfiles[i].chksum, list) != 0)
1482 outputfiles[i].chksum_size = -1;
1489 void write_checkpoint(const char *fn, gmx_bool bNumberAndKeep,
1490 FILE *fplog, t_commrec *cr,
1491 int eIntegrator, int simulation_part,
1492 gmx_bool bExpanded, int elamstats,
1493 gmx_int64_t step, double t, t_state *state)
1503 char *fntemp; /* the temporary checkpoint file name */
1505 char timebuf[STRLEN];
1506 int nppnodes, npmenodes;
1507 char buf[1024], suffix[5+STEPSTRSIZE], sbuf[STEPSTRSIZE];
1508 gmx_file_position_t *outputfiles;
1511 int flags_eks, flags_enh, flags_dfh;
1514 if (DOMAINDECOMP(cr))
1516 nppnodes = cr->dd->nnodes;
1517 npmenodes = cr->npmenodes;
1525 #ifndef GMX_NO_RENAME
1526 /* make the new temporary filename */
1527 snew(fntemp, strlen(fn)+5+STEPSTRSIZE);
1529 fntemp[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1530 sprintf(suffix, "_%s%s", "step", gmx_step_str(step, sbuf));
1531 strcat(fntemp, suffix);
1532 strcat(fntemp, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1534 /* if we can't rename, we just overwrite the cpt file.
1535 * dangerous if interrupted.
1537 snew(fntemp, strlen(fn));
1541 gmx_ctime_r(&now, timebuf, STRLEN);
1545 fprintf(fplog, "Writing checkpoint, step %s at %s\n\n",
1546 gmx_step_str(step, buf), timebuf);
1549 /* Get offsets for open files */
1550 gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1552 fp = gmx_fio_open(fntemp, "w");
1554 if (state->ekinstate.bUpToDate)
1557 ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) |
1558 (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) |
1559 (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1567 if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1569 flags_enh |= (1<<eenhENERGY_N);
1570 if (state->enerhist.nsum > 0)
1572 flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1573 (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSUM));
1575 if (state->enerhist.nsum_sim > 0)
1577 flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSTEPS_SIM) |
1578 (1<<eenhENERGY_NSUM_SIM));
1580 if (state->enerhist.dht)
1582 flags_enh |= ( (1<< eenhENERGY_DELTA_H_NN) |
1583 (1<< eenhENERGY_DELTA_H_LIST) |
1584 (1<< eenhENERGY_DELTA_H_STARTTIME) |
1585 (1<< eenhENERGY_DELTA_H_STARTLAMBDA) );
1591 flags_dfh = ((1<<edfhBEQUIL) | (1<<edfhNATLAMBDA) | (1<<edfhSUMWEIGHTS) | (1<<edfhSUMDG) |
1592 (1<<edfhTIJ) | (1<<edfhTIJEMP));
1595 flags_dfh |= ((1<<edfhWLDELTA) | (1<<edfhWLHISTO));
1597 if ((elamstats == elamstatsMINVAR) || (elamstats == elamstatsBARKER) || (elamstats == elamstatsMETROPOLIS))
1599 flags_dfh |= ((1<<edfhACCUMP) | (1<<edfhACCUMM) | (1<<edfhACCUMP2) | (1<<edfhACCUMM2)
1600 | (1<<edfhSUMMINVAR) | (1<<edfhSUMVAR));
1608 /* We can check many more things now (CPU, acceleration, etc), but
1609 * it is highly unlikely to have two separate builds with exactly
1610 * the same version, user, time, and build host!
1613 version = gmx_strdup(gmx_version());
1614 btime = gmx_strdup(BUILD_TIME);
1615 buser = gmx_strdup(BUILD_USER);
1616 bhost = gmx_strdup(BUILD_HOST);
1618 double_prec = GMX_CPT_BUILD_DP;
1619 fprog = gmx_strdup(Program());
1621 ftime = &(timebuf[0]);
1623 do_cpt_header(gmx_fio_getxdr(fp), FALSE, &file_version,
1624 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1625 &eIntegrator, &simulation_part, &step, &t, &nppnodes,
1626 DOMAINDECOMP(cr) ? cr->dd->nc : NULL, &npmenodes,
1627 &state->natoms, &state->ngtc, &state->nnhpres,
1628 &state->nhchainlength, &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
1629 &state->edsamstate.nED, &state->swapstate.eSwapCoords,
1638 if ((do_cpt_state(gmx_fio_getxdr(fp), FALSE, state->flags, state, NULL) < 0) ||
1639 (do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL) < 0) ||
1640 (do_cpt_enerhist(gmx_fio_getxdr(fp), FALSE, flags_enh, &state->enerhist, NULL) < 0) ||
1641 (do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL) < 0) ||
1642 (do_cpt_EDstate(gmx_fio_getxdr(fp), FALSE, &state->edsamstate, NULL) < 0) ||
1643 (do_cpt_swapstate(gmx_fio_getxdr(fp), FALSE, &state->swapstate, NULL) < 0) ||
1644 (do_cpt_files(gmx_fio_getxdr(fp), FALSE, &outputfiles, &noutputfiles, NULL,
1647 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1650 do_cpt_footer(gmx_fio_getxdr(fp), file_version);
1652 /* we really, REALLY, want to make sure to physically write the checkpoint,
1653 and all the files it depends on, out to disk. Because we've
1654 opened the checkpoint with gmx_fio_open(), it's in our list
1656 ret = gmx_fio_all_output_fsync();
1662 "Cannot fsync '%s'; maybe you are out of disk space?",
1663 gmx_fio_getname(ret));
1665 if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV) == NULL)
1675 if (gmx_fio_close(fp) != 0)
1677 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1680 /* we don't move the checkpoint if the user specified they didn't want it,
1681 or if the fsyncs failed */
1682 #ifndef GMX_NO_RENAME
1683 if (!bNumberAndKeep && !ret)
1687 /* Rename the previous checkpoint file */
1689 buf[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1690 strcat(buf, "_prev");
1691 strcat(buf, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1693 /* we copy here so that if something goes wrong between now and
1694 * the rename below, there's always a state.cpt.
1695 * If renames are atomic (such as in POSIX systems),
1696 * this copying should be unneccesary.
1698 gmx_file_copy(fn, buf, FALSE);
1699 /* We don't really care if this fails:
1700 * there's already a new checkpoint.
1703 gmx_file_rename(fn, buf);
1706 if (gmx_file_rename(fntemp, fn) != 0)
1708 gmx_file("Cannot rename checkpoint file; maybe you are out of disk space?");
1711 #endif /* GMX_NO_RENAME */
1717 /*code for alternate checkpointing scheme. moved from top of loop over
1719 fcRequestCheckPoint();
1720 if (fcCheckPointParallel( cr->nodeid, NULL, 0) == 0)
1722 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", step );
1724 #endif /* end GMX_FAHCORE block */
1727 static void print_flag_mismatch(FILE *fplog, int sflags, int fflags)
1731 fprintf(fplog, "\nState entry mismatch between the simulation and the checkpoint file\n");
1732 fprintf(fplog, "Entries which are not present in the checkpoint file will not be updated\n");
1733 fprintf(fplog, " %24s %11s %11s\n", "", "simulation", "checkpoint");
1734 for (i = 0; i < estNR; i++)
1736 if ((sflags & (1<<i)) || (fflags & (1<<i)))
1738 fprintf(fplog, " %24s %11s %11s\n",
1740 (sflags & (1<<i)) ? " present " : "not present",
1741 (fflags & (1<<i)) ? " present " : "not present");
1746 static void check_int(FILE *fplog, const char *type, int p, int f, gmx_bool *mm)
1748 FILE *fp = fplog ? fplog : stderr;
1752 fprintf(fp, " %s mismatch,\n", type);
1753 fprintf(fp, " current program: %d\n", p);
1754 fprintf(fp, " checkpoint file: %d\n", f);
1760 static void check_string(FILE *fplog, const char *type, const char *p,
1761 const char *f, gmx_bool *mm)
1763 FILE *fp = fplog ? fplog : stderr;
1765 if (strcmp(p, f) != 0)
1767 fprintf(fp, " %s mismatch,\n", type);
1768 fprintf(fp, " current program: %s\n", p);
1769 fprintf(fp, " checkpoint file: %s\n", f);
1775 static void check_match(FILE *fplog,
1777 char *btime, char *buser, char *bhost, int double_prec,
1779 t_commrec *cr, int npp_f, int npme_f,
1780 ivec dd_nc, ivec dd_nc_f)
1787 check_string(fplog, "Version", gmx_version(), version, &mm);
1788 check_string(fplog, "Build time", BUILD_TIME, btime, &mm);
1789 check_string(fplog, "Build user", BUILD_USER, buser, &mm);
1790 check_string(fplog, "Build host", BUILD_HOST, bhost, &mm);
1791 check_int (fplog, "Double prec.", GMX_CPT_BUILD_DP, double_prec, &mm);
1792 check_string(fplog, "Program name", Program(), fprog, &mm);
1794 check_int (fplog, "#nodes", cr->nnodes, npp_f+npme_f, &mm);
1797 check_int (fplog, "#PME-nodes", cr->npmenodes, npme_f, &mm);
1800 if (cr->npmenodes >= 0)
1802 npp -= cr->npmenodes;
1806 check_int (fplog, "#DD-cells[x]", dd_nc[XX], dd_nc_f[XX], &mm);
1807 check_int (fplog, "#DD-cells[y]", dd_nc[YY], dd_nc_f[YY], &mm);
1808 check_int (fplog, "#DD-cells[z]", dd_nc[ZZ], dd_nc_f[ZZ], &mm);
1815 "Gromacs binary or parallel settings not identical to previous run.\n"
1816 "Continuation is exact, but is not guaranteed to be binary identical%s.\n\n",
1817 fplog ? ",\n see the log file for details" : "");
1822 "Gromacs binary or parallel settings not identical to previous run.\n"
1823 "Continuation is exact, but is not guaranteed to be binary identical.\n\n");
1828 static void read_checkpoint(const char *fn, FILE **pfplog,
1829 t_commrec *cr, ivec dd_nc,
1830 int eIntegrator, int *init_fep_state, gmx_int64_t *step, double *t,
1831 t_state *state, gmx_bool *bReadEkin,
1832 int *simulation_part,
1833 gmx_bool bAppendOutputFiles, gmx_bool bForceAppend)
1838 char *version, *btime, *buser, *bhost, *fprog, *ftime;
1840 char buf[STEPSTRSIZE];
1841 int eIntegrator_f, nppnodes_f, npmenodes_f;
1843 int natoms, ngtc, nnhpres, nhchainlength, nlambda, fflags, flags_eks, flags_enh, flags_dfh;
1846 gmx_file_position_t *outputfiles;
1848 t_fileio *chksum_file;
1849 FILE * fplog = *pfplog;
1850 unsigned char digest[16];
1851 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1852 struct flock fl; /* don't initialize here: the struct order is OS
1856 const char *int_warn =
1857 "WARNING: The checkpoint file was generated with integrator %s,\n"
1858 " while the simulation uses integrator %s\n\n";
1860 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1861 fl.l_type = F_WRLCK;
1862 fl.l_whence = SEEK_SET;
1868 fp = gmx_fio_open(fn, "r");
1869 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
1870 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1871 &eIntegrator_f, simulation_part, step, t,
1872 &nppnodes_f, dd_nc_f, &npmenodes_f,
1873 &natoms, &ngtc, &nnhpres, &nhchainlength, &nlambda,
1874 &fflags, &flags_eks, &flags_enh, &flags_dfh,
1875 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
1877 if (bAppendOutputFiles &&
1878 file_version >= 13 && double_prec != GMX_CPT_BUILD_DP)
1880 gmx_fatal(FARGS, "Output file appending requested, but the code and checkpoint file precision (single/double) don't match");
1883 if (cr == NULL || MASTER(cr))
1885 fprintf(stderr, "\nReading checkpoint file %s generated: %s\n\n",
1889 /* This will not be written if we do appending, since fplog is still NULL then */
1892 fprintf(fplog, "\n");
1893 fprintf(fplog, "Reading checkpoint file %s\n", fn);
1894 fprintf(fplog, " file generated by: %s\n", fprog);
1895 fprintf(fplog, " file generated at: %s\n", ftime);
1896 fprintf(fplog, " GROMACS build time: %s\n", btime);
1897 fprintf(fplog, " GROMACS build user: %s\n", buser);
1898 fprintf(fplog, " GROMACS build host: %s\n", bhost);
1899 fprintf(fplog, " GROMACS double prec.: %d\n", double_prec);
1900 fprintf(fplog, " simulation part #: %d\n", *simulation_part);
1901 fprintf(fplog, " step: %s\n", gmx_step_str(*step, buf));
1902 fprintf(fplog, " time: %f\n", *t);
1903 fprintf(fplog, "\n");
1906 if (natoms != state->natoms)
1908 gmx_fatal(FARGS, "Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms", natoms, state->natoms);
1910 if (ngtc != state->ngtc)
1912 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);
1914 if (nnhpres != state->nnhpres)
1916 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);
1919 if (nlambda != state->dfhist.nlambda)
1921 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);
1924 init_gtc_state(state, state->ngtc, state->nnhpres, nhchainlength); /* need to keep this here to keep the tpr format working */
1925 /* write over whatever was read; we use the number of Nose-Hoover chains from the checkpoint */
1927 if (eIntegrator_f != eIntegrator)
1931 fprintf(stderr, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1933 if (bAppendOutputFiles)
1936 "Output file appending requested, but input/checkpoint integrators do not match.\n"
1937 "Stopping the run to prevent you from ruining all your data...\n"
1938 "If you _really_ know what you are doing, try with the -noappend option.\n");
1942 fprintf(fplog, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1950 else if (cr->nnodes == nppnodes_f + npmenodes_f)
1952 if (cr->npmenodes < 0)
1954 cr->npmenodes = npmenodes_f;
1956 int nppnodes = cr->nnodes - cr->npmenodes;
1957 if (nppnodes == nppnodes_f)
1959 for (d = 0; d < DIM; d++)
1963 dd_nc[d] = dd_nc_f[d];
1969 if (fflags != state->flags)
1974 if (bAppendOutputFiles)
1977 "Output file appending requested, but input and checkpoint states are not identical.\n"
1978 "Stopping the run to prevent you from ruining all your data...\n"
1979 "You can try with the -noappend option, and get more info in the log file.\n");
1982 if (getenv("GMX_ALLOW_CPT_MISMATCH") == NULL)
1984 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");
1989 "WARNING: The checkpoint state entries do not match the simulation,\n"
1990 " see the log file for details\n\n");
1996 print_flag_mismatch(fplog, state->flags, fflags);
2003 check_match(fplog, version, btime, buser, bhost, double_prec, fprog,
2004 cr, nppnodes_f, npmenodes_f, dd_nc, dd_nc_f);
2007 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, fflags, state, NULL);
2008 *init_fep_state = state->fep_state; /* there should be a better way to do this than setting it here.
2009 Investigate for 5.0. */
2014 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2019 *bReadEkin = ((flags_eks & (1<<eeksEKINH)) || (flags_eks & (1<<eeksEKINF)) || (flags_eks & (1<<eeksEKINO)) ||
2020 ((flags_eks & (1<<eeksEKINSCALEF)) | (flags_eks & (1<<eeksEKINSCALEH)) | (flags_eks & (1<<eeksVSCALE))));
2022 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2023 flags_enh, &state->enerhist, NULL);
2029 if (file_version < 6)
2031 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.";
2033 fprintf(stderr, "\nWARNING: %s\n\n", warn);
2036 fprintf(fplog, "\nWARNING: %s\n\n", warn);
2038 state->enerhist.nsum = *step;
2039 state->enerhist.nsum_sim = *step;
2042 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2048 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2054 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2060 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, NULL, file_version);
2066 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2071 if (gmx_fio_close(fp) != 0)
2073 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2082 /* If the user wants to append to output files,
2083 * we use the file pointer positions of the output files stored
2084 * in the checkpoint file and truncate the files such that any frames
2085 * written after the checkpoint time are removed.
2086 * All files are md5sum checked such that we can be sure that
2087 * we do not truncate other (maybe imprortant) files.
2089 if (bAppendOutputFiles)
2091 if (fn2ftp(outputfiles[0].filename) != efLOG)
2093 /* make sure first file is log file so that it is OK to use it for
2096 gmx_fatal(FARGS, "The first output file should always be the log "
2097 "file but instead is: %s. Cannot do appending because of this condition.", outputfiles[0].filename);
2099 for (i = 0; i < nfiles; i++)
2101 if (outputfiles[i].offset < 0)
2103 gmx_fatal(FARGS, "The original run wrote a file called '%s' which "
2104 "is larger than 2 GB, but mdrun did not support large file"
2105 " offsets. Can not append. Run mdrun with -noappend",
2106 outputfiles[i].filename);
2109 chksum_file = gmx_fio_open(outputfiles[i].filename, "a");
2112 chksum_file = gmx_fio_open(outputfiles[i].filename, "r+");
2117 /* Note that there are systems where the lock operation
2118 * will succeed, but a second process can also lock the file.
2119 * We should probably try to detect this.
2121 #if defined __native_client__
2125 #elif defined GMX_NATIVE_WINDOWS
2126 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX) == -1)
2128 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl) == -1)
2131 if (errno == ENOSYS)
2135 gmx_fatal(FARGS, "File locking is not supported on this system. Use -noappend or specify -append explicitly to append anyhow.");
2139 fprintf(stderr, "\nNOTE: File locking is not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2142 fprintf(fplog, "\nNOTE: File locking not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2146 else if (errno == EACCES || errno == EAGAIN)
2148 gmx_fatal(FARGS, "Failed to lock: %s. Already running "
2149 "simulation?", outputfiles[i].filename);
2153 gmx_fatal(FARGS, "Failed to lock: %s. %s.",
2154 outputfiles[i].filename, strerror(errno));
2159 /* compute md5 chksum */
2160 if (outputfiles[i].chksum_size != -1)
2162 if (gmx_fio_get_file_md5(chksum_file, outputfiles[i].offset,
2163 digest) != outputfiles[i].chksum_size) /*at the end of the call the file position is at the end of the file*/
2165 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.",
2166 outputfiles[i].chksum_size,
2167 outputfiles[i].filename);
2170 if (i == 0) /*log file needs to be seeked in case we need to truncate (other files are truncated below)*/
2172 if (gmx_fio_seek(chksum_file, outputfiles[i].offset))
2174 gmx_fatal(FARGS, "Seek error! Failed to truncate log-file: %s.", strerror(errno));
2179 if (i == 0) /*open log file here - so that lock is never lifted
2180 after chksum is calculated */
2182 *pfplog = gmx_fio_getfp(chksum_file);
2186 gmx_fio_close(chksum_file);
2189 /* compare md5 chksum */
2190 if (outputfiles[i].chksum_size != -1 &&
2191 memcmp(digest, outputfiles[i].chksum, 16) != 0)
2195 fprintf(debug, "chksum for %s: ", outputfiles[i].filename);
2196 for (j = 0; j < 16; j++)
2198 fprintf(debug, "%02x", digest[j]);
2200 fprintf(debug, "\n");
2202 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.",
2203 outputfiles[i].filename);
2208 if (i != 0) /*log file is already seeked to correct position */
2210 #ifdef GMX_NATIVE_WINDOWS
2211 rc = gmx_wintruncate(outputfiles[i].filename, outputfiles[i].offset);
2213 rc = 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);
2227 void load_checkpoint(const char *fn, FILE **fplog,
2228 t_commrec *cr, ivec dd_nc,
2229 t_inputrec *ir, t_state *state,
2230 gmx_bool *bReadEkin,
2231 gmx_bool bAppend, gmx_bool bForceAppend)
2238 /* Read the state from the checkpoint file */
2239 read_checkpoint(fn, fplog,
2241 ir->eI, &(ir->fepvals->init_fep_state), &step, &t, state, bReadEkin,
2242 &ir->simulation_part, bAppend, bForceAppend);
2246 gmx_bcast(sizeof(cr->npmenodes), &cr->npmenodes, cr);
2247 gmx_bcast(DIM*sizeof(dd_nc[0]), dd_nc, cr);
2248 gmx_bcast(sizeof(step), &step, cr);
2249 gmx_bcast(sizeof(*bReadEkin), bReadEkin, cr);
2251 ir->bContinuation = TRUE;
2252 if (ir->nsteps >= 0)
2254 ir->nsteps += ir->init_step - step;
2256 ir->init_step = step;
2257 ir->simulation_part += 1;
2260 static void read_checkpoint_data(t_fileio *fp, int *simulation_part,
2261 gmx_int64_t *step, double *t, t_state *state,
2262 int *nfiles, gmx_file_position_t **outputfiles)
2265 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2270 int flags_eks, flags_enh, flags_dfh;
2272 gmx_file_position_t *files_loc = NULL;
2275 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2276 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2277 &eIntegrator, simulation_part, step, t, &nppnodes, dd_nc, &npme,
2278 &state->natoms, &state->ngtc, &state->nnhpres, &state->nhchainlength,
2279 &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
2280 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
2282 do_cpt_state(gmx_fio_getxdr(fp), TRUE, state->flags, state, NULL);
2287 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2292 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2293 flags_enh, &state->enerhist, NULL);
2298 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2304 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2310 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2316 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE,
2317 outputfiles != NULL ? outputfiles : &files_loc,
2318 outputfiles != NULL ? nfiles : &nfiles_loc,
2319 NULL, file_version);
2320 if (files_loc != NULL)
2330 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2344 read_checkpoint_state(const char *fn, int *simulation_part,
2345 gmx_int64_t *step, double *t, t_state *state)
2349 fp = gmx_fio_open(fn, "r");
2350 read_checkpoint_data(fp, simulation_part, step, t, state, NULL, NULL);
2351 if (gmx_fio_close(fp) != 0)
2353 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2357 void read_checkpoint_trxframe(t_fileio *fp, t_trxframe *fr)
2359 /* This next line is nasty because the sub-structures of t_state
2360 * cannot be assumed to be zeroed (or even initialized in ways the
2361 * rest of the code might assume). Using snew would be better, but
2362 * this will all go away for 5.0. */
2364 int simulation_part;
2368 init_state(&state, 0, 0, 0, 0, 0);
2370 read_checkpoint_data(fp, &simulation_part, &step, &t, &state, NULL, NULL);
2372 fr->natoms = state.natoms;
2375 fr->step = gmx_int64_to_int(step,
2376 "conversion of checkpoint to trajectory");
2380 fr->lambda = state.lambda[efptFEP];
2381 fr->fep_state = state.fep_state;
2383 fr->bX = (state.flags & (1<<estX));
2389 fr->bV = (state.flags & (1<<estV));
2396 fr->bBox = (state.flags & (1<<estBOX));
2399 copy_mat(state.box, fr->box);
2404 void list_checkpoint(const char *fn, FILE *out)
2408 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2410 int eIntegrator, simulation_part, nppnodes, npme;
2415 int flags_eks, flags_enh, flags_dfh;
2417 gmx_file_position_t *outputfiles;
2420 init_state(&state, -1, -1, -1, -1, 0);
2422 fp = gmx_fio_open(fn, "r");
2423 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2424 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2425 &eIntegrator, &simulation_part, &step, &t, &nppnodes, dd_nc, &npme,
2426 &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
2427 &(state.dfhist.nlambda), &state.flags,
2428 &flags_eks, &flags_enh, &flags_dfh, &state.edsamstate.nED,
2429 &state.swapstate.eSwapCoords, out);
2430 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, state.flags, &state, out);
2435 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state.ekinstate, out);
2440 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2441 flags_enh, &state.enerhist, out);
2445 ret = do_cpt_df_hist(gmx_fio_getxdr(fp),
2446 flags_dfh, &state.dfhist, out);
2451 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state.edsamstate, out);
2456 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state.swapstate, out);
2461 do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, out, file_version);
2466 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2473 if (gmx_fio_close(fp) != 0)
2475 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2482 static gmx_bool exist_output_file(const char *fnm_cp, int nfile, const t_filenm fnm[])
2486 /* Check if the output file name stored in the checkpoint file
2487 * is one of the output file names of mdrun.
2491 !(is_output(&fnm[i]) && strcmp(fnm_cp, fnm[i].fns[0]) == 0))
2496 return (i < nfile && gmx_fexist(fnm_cp));
2499 /* This routine cannot print tons of data, since it is called before the log file is opened. */
2500 gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
2501 gmx_int64_t *cpt_step, t_commrec *cr,
2502 gmx_bool bAppendReq,
2503 int nfile, const t_filenm fnm[],
2504 const char *part_suffix, gmx_bool *bAddPart)
2507 gmx_int64_t step = 0;
2509 /* This next line is nasty because the sub-structures of t_state
2510 * cannot be assumed to be zeroed (or even initialized in ways the
2511 * rest of the code might assume). Using snew would be better, but
2512 * this will all go away for 5.0. */
2515 gmx_file_position_t *outputfiles;
2518 char *fn, suf_up[STRLEN];
2524 if (!gmx_fexist(filename) || (!(fp = gmx_fio_open(filename, "r")) ))
2526 *simulation_part = 0;
2530 init_state(&state, 0, 0, 0, 0, 0);
2532 read_checkpoint_data(fp, simulation_part, &step, &t, &state,
2533 &nfiles, &outputfiles);
2534 if (gmx_fio_close(fp) != 0)
2536 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2543 for (f = 0; f < nfiles; f++)
2545 if (exist_output_file(outputfiles[f].filename, nfile, fnm))
2550 if (nexist == nfiles)
2552 bAppend = bAppendReq;
2554 else if (nexist > 0)
2557 "Output file appending has been requested,\n"
2558 "but some output files listed in the checkpoint file %s\n"
2559 "are not present or are named differently by the current program:\n",
2561 fprintf(stderr, "output files present:");
2562 for (f = 0; f < nfiles; f++)
2564 if (exist_output_file(outputfiles[f].filename,
2567 fprintf(stderr, " %s", outputfiles[f].filename);
2570 fprintf(stderr, "\n");
2571 fprintf(stderr, "output files not present or named differently:");
2572 for (f = 0; f < nfiles; f++)
2574 if (!exist_output_file(outputfiles[f].filename,
2577 fprintf(stderr, " %s", outputfiles[f].filename);
2580 fprintf(stderr, "\n");
2582 gmx_fatal(FARGS, "File appending requested, but only %d of the %d output files are present", nexist, nfiles);
2590 gmx_fatal(FARGS, "File appending requested, but no output file information is stored in the checkpoint file");
2592 fn = outputfiles[0].filename;
2593 if (strlen(fn) < 4 ||
2594 gmx_strcasecmp(fn+strlen(fn)-4, ftp2ext(efLOG)) == 0)
2596 gmx_fatal(FARGS, "File appending requested, but the log file is not the first file listed in the checkpoint file");
2598 /* Set bAddPart to whether the suffix string '.part' is present
2599 * in the log file name.
2601 strcpy(suf_up, part_suffix);
2603 *bAddPart = (strstr(fn, part_suffix) != NULL ||
2604 strstr(fn, suf_up) != NULL);
2612 gmx_bcast(sizeof(*simulation_part), simulation_part, cr);
2614 if (*simulation_part > 0 && bAppendReq)
2616 gmx_bcast(sizeof(bAppend), &bAppend, cr);
2617 gmx_bcast(sizeof(*bAddPart), bAddPart, cr);
2620 if (NULL != cpt_step)