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. */
50 #ifdef HAVE_SYS_TIME_H
58 #ifdef GMX_NATIVE_WINDOWS
61 #include <sys/locking.h>
67 #include "types/commrec.h"
68 #include "gromacs/utility/smalloc.h"
70 #include "gromacs/math/vec.h"
72 #include "checkpoint.h"
73 #include "gromacs/utility/cstringutil.h"
75 #include "gromacs/fileio/filenm.h"
76 #include "gromacs/utility/futil.h"
77 #include "gromacs/fileio/gmxfio.h"
78 #include "gromacs/fileio/xdrf.h"
79 #include "gromacs/fileio/xdr_datatype.h"
80 #include "gromacs/utility/basenetwork.h"
81 #include "gromacs/utility/baseversion.h"
82 #include "gromacs/utility/fatalerror.h"
84 #include "buildinfo.h"
90 #define CPT_MAGIC1 171817
91 #define CPT_MAGIC2 171819
92 #define CPTSTRLEN 1024
95 #define GMX_CPT_BUILD_DP 1
97 #define GMX_CPT_BUILD_DP 0
100 /* cpt_version should normally only be changed
101 * when the header of footer format changes.
102 * The state data format itself is backward and forward compatible.
103 * But old code can not read a new entry that is present in the file
104 * (but can read a new format when new entries are not present).
106 static const int cpt_version = 16;
109 const char *est_names[estNR] =
112 "box", "box-rel", "box-v", "pres_prev",
113 "nosehoover-xi", "thermostat-integral",
114 "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
115 "disre_initf", "disre_rm3tav",
116 "orire_initf", "orire_Dtav",
117 "svir_prev", "nosehoover-vxi", "v_eta", "vol0", "nhpres_xi", "nhpres_vxi", "fvir_prev", "fep_state", "MC-rng", "MC-rng-i"
121 eeksEKIN_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksEKINF, eeksEKINO, eeksEKINSCALEF, eeksEKINSCALEH, eeksVSCALE, eeksEKINTOTAL, eeksNR
124 const char *eeks_names[eeksNR] =
126 "Ekin_n", "Ekinh", "dEkindlambda", "mv_cos",
127 "Ekinf", "Ekinh_old", "EkinScaleF_NHC", "EkinScaleH_NHC", "Vscale_NHC", "Ekin_Total"
131 eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
132 eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
133 eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM,
134 eenhENERGY_DELTA_H_NN,
135 eenhENERGY_DELTA_H_LIST,
136 eenhENERGY_DELTA_H_STARTTIME,
137 eenhENERGY_DELTA_H_STARTLAMBDA,
141 const char *eenh_names[eenhNR] =
143 "energy_n", "energy_aver", "energy_sum", "energy_nsum",
144 "energy_sum_sim", "energy_nsum_sim",
145 "energy_nsteps", "energy_nsteps_sim",
147 "energy_delta_h_list",
148 "energy_delta_h_start_time",
149 "energy_delta_h_start_lambda"
152 /* free energy history variables -- need to be preserved over checkpoint */
154 edfhBEQUIL, edfhNATLAMBDA, edfhWLHISTO, edfhWLDELTA, edfhSUMWEIGHTS, edfhSUMDG, edfhSUMMINVAR, edfhSUMVAR,
155 edfhACCUMP, edfhACCUMM, edfhACCUMP2, edfhACCUMM2, edfhTIJ, edfhTIJEMP, edfhNR
157 /* free energy history variable names */
158 const char *edfh_names[edfhNR] =
160 "bEquilibrated", "N_at_state", "Wang-Landau Histogram", "Wang-Landau Delta", "Weights", "Free Energies", "minvar", "variance",
161 "accumulated_plus", "accumulated_minus", "accumulated_plus_2", "accumulated_minus_2", "Tij", "Tij_empirical"
164 #ifdef GMX_NATIVE_WINDOWS
166 gmx_wintruncate(const char *filename, __int64 size)
169 /*we do this elsewhere*/
175 fp = fopen(filename, "rb+");
182 return _chsize_s( fileno(fp), size);
189 ecprREAL, ecprRVEC, ecprMATRIX
193 cptpEST, cptpEEKS, cptpEENH, cptpEDFH
195 /* enums for the different components of checkpoint variables, replacing the hard coded ones.
196 cptpEST - state variables.
197 cptpEEKS - Kinetic energy state variables.
198 cptpEENH - Energy history state variables.
199 cptpEDFH - free energy history variables.
203 static const char *st_names(int cptp, int ecpt)
207 case cptpEST: return est_names [ecpt]; break;
208 case cptpEEKS: return eeks_names[ecpt]; break;
209 case cptpEENH: return eenh_names[ecpt]; break;
210 case cptpEDFH: return edfh_names[ecpt]; break;
216 static void cp_warning(FILE *fp)
218 fprintf(fp, "\nWARNING: Checkpoint file is corrupted or truncated\n\n");
221 static void cp_error()
223 gmx_fatal(FARGS, "Checkpoint file corrupted/truncated, or maybe you are out of disk space?");
226 static void do_cpt_string_err(XDR *xd, gmx_bool bRead, const char *desc, char **s, FILE *list)
234 res = xdr_string(xd, s, CPTSTRLEN);
241 fprintf(list, "%s = %s\n", desc, *s);
246 static int do_cpt_int(XDR *xd, const char *desc, int *i, FILE *list)
250 res = xdr_int(xd, i);
257 fprintf(list, "%s = %d\n", desc, *i);
262 static int do_cpt_u_chars(XDR *xd, const char *desc, int n, unsigned char *i, FILE *list)
268 fprintf(list, "%s = ", desc);
270 for (j = 0; j < n && res; j++)
272 res &= xdr_u_char(xd, &i[j]);
275 fprintf(list, "%02x", i[j]);
290 static void do_cpt_int_err(XDR *xd, const char *desc, int *i, FILE *list)
292 if (do_cpt_int(xd, desc, i, list) < 0)
298 static void do_cpt_step_err(XDR *xd, const char *desc, gmx_int64_t *i, FILE *list)
301 char buf[STEPSTRSIZE];
303 res = xdr_int64(xd, i);
310 fprintf(list, "%s = %s\n", desc, gmx_step_str(*i, buf));
314 static void do_cpt_double_err(XDR *xd, const char *desc, double *f, FILE *list)
318 res = xdr_double(xd, f);
325 fprintf(list, "%s = %f\n", desc, *f);
329 static void do_cpt_real_err(XDR *xd, real *f)
334 res = xdr_double(xd, f);
336 res = xdr_float(xd, f);
344 static void do_cpt_n_rvecs_err(XDR *xd, const char *desc, int n, rvec f[], FILE *list)
348 for (i = 0; i < n; i++)
350 for (j = 0; j < DIM; j++)
352 do_cpt_real_err(xd, &f[i][j]);
358 pr_rvecs(list, 0, desc, f, n);
362 /* If nval >= 0, nval is used; on read this should match the passed value.
363 * If nval n<0, *nptr is used; on read the value is stored in nptr
365 static int do_cpte_reals_low(XDR *xd, int cptp, int ecpt, int sflags,
366 int nval, int *nptr, real **v,
367 FILE *list, int erealtype)
371 int dtc = xdr_datatype_float;
373 int dtc = xdr_datatype_double;
375 real *vp, *va = NULL;
390 gmx_incons("*ntpr=NULL in do_cpte_reals_low");
395 res = xdr_int(xd, &nf);
406 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), nval, nf);
415 res = xdr_int(xd, &dt);
422 fprintf(stderr, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
423 st_names(cptp, ecpt), xdr_datatype_names[dtc],
424 xdr_datatype_names[dt]);
426 if (list || !(sflags & (1<<ecpt)))
439 if (dt == xdr_datatype_float)
441 if (dtc == xdr_datatype_float)
449 res = xdr_vector(xd, (char *)vf, nf,
450 (unsigned int)sizeof(float), (xdrproc_t)xdr_float);
455 if (dtc != xdr_datatype_float)
457 for (i = 0; i < nf; i++)
466 if (dtc == xdr_datatype_double)
468 /* cppcheck-suppress invalidPointerCast
469 * Only executed if real is anyhow double */
476 res = xdr_vector(xd, (char *)vd, nf,
477 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
482 if (dtc != xdr_datatype_double)
484 for (i = 0; i < nf; i++)
497 pr_reals(list, 0, st_names(cptp, ecpt), vp, nf);
500 pr_rvecs(list, 0, st_names(cptp, ecpt), (rvec *)vp, nf/3);
503 gmx_incons("Unknown checkpoint real type");
515 /* This function stores n along with the reals for reading,
516 * but on reading it assumes that n matches the value in the checkpoint file,
517 * a fatal error is generated when this is not the case.
519 static int do_cpte_reals(XDR *xd, int cptp, int ecpt, int sflags,
520 int n, real **v, FILE *list)
522 return do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, v, list, ecprREAL);
525 /* This function does the same as do_cpte_reals,
526 * except that on reading it ignores the passed value of *n
527 * and stored the value read from the checkpoint file in *n.
529 static int do_cpte_n_reals(XDR *xd, int cptp, int ecpt, int sflags,
530 int *n, real **v, FILE *list)
532 return do_cpte_reals_low(xd, cptp, ecpt, sflags, -1, n, v, list, ecprREAL);
535 static int do_cpte_real(XDR *xd, int cptp, int ecpt, int sflags,
540 return do_cpte_reals_low(xd, cptp, ecpt, sflags, 1, NULL, &r, list, ecprREAL);
543 static int do_cpte_ints(XDR *xd, int cptp, int ecpt, int sflags,
544 int n, int **v, FILE *list)
547 int dtc = xdr_datatype_int;
552 res = xdr_int(xd, &nf);
557 if (list == NULL && v != NULL && nf != n)
559 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
562 res = xdr_int(xd, &dt);
569 gmx_fatal(FARGS, "Type mismatch for state entry %s, code type is %s, file type is %s\n",
570 st_names(cptp, ecpt), xdr_datatype_names[dtc],
571 xdr_datatype_names[dt]);
573 if (list || !(sflags & (1<<ecpt)) || v == NULL)
586 res = xdr_vector(xd, (char *)vp, nf,
587 (unsigned int)sizeof(int), (xdrproc_t)xdr_int);
594 pr_ivec(list, 0, st_names(cptp, ecpt), vp, nf, TRUE);
604 static int do_cpte_int(XDR *xd, int cptp, int ecpt, int sflags,
607 return do_cpte_ints(xd, cptp, ecpt, sflags, 1, &i, list);
610 static int do_cpte_doubles(XDR *xd, int cptp, int ecpt, int sflags,
611 int n, double **v, FILE *list)
614 int dtc = xdr_datatype_double;
615 double *vp, *va = NULL;
619 res = xdr_int(xd, &nf);
624 if (list == NULL && nf != n)
626 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
629 res = xdr_int(xd, &dt);
636 gmx_fatal(FARGS, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
637 st_names(cptp, ecpt), xdr_datatype_names[dtc],
638 xdr_datatype_names[dt]);
640 if (list || !(sflags & (1<<ecpt)))
653 res = xdr_vector(xd, (char *)vp, nf,
654 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
661 pr_doubles(list, 0, st_names(cptp, ecpt), vp, nf);
671 static int do_cpte_double(XDR *xd, int cptp, int ecpt, int sflags,
672 double *r, FILE *list)
674 return do_cpte_doubles(xd, cptp, ecpt, sflags, 1, &r, list);
678 static int do_cpte_rvecs(XDR *xd, int cptp, int ecpt, int sflags,
679 int n, rvec **v, FILE *list)
683 return do_cpte_reals_low(xd, cptp, ecpt, sflags,
684 n*DIM, NULL, (real **)v, list, ecprRVEC);
687 static int do_cpte_matrix(XDR *xd, int cptp, int ecpt, int sflags,
688 matrix v, FILE *list)
693 vr = (real *)&(v[0][0]);
694 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
695 DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
697 if (list && ret == 0)
699 pr_rvecs(list, 0, st_names(cptp, ecpt), v, DIM);
706 static int do_cpte_nmatrix(XDR *xd, int cptp, int ecpt, int sflags,
707 int n, real **v, FILE *list)
712 char name[CPTSTRLEN];
719 for (i = 0; i < n; i++)
722 reti = do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, &(v[i]), NULL, ecprREAL);
723 if (list && reti == 0)
725 sprintf(name, "%s[%d]", st_names(cptp, ecpt), i);
726 pr_reals(list, 0, name, v[i], n);
736 static int do_cpte_matrices(XDR *xd, int cptp, int ecpt, int sflags,
737 int n, matrix **v, FILE *list)
740 matrix *vp, *va = NULL;
746 res = xdr_int(xd, &nf);
751 if (list == NULL && nf != n)
753 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
755 if (list || !(sflags & (1<<ecpt)))
768 snew(vr, nf*DIM*DIM);
769 for (i = 0; i < nf; i++)
771 for (j = 0; j < DIM; j++)
773 for (k = 0; k < DIM; k++)
775 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
779 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
780 nf*DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
781 for (i = 0; i < nf; i++)
783 for (j = 0; j < DIM; j++)
785 for (k = 0; k < DIM; k++)
787 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
793 if (list && ret == 0)
795 for (i = 0; i < nf; i++)
797 pr_rvecs(list, 0, st_names(cptp, ecpt), vp[i], DIM);
808 static void do_cpt_header(XDR *xd, gmx_bool bRead, int *file_version,
809 char **version, char **btime, char **buser, char **bhost,
811 char **fprog, char **ftime,
812 int *eIntegrator, int *simulation_part,
813 gmx_int64_t *step, double *t,
814 int *nnodes, int *dd_nc, int *npme,
815 int *natoms, int *ngtc, int *nnhpres, int *nhchainlength,
816 int *nlambda, int *flags_state,
817 int *flags_eks, int *flags_enh, int *flags_dfh,
818 int *nED, int *eSwapCoords,
835 res = xdr_int(xd, &magic);
838 gmx_fatal(FARGS, "The checkpoint file is empty/corrupted, or maybe you are out of disk space?");
840 if (magic != CPT_MAGIC1)
842 gmx_fatal(FARGS, "Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
843 "The checkpoint file is corrupted or not a checkpoint file",
849 gmx_gethostname(fhost, 255);
851 do_cpt_string_err(xd, bRead, "GROMACS version", version, list);
852 do_cpt_string_err(xd, bRead, "GROMACS build time", btime, list);
853 do_cpt_string_err(xd, bRead, "GROMACS build user", buser, list);
854 do_cpt_string_err(xd, bRead, "GROMACS build host", bhost, list);
855 do_cpt_string_err(xd, bRead, "generating program", fprog, list);
856 do_cpt_string_err(xd, bRead, "generation time", ftime, list);
857 *file_version = cpt_version;
858 do_cpt_int_err(xd, "checkpoint file version", file_version, list);
859 if (*file_version > cpt_version)
861 gmx_fatal(FARGS, "Attempting to read a checkpoint file of version %d with code of version %d\n", *file_version, cpt_version);
863 if (*file_version >= 13)
865 do_cpt_int_err(xd, "GROMACS double precision", double_prec, list);
871 if (*file_version >= 12)
873 do_cpt_string_err(xd, bRead, "generating host", &fhost, list);
879 do_cpt_int_err(xd, "#atoms", natoms, list);
880 do_cpt_int_err(xd, "#T-coupling groups", ngtc, list);
881 if (*file_version >= 10)
883 do_cpt_int_err(xd, "#Nose-Hoover T-chains", nhchainlength, list);
889 if (*file_version >= 11)
891 do_cpt_int_err(xd, "#Nose-Hoover T-chains for barostat ", nnhpres, list);
897 if (*file_version >= 14)
899 do_cpt_int_err(xd, "# of total lambda states ", nlambda, list);
905 do_cpt_int_err(xd, "integrator", eIntegrator, list);
906 if (*file_version >= 3)
908 do_cpt_int_err(xd, "simulation part #", simulation_part, list);
912 *simulation_part = 1;
914 if (*file_version >= 5)
916 do_cpt_step_err(xd, "step", step, list);
920 do_cpt_int_err(xd, "step", &idum, list);
923 do_cpt_double_err(xd, "t", t, list);
924 do_cpt_int_err(xd, "#PP-nodes", nnodes, list);
926 do_cpt_int_err(xd, "dd_nc[x]", dd_nc ? &(dd_nc[0]) : &idum, list);
927 do_cpt_int_err(xd, "dd_nc[y]", dd_nc ? &(dd_nc[1]) : &idum, list);
928 do_cpt_int_err(xd, "dd_nc[z]", dd_nc ? &(dd_nc[2]) : &idum, list);
929 do_cpt_int_err(xd, "#PME-only nodes", npme, list);
930 do_cpt_int_err(xd, "state flags", flags_state, list);
931 if (*file_version >= 4)
933 do_cpt_int_err(xd, "ekin data flags", flags_eks, list);
934 do_cpt_int_err(xd, "energy history flags", flags_enh, list);
939 *flags_enh = (*flags_state >> (estORIRE_DTAV+1));
940 *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
941 (1<<(estORIRE_DTAV+2)) |
942 (1<<(estORIRE_DTAV+3))));
944 if (*file_version >= 14)
946 do_cpt_int_err(xd, "df history flags", flags_dfh, list);
953 if (*file_version >= 15)
955 do_cpt_int_err(xd, "ED data sets", nED, list);
961 if (*file_version >= 16)
963 do_cpt_int_err(xd, "swap", eSwapCoords, list);
967 static int do_cpt_footer(XDR *xd, int file_version)
972 if (file_version >= 2)
975 res = xdr_int(xd, &magic);
980 if (magic != CPT_MAGIC2)
989 static int do_cpt_state(XDR *xd, gmx_bool bRead,
990 int fflags, t_state *state,
1000 nnht = state->nhchainlength*state->ngtc;
1001 nnhtp = state->nhchainlength*state->nnhpres;
1003 if (bRead) /* we need to allocate space for dfhist if we are reading */
1005 init_df_history(&state->dfhist, state->dfhist.nlambda);
1008 sflags = state->flags;
1009 for (i = 0; (i < estNR && ret == 0); i++)
1011 if (fflags & (1<<i))
1015 case estLAMBDA: ret = do_cpte_reals(xd, cptpEST, i, sflags, efptNR, &(state->lambda), list); break;
1016 case estFEPSTATE: ret = do_cpte_int (xd, cptpEST, i, sflags, &state->fep_state, list); break;
1017 case estBOX: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box, list); break;
1018 case estBOX_REL: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box_rel, list); break;
1019 case estBOXV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->boxv, list); break;
1020 case estPRES_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->pres_prev, list); break;
1021 case estSVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->svir_prev, list); break;
1022 case estFVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->fvir_prev, list); break;
1023 case estNH_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_xi, list); break;
1024 case estNH_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_vxi, list); break;
1025 case estNHPRES_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_xi, list); break;
1026 case estNHPRES_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_vxi, list); break;
1027 case estTC_INT: ret = do_cpte_doubles(xd, cptpEST, i, sflags, state->ngtc, &state->therm_integral, list); break;
1028 case estVETA: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->veta, list); break;
1029 case estVOL0: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->vol0, list); break;
1030 case estX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->x, list); break;
1031 case estV: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->v, list); break;
1032 case estSDX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->sd_X, list); break;
1033 /* The RNG entries are no longer written,
1034 * the next 4 lines are only for reading old files.
1036 case estLD_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1037 case estLD_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1038 case estMC_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1039 case estMC_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1040 case estDISRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.disre_initf, list); break;
1041 case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.ndisrepairs, &state->hist.disre_rm3tav, list); break;
1042 case estORIRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.orire_initf, list); break;
1043 case estORIRE_DTAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.norire_Dtav, &state->hist.orire_Dtav, list); break;
1045 gmx_fatal(FARGS, "Unknown state entry %d\n"
1046 "You are probably reading a new checkpoint file with old code", i);
1054 static int do_cpt_ekinstate(XDR *xd, int fflags, ekinstate_t *ekins,
1062 for (i = 0; (i < eeksNR && ret == 0); i++)
1064 if (fflags & (1<<i))
1069 case eeksEKIN_N: ret = do_cpte_int(xd, cptpEEKS, i, fflags, &ekins->ekin_n, list); break;
1070 case eeksEKINH: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh, list); break;
1071 case eeksEKINF: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinf, list); break;
1072 case eeksEKINO: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh_old, list); break;
1073 case eeksEKINTOTAL: ret = do_cpte_matrix(xd, cptpEEKS, i, fflags, ekins->ekin_total, list); break;
1074 case eeksEKINSCALEF: ret = do_cpte_doubles(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinscalef_nhc, list); break;
1075 case eeksVSCALE: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->vscale_nhc, list); break;
1076 case eeksEKINSCALEH: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->ekinscaleh_nhc, list); break;
1077 case eeksDEKINDL: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->dekindl, list); break;
1078 case eeksMVCOS: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->mvcos, list); break;
1080 gmx_fatal(FARGS, "Unknown ekin data state entry %d\n"
1081 "You are probably reading a new checkpoint file with old code", i);
1090 static int do_cpt_swapstate(XDR *xd, gmx_bool bRead, swapstate_t *swapstate, FILE *list)
1094 int swap_cpt_version = 1;
1097 if (eswapNO == swapstate->eSwapCoords)
1102 swapstate->bFromCpt = bRead;
1104 do_cpt_int_err(xd, "swap checkpoint version", &swap_cpt_version, list);
1105 do_cpt_int_err(xd, "swap coupling steps", &swapstate->nAverage, list);
1107 /* When reading, init_swapcoords has not been called yet,
1108 * so we have to allocate memory first. */
1110 for (ic = 0; ic < eCompNR; ic++)
1112 for (ii = 0; ii < eIonNR; ii++)
1116 do_cpt_int_err(xd, "swap requested atoms", &swapstate->nat_req[ic][ii], list);
1120 do_cpt_int_err(xd, "swap requested atoms p", swapstate->nat_req_p[ic][ii], list);
1125 do_cpt_int_err(xd, "swap influx netto", &swapstate->inflow_netto[ic][ii], list);
1129 do_cpt_int_err(xd, "swap influx netto p", swapstate->inflow_netto_p[ic][ii], list);
1132 if (bRead && (NULL == swapstate->nat_past[ic][ii]) )
1134 snew(swapstate->nat_past[ic][ii], swapstate->nAverage);
1137 for (j = 0; j < swapstate->nAverage; j++)
1141 do_cpt_int_err(xd, "swap past atom counts", &swapstate->nat_past[ic][ii][j], list);
1145 do_cpt_int_err(xd, "swap past atom counts p", &swapstate->nat_past_p[ic][ii][j], list);
1151 /* Ion flux per channel */
1152 for (ic = 0; ic < eChanNR; ic++)
1154 for (ii = 0; ii < eIonNR; ii++)
1158 do_cpt_int_err(xd, "channel flux", &swapstate->fluxfromAtoB[ic][ii], list);
1162 do_cpt_int_err(xd, "channel flux p", swapstate->fluxfromAtoB_p[ic][ii], list);
1167 /* Ion flux leakage */
1170 snew(swapstate->fluxleak, 1);
1172 do_cpt_int_err(xd, "flux leakage", swapstate->fluxleak, list);
1175 do_cpt_int_err(xd, "number of ions", &swapstate->nions, list);
1179 snew(swapstate->channel_label, swapstate->nions);
1180 snew(swapstate->comp_from, swapstate->nions);
1183 do_cpt_u_chars(xd, "channel history", swapstate->nions, swapstate->channel_label, list);
1184 do_cpt_u_chars(xd, "domain history", swapstate->nions, swapstate->comp_from, list);
1186 /* Save the last known whole positions to checkpoint
1187 * file to be able to also make multimeric channels whole in PBC */
1188 do_cpt_int_err(xd, "Ch0 atoms", &swapstate->nat[eChan0], list);
1189 do_cpt_int_err(xd, "Ch1 atoms", &swapstate->nat[eChan1], list);
1192 snew(swapstate->xc_old_whole[eChan0], swapstate->nat[eChan0]);
1193 snew(swapstate->xc_old_whole[eChan1], swapstate->nat[eChan1]);
1194 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], swapstate->xc_old_whole[eChan0], list);
1195 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], swapstate->xc_old_whole[eChan1], list);
1199 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], *swapstate->xc_old_whole_p[eChan0], list);
1200 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], *swapstate->xc_old_whole_p[eChan1], list);
1207 static int do_cpt_enerhist(XDR *xd, gmx_bool bRead,
1208 int fflags, energyhistory_t *enerhist,
1219 enerhist->nsteps = 0;
1221 enerhist->nsteps_sim = 0;
1222 enerhist->nsum_sim = 0;
1223 enerhist->dht = NULL;
1225 if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
1227 snew(enerhist->dht, 1);
1228 enerhist->dht->ndh = NULL;
1229 enerhist->dht->dh = NULL;
1230 enerhist->dht->start_lambda_set = FALSE;
1234 for (i = 0; (i < eenhNR && ret == 0); i++)
1236 if (fflags & (1<<i))
1240 case eenhENERGY_N: ret = do_cpte_int(xd, cptpEENH, i, fflags, &enerhist->nener, list); break;
1241 case eenhENERGY_AVER: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_ave, list); break;
1242 case eenhENERGY_SUM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum, list); break;
1243 case eenhENERGY_NSUM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum, list); break;
1244 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum_sim, list); break;
1245 case eenhENERGY_NSUM_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum_sim, list); break;
1246 case eenhENERGY_NSTEPS: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps, list); break;
1247 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps_sim, list); break;
1248 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd, eenh_names[i], &(enerhist->dht->nndh), list);
1249 if (bRead) /* now allocate memory for it */
1251 snew(enerhist->dht->dh, enerhist->dht->nndh);
1252 snew(enerhist->dht->ndh, enerhist->dht->nndh);
1253 for (j = 0; j < enerhist->dht->nndh; j++)
1255 enerhist->dht->ndh[j] = 0;
1256 enerhist->dht->dh[j] = NULL;
1260 case eenhENERGY_DELTA_H_LIST:
1261 for (j = 0; j < enerhist->dht->nndh; j++)
1263 ret = do_cpte_n_reals(xd, cptpEENH, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list);
1266 case eenhENERGY_DELTA_H_STARTTIME:
1267 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_time), list); break;
1268 case eenhENERGY_DELTA_H_STARTLAMBDA:
1269 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_lambda), list); break;
1271 gmx_fatal(FARGS, "Unknown energy history entry %d\n"
1272 "You are probably reading a new checkpoint file with old code", i);
1277 if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
1279 /* Assume we have an old file format and copy sum to sum_sim */
1280 srenew(enerhist->ener_sum_sim, enerhist->nener);
1281 for (i = 0; i < enerhist->nener; i++)
1283 enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1287 if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1288 !(fflags & (1<<eenhENERGY_NSTEPS)))
1290 /* Assume we have an old file format and copy nsum to nsteps */
1291 enerhist->nsteps = enerhist->nsum;
1293 if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1294 !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1296 /* Assume we have an old file format and copy nsum to nsteps */
1297 enerhist->nsteps_sim = enerhist->nsum_sim;
1303 static int do_cpt_df_hist(XDR *xd, int fflags, df_history_t *dfhist, FILE *list)
1308 nlambda = dfhist->nlambda;
1311 for (i = 0; (i < edfhNR && ret == 0); i++)
1313 if (fflags & (1<<i))
1317 case edfhBEQUIL: ret = do_cpte_int(xd, cptpEDFH, i, fflags, &dfhist->bEquil, list); break;
1318 case edfhNATLAMBDA: ret = do_cpte_ints(xd, cptpEDFH, i, fflags, nlambda, &dfhist->n_at_lam, list); break;
1319 case edfhWLHISTO: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->wl_histo, list); break;
1320 case edfhWLDELTA: ret = do_cpte_real(xd, cptpEDFH, i, fflags, &dfhist->wl_delta, list); break;
1321 case edfhSUMWEIGHTS: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_weights, list); break;
1322 case edfhSUMDG: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_dg, list); break;
1323 case edfhSUMMINVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_minvar, list); break;
1324 case edfhSUMVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_variance, list); break;
1325 case edfhACCUMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p, list); break;
1326 case edfhACCUMM: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m, list); break;
1327 case edfhACCUMP2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p2, list); break;
1328 case edfhACCUMM2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m2, list); break;
1329 case edfhTIJ: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij, list); break;
1330 case edfhTIJEMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij_empirical, list); break;
1333 gmx_fatal(FARGS, "Unknown df history entry %d\n"
1334 "You are probably reading a new checkpoint file with old code", i);
1343 /* This function stores the last whole configuration of the reference and
1344 * average structure in the .cpt file
1346 static int do_cpt_EDstate(XDR *xd, gmx_bool bRead,
1347 edsamstate_t *EDstate, FILE *list)
1354 EDstate->bFromCpt = bRead;
1356 if (EDstate->nED <= 0)
1361 /* When reading, init_edsam has not been called yet,
1362 * so we have to allocate memory first. */
1365 snew(EDstate->nref, EDstate->nED);
1366 snew(EDstate->old_sref, EDstate->nED);
1367 snew(EDstate->nav, EDstate->nED);
1368 snew(EDstate->old_sav, EDstate->nED);
1371 /* Read/write the last whole conformation of SREF and SAV for each ED dataset (usually only one) */
1372 for (i = 0; i < EDstate->nED; i++)
1374 /* Reference structure SREF */
1375 sprintf(buf, "ED%d # of atoms in reference structure", i+1);
1376 do_cpt_int_err(xd, buf, &EDstate->nref[i], list);
1377 sprintf(buf, "ED%d x_ref", i+1);
1380 snew(EDstate->old_sref[i], EDstate->nref[i]);
1381 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref[i], list);
1385 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref_p[i], list);
1388 /* Average structure SAV */
1389 sprintf(buf, "ED%d # of atoms in average structure", i+1);
1390 do_cpt_int_err(xd, buf, &EDstate->nav[i], list);
1391 sprintf(buf, "ED%d x_av", i+1);
1394 snew(EDstate->old_sav[i], EDstate->nav[i]);
1395 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav[i], list);
1399 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav_p[i], list);
1407 static int do_cpt_files(XDR *xd, gmx_bool bRead,
1408 gmx_file_position_t **p_outputfiles, int *nfiles,
1409 FILE *list, int file_version)
1413 gmx_off_t mask = 0xFFFFFFFFL;
1414 int offset_high, offset_low;
1416 gmx_file_position_t *outputfiles;
1418 if (do_cpt_int(xd, "number of output files", nfiles, list) != 0)
1425 snew(*p_outputfiles, *nfiles);
1428 outputfiles = *p_outputfiles;
1430 for (i = 0; i < *nfiles; i++)
1432 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1435 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1436 strncpy(outputfiles[i].filename, buf, CPTSTRLEN-1);
1442 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1446 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1450 outputfiles[i].offset = ( ((gmx_off_t) offset_high) << 32 ) | ( (gmx_off_t) offset_low & mask );
1454 buf = outputfiles[i].filename;
1455 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1457 offset = outputfiles[i].offset;
1465 offset_low = (int) (offset & mask);
1466 offset_high = (int) ((offset >> 32) & mask);
1468 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1472 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1477 if (file_version >= 8)
1479 if (do_cpt_int(xd, "file_checksum_size", &(outputfiles[i].chksum_size),
1484 if (do_cpt_u_chars(xd, "file_checksum", 16, outputfiles[i].chksum, list) != 0)
1491 outputfiles[i].chksum_size = -1;
1498 void write_checkpoint(const char *fn, gmx_bool bNumberAndKeep,
1499 FILE *fplog, t_commrec *cr,
1500 int eIntegrator, int simulation_part,
1501 gmx_bool bExpanded, int elamstats,
1502 gmx_int64_t step, double t, t_state *state)
1512 char *fntemp; /* the temporary checkpoint file name */
1514 char timebuf[STRLEN];
1515 int nppnodes, npmenodes, flag_64bit;
1516 char buf[1024], suffix[5+STEPSTRSIZE], sbuf[STEPSTRSIZE];
1517 gmx_file_position_t *outputfiles;
1520 int flags_eks, flags_enh, flags_dfh, i;
1523 if (DOMAINDECOMP(cr))
1525 nppnodes = cr->dd->nnodes;
1526 npmenodes = cr->npmenodes;
1534 #ifndef GMX_NO_RENAME
1535 /* make the new temporary filename */
1536 snew(fntemp, strlen(fn)+5+STEPSTRSIZE);
1538 fntemp[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1539 sprintf(suffix, "_%s%s", "step", gmx_step_str(step, sbuf));
1540 strcat(fntemp, suffix);
1541 strcat(fntemp, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1543 /* if we can't rename, we just overwrite the cpt file.
1544 * dangerous if interrupted.
1546 snew(fntemp, strlen(fn));
1550 gmx_ctime_r(&now, timebuf, STRLEN);
1554 fprintf(fplog, "Writing checkpoint, step %s at %s\n\n",
1555 gmx_step_str(step, buf), timebuf);
1558 /* Get offsets for open files */
1559 gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1561 fp = gmx_fio_open(fntemp, "w");
1563 if (state->ekinstate.bUpToDate)
1566 ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) |
1567 (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) |
1568 (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1576 if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1578 flags_enh |= (1<<eenhENERGY_N);
1579 if (state->enerhist.nsum > 0)
1581 flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1582 (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSUM));
1584 if (state->enerhist.nsum_sim > 0)
1586 flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSTEPS_SIM) |
1587 (1<<eenhENERGY_NSUM_SIM));
1589 if (state->enerhist.dht)
1591 flags_enh |= ( (1<< eenhENERGY_DELTA_H_NN) |
1592 (1<< eenhENERGY_DELTA_H_LIST) |
1593 (1<< eenhENERGY_DELTA_H_STARTTIME) |
1594 (1<< eenhENERGY_DELTA_H_STARTLAMBDA) );
1600 flags_dfh = ((1<<edfhBEQUIL) | (1<<edfhNATLAMBDA) | (1<<edfhSUMWEIGHTS) | (1<<edfhSUMDG) |
1601 (1<<edfhTIJ) | (1<<edfhTIJEMP));
1604 flags_dfh |= ((1<<edfhWLDELTA) | (1<<edfhWLHISTO));
1606 if ((elamstats == elamstatsMINVAR) || (elamstats == elamstatsBARKER) || (elamstats == elamstatsMETROPOLIS))
1608 flags_dfh |= ((1<<edfhACCUMP) | (1<<edfhACCUMM) | (1<<edfhACCUMP2) | (1<<edfhACCUMM2)
1609 | (1<<edfhSUMMINVAR) | (1<<edfhSUMVAR));
1617 /* We can check many more things now (CPU, acceleration, etc), but
1618 * it is highly unlikely to have two separate builds with exactly
1619 * the same version, user, time, and build host!
1622 version = gmx_strdup(gmx_version());
1623 btime = gmx_strdup(BUILD_TIME);
1624 buser = gmx_strdup(BUILD_USER);
1625 bhost = gmx_strdup(BUILD_HOST);
1627 double_prec = GMX_CPT_BUILD_DP;
1628 fprog = gmx_strdup(Program());
1630 ftime = &(timebuf[0]);
1632 do_cpt_header(gmx_fio_getxdr(fp), FALSE, &file_version,
1633 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1634 &eIntegrator, &simulation_part, &step, &t, &nppnodes,
1635 DOMAINDECOMP(cr) ? cr->dd->nc : NULL, &npmenodes,
1636 &state->natoms, &state->ngtc, &state->nnhpres,
1637 &state->nhchainlength, &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
1638 &state->edsamstate.nED, &state->swapstate.eSwapCoords,
1647 if ((do_cpt_state(gmx_fio_getxdr(fp), FALSE, state->flags, state, NULL) < 0) ||
1648 (do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL) < 0) ||
1649 (do_cpt_enerhist(gmx_fio_getxdr(fp), FALSE, flags_enh, &state->enerhist, NULL) < 0) ||
1650 (do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL) < 0) ||
1651 (do_cpt_EDstate(gmx_fio_getxdr(fp), FALSE, &state->edsamstate, NULL) < 0) ||
1652 (do_cpt_swapstate(gmx_fio_getxdr(fp), FALSE, &state->swapstate, NULL) < 0) ||
1653 (do_cpt_files(gmx_fio_getxdr(fp), FALSE, &outputfiles, &noutputfiles, NULL,
1656 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1659 do_cpt_footer(gmx_fio_getxdr(fp), file_version);
1661 /* we really, REALLY, want to make sure to physically write the checkpoint,
1662 and all the files it depends on, out to disk. Because we've
1663 opened the checkpoint with gmx_fio_open(), it's in our list
1665 ret = gmx_fio_all_output_fsync();
1671 "Cannot fsync '%s'; maybe you are out of disk space?",
1672 gmx_fio_getname(ret));
1674 if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV) == NULL)
1684 if (gmx_fio_close(fp) != 0)
1686 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1689 /* we don't move the checkpoint if the user specified they didn't want it,
1690 or if the fsyncs failed */
1691 #ifndef GMX_NO_RENAME
1692 if (!bNumberAndKeep && !ret)
1696 /* Rename the previous checkpoint file */
1698 buf[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1699 strcat(buf, "_prev");
1700 strcat(buf, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1702 /* we copy here so that if something goes wrong between now and
1703 * the rename below, there's always a state.cpt.
1704 * If renames are atomic (such as in POSIX systems),
1705 * this copying should be unneccesary.
1707 gmx_file_copy(fn, buf, FALSE);
1708 /* We don't really care if this fails:
1709 * there's already a new checkpoint.
1712 gmx_file_rename(fn, buf);
1715 if (gmx_file_rename(fntemp, fn) != 0)
1717 gmx_file("Cannot rename checkpoint file; maybe you are out of disk space?");
1720 #endif /* GMX_NO_RENAME */
1726 /*code for alternate checkpointing scheme. moved from top of loop over
1728 fcRequestCheckPoint();
1729 if (fcCheckPointParallel( cr->nodeid, NULL, 0) == 0)
1731 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", step );
1733 #endif /* end GMX_FAHCORE block */
1736 static void print_flag_mismatch(FILE *fplog, int sflags, int fflags)
1740 fprintf(fplog, "\nState entry mismatch between the simulation and the checkpoint file\n");
1741 fprintf(fplog, "Entries which are not present in the checkpoint file will not be updated\n");
1742 fprintf(fplog, " %24s %11s %11s\n", "", "simulation", "checkpoint");
1743 for (i = 0; i < estNR; i++)
1745 if ((sflags & (1<<i)) || (fflags & (1<<i)))
1747 fprintf(fplog, " %24s %11s %11s\n",
1749 (sflags & (1<<i)) ? " present " : "not present",
1750 (fflags & (1<<i)) ? " present " : "not present");
1755 static void check_int(FILE *fplog, const char *type, int p, int f, gmx_bool *mm)
1757 FILE *fp = fplog ? fplog : stderr;
1761 fprintf(fp, " %s mismatch,\n", type);
1762 fprintf(fp, " current program: %d\n", p);
1763 fprintf(fp, " checkpoint file: %d\n", f);
1769 static void check_string(FILE *fplog, const char *type, const char *p,
1770 const char *f, gmx_bool *mm)
1772 FILE *fp = fplog ? fplog : stderr;
1774 if (strcmp(p, f) != 0)
1776 fprintf(fp, " %s mismatch,\n", type);
1777 fprintf(fp, " current program: %s\n", p);
1778 fprintf(fp, " checkpoint file: %s\n", f);
1784 static void check_match(FILE *fplog,
1786 char *btime, char *buser, char *bhost, int double_prec,
1788 t_commrec *cr, int npp_f, int npme_f,
1789 ivec dd_nc, ivec dd_nc_f)
1796 check_string(fplog, "Version", gmx_version(), version, &mm);
1797 check_string(fplog, "Build time", BUILD_TIME, btime, &mm);
1798 check_string(fplog, "Build user", BUILD_USER, buser, &mm);
1799 check_string(fplog, "Build host", BUILD_HOST, bhost, &mm);
1800 check_int (fplog, "Double prec.", GMX_CPT_BUILD_DP, double_prec, &mm);
1801 check_string(fplog, "Program name", Program(), fprog, &mm);
1803 check_int (fplog, "#nodes", cr->nnodes, npp_f+npme_f, &mm);
1806 check_int (fplog, "#PME-nodes", cr->npmenodes, npme_f, &mm);
1809 if (cr->npmenodes >= 0)
1811 npp -= cr->npmenodes;
1815 check_int (fplog, "#DD-cells[x]", dd_nc[XX], dd_nc_f[XX], &mm);
1816 check_int (fplog, "#DD-cells[y]", dd_nc[YY], dd_nc_f[YY], &mm);
1817 check_int (fplog, "#DD-cells[z]", dd_nc[ZZ], dd_nc_f[ZZ], &mm);
1824 "Gromacs binary or parallel settings not identical to previous run.\n"
1825 "Continuation is exact, but is not guaranteed to be binary identical%s.\n\n",
1826 fplog ? ",\n see the log file for details" : "");
1831 "Gromacs binary or parallel settings not identical to previous run.\n"
1832 "Continuation is exact, but is not guaranteed to be binary identical.\n\n");
1837 static void read_checkpoint(const char *fn, FILE **pfplog,
1838 t_commrec *cr, ivec dd_nc,
1839 int eIntegrator, int *init_fep_state, gmx_int64_t *step, double *t,
1840 t_state *state, gmx_bool *bReadEkin,
1841 int *simulation_part,
1842 gmx_bool bAppendOutputFiles, gmx_bool bForceAppend)
1847 char *version, *btime, *buser, *bhost, *fprog, *ftime;
1849 char filename[STRLEN], buf[STEPSTRSIZE];
1850 int nppnodes, eIntegrator_f, nppnodes_f, npmenodes_f;
1852 int natoms, ngtc, nnhpres, nhchainlength, nlambda, fflags, flags_eks, flags_enh, flags_dfh;
1855 gmx_file_position_t *outputfiles;
1857 t_fileio *chksum_file;
1858 FILE * fplog = *pfplog;
1859 unsigned char digest[16];
1860 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1861 struct flock fl; /* don't initialize here: the struct order is OS
1865 const char *int_warn =
1866 "WARNING: The checkpoint file was generated with integrator %s,\n"
1867 " while the simulation uses integrator %s\n\n";
1869 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1870 fl.l_type = F_WRLCK;
1871 fl.l_whence = SEEK_SET;
1877 fp = gmx_fio_open(fn, "r");
1878 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
1879 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1880 &eIntegrator_f, simulation_part, step, t,
1881 &nppnodes_f, dd_nc_f, &npmenodes_f,
1882 &natoms, &ngtc, &nnhpres, &nhchainlength, &nlambda,
1883 &fflags, &flags_eks, &flags_enh, &flags_dfh,
1884 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
1886 if (bAppendOutputFiles &&
1887 file_version >= 13 && double_prec != GMX_CPT_BUILD_DP)
1889 gmx_fatal(FARGS, "Output file appending requested, but the code and checkpoint file precision (single/double) don't match");
1892 if (cr == NULL || MASTER(cr))
1894 fprintf(stderr, "\nReading checkpoint file %s generated: %s\n\n",
1898 /* This will not be written if we do appending, since fplog is still NULL then */
1901 fprintf(fplog, "\n");
1902 fprintf(fplog, "Reading checkpoint file %s\n", fn);
1903 fprintf(fplog, " file generated by: %s\n", fprog);
1904 fprintf(fplog, " file generated at: %s\n", ftime);
1905 fprintf(fplog, " GROMACS build time: %s\n", btime);
1906 fprintf(fplog, " GROMACS build user: %s\n", buser);
1907 fprintf(fplog, " GROMACS build host: %s\n", bhost);
1908 fprintf(fplog, " GROMACS double prec.: %d\n", double_prec);
1909 fprintf(fplog, " simulation part #: %d\n", *simulation_part);
1910 fprintf(fplog, " step: %s\n", gmx_step_str(*step, buf));
1911 fprintf(fplog, " time: %f\n", *t);
1912 fprintf(fplog, "\n");
1915 if (natoms != state->natoms)
1917 gmx_fatal(FARGS, "Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms", natoms, state->natoms);
1919 if (ngtc != state->ngtc)
1921 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);
1923 if (nnhpres != state->nnhpres)
1925 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);
1928 if (nlambda != state->dfhist.nlambda)
1930 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);
1933 init_gtc_state(state, state->ngtc, state->nnhpres, nhchainlength); /* need to keep this here to keep the tpr format working */
1934 /* write over whatever was read; we use the number of Nose-Hoover chains from the checkpoint */
1936 if (eIntegrator_f != eIntegrator)
1940 fprintf(stderr, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1942 if (bAppendOutputFiles)
1945 "Output file appending requested, but input/checkpoint integrators do not match.\n"
1946 "Stopping the run to prevent you from ruining all your data...\n"
1947 "If you _really_ know what you are doing, try with the -noappend option.\n");
1951 fprintf(fplog, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1960 else if (cr->nnodes == nppnodes_f + npmenodes_f)
1962 if (cr->npmenodes < 0)
1964 cr->npmenodes = npmenodes_f;
1966 nppnodes = cr->nnodes - cr->npmenodes;
1967 if (nppnodes == nppnodes_f)
1969 for (d = 0; d < DIM; d++)
1973 dd_nc[d] = dd_nc_f[d];
1980 /* The number of PP nodes has not been set yet */
1984 if (fflags != state->flags)
1989 if (bAppendOutputFiles)
1992 "Output file appending requested, but input and checkpoint states are not identical.\n"
1993 "Stopping the run to prevent you from ruining all your data...\n"
1994 "You can try with the -noappend option, and get more info in the log file.\n");
1997 if (getenv("GMX_ALLOW_CPT_MISMATCH") == NULL)
1999 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");
2004 "WARNING: The checkpoint state entries do not match the simulation,\n"
2005 " see the log file for details\n\n");
2011 print_flag_mismatch(fplog, state->flags, fflags);
2018 check_match(fplog, version, btime, buser, bhost, double_prec, fprog,
2019 cr, nppnodes_f, npmenodes_f, dd_nc, dd_nc_f);
2022 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, fflags, state, NULL);
2023 *init_fep_state = state->fep_state; /* there should be a better way to do this than setting it here.
2024 Investigate for 5.0. */
2029 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2034 *bReadEkin = ((flags_eks & (1<<eeksEKINH)) || (flags_eks & (1<<eeksEKINF)) || (flags_eks & (1<<eeksEKINO)) ||
2035 ((flags_eks & (1<<eeksEKINSCALEF)) | (flags_eks & (1<<eeksEKINSCALEH)) | (flags_eks & (1<<eeksVSCALE))));
2037 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2038 flags_enh, &state->enerhist, NULL);
2044 if (file_version < 6)
2046 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.";
2048 fprintf(stderr, "\nWARNING: %s\n\n", warn);
2051 fprintf(fplog, "\nWARNING: %s\n\n", warn);
2053 state->enerhist.nsum = *step;
2054 state->enerhist.nsum_sim = *step;
2057 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2063 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2069 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2075 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, NULL, file_version);
2081 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2086 if (gmx_fio_close(fp) != 0)
2088 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2097 /* If the user wants to append to output files,
2098 * we use the file pointer positions of the output files stored
2099 * in the checkpoint file and truncate the files such that any frames
2100 * written after the checkpoint time are removed.
2101 * All files are md5sum checked such that we can be sure that
2102 * we do not truncate other (maybe imprortant) files.
2104 if (bAppendOutputFiles)
2106 if (fn2ftp(outputfiles[0].filename) != efLOG)
2108 /* make sure first file is log file so that it is OK to use it for
2111 gmx_fatal(FARGS, "The first output file should always be the log "
2112 "file but instead is: %s. Cannot do appending because of this condition.", outputfiles[0].filename);
2114 for (i = 0; i < nfiles; i++)
2116 if (outputfiles[i].offset < 0)
2118 gmx_fatal(FARGS, "The original run wrote a file called '%s' which "
2119 "is larger than 2 GB, but mdrun did not support large file"
2120 " offsets. Can not append. Run mdrun with -noappend",
2121 outputfiles[i].filename);
2124 chksum_file = gmx_fio_open(outputfiles[i].filename, "a");
2127 chksum_file = gmx_fio_open(outputfiles[i].filename, "r+");
2132 /* Note that there are systems where the lock operation
2133 * will succeed, but a second process can also lock the file.
2134 * We should probably try to detect this.
2136 #if defined __native_client__
2140 #elif defined GMX_NATIVE_WINDOWS
2141 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX) == -1)
2143 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl) == -1)
2146 if (errno == ENOSYS)
2150 gmx_fatal(FARGS, "File locking is not supported on this system. Use -noappend or specify -append explicitly to append anyhow.");
2154 fprintf(stderr, "\nNOTE: File locking is not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2157 fprintf(fplog, "\nNOTE: File locking not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2161 else if (errno == EACCES || errno == EAGAIN)
2163 gmx_fatal(FARGS, "Failed to lock: %s. Already running "
2164 "simulation?", outputfiles[i].filename);
2168 gmx_fatal(FARGS, "Failed to lock: %s. %s.",
2169 outputfiles[i].filename, strerror(errno));
2174 /* compute md5 chksum */
2175 if (outputfiles[i].chksum_size != -1)
2177 if (gmx_fio_get_file_md5(chksum_file, outputfiles[i].offset,
2178 digest) != outputfiles[i].chksum_size) /*at the end of the call the file position is at the end of the file*/
2180 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.",
2181 outputfiles[i].chksum_size,
2182 outputfiles[i].filename);
2185 if (i == 0) /*log file needs to be seeked in case we need to truncate (other files are truncated below)*/
2187 if (gmx_fio_seek(chksum_file, outputfiles[i].offset))
2189 gmx_fatal(FARGS, "Seek error! Failed to truncate log-file: %s.", strerror(errno));
2194 if (i == 0) /*open log file here - so that lock is never lifted
2195 after chksum is calculated */
2197 *pfplog = gmx_fio_getfp(chksum_file);
2201 gmx_fio_close(chksum_file);
2204 /* compare md5 chksum */
2205 if (outputfiles[i].chksum_size != -1 &&
2206 memcmp(digest, outputfiles[i].chksum, 16) != 0)
2210 fprintf(debug, "chksum for %s: ", outputfiles[i].filename);
2211 for (j = 0; j < 16; j++)
2213 fprintf(debug, "%02x", digest[j]);
2215 fprintf(debug, "\n");
2217 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.",
2218 outputfiles[i].filename);
2223 if (i != 0) /*log file is already seeked to correct position */
2225 #ifdef GMX_NATIVE_WINDOWS
2226 rc = gmx_wintruncate(outputfiles[i].filename, outputfiles[i].offset);
2228 rc = truncate(outputfiles[i].filename, outputfiles[i].offset);
2232 gmx_fatal(FARGS, "Truncation of file %s failed. Cannot do appending because of this failure.", outputfiles[i].filename);
2242 void load_checkpoint(const char *fn, FILE **fplog,
2243 t_commrec *cr, ivec dd_nc,
2244 t_inputrec *ir, t_state *state,
2245 gmx_bool *bReadEkin,
2246 gmx_bool bAppend, gmx_bool bForceAppend)
2253 /* Read the state from the checkpoint file */
2254 read_checkpoint(fn, fplog,
2256 ir->eI, &(ir->fepvals->init_fep_state), &step, &t, state, bReadEkin,
2257 &ir->simulation_part, bAppend, bForceAppend);
2261 gmx_bcast(sizeof(cr->npmenodes), &cr->npmenodes, cr);
2262 gmx_bcast(DIM*sizeof(dd_nc[0]), dd_nc, cr);
2263 gmx_bcast(sizeof(step), &step, cr);
2264 gmx_bcast(sizeof(*bReadEkin), bReadEkin, cr);
2266 ir->bContinuation = TRUE;
2267 if (ir->nsteps >= 0)
2269 ir->nsteps += ir->init_step - step;
2271 ir->init_step = step;
2272 ir->simulation_part += 1;
2275 static void read_checkpoint_data(t_fileio *fp, int *simulation_part,
2276 gmx_int64_t *step, double *t, t_state *state,
2277 int *nfiles, gmx_file_position_t **outputfiles)
2280 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2285 int flags_eks, flags_enh, flags_dfh;
2287 gmx_file_position_t *files_loc = NULL;
2290 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2291 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2292 &eIntegrator, simulation_part, step, t, &nppnodes, dd_nc, &npme,
2293 &state->natoms, &state->ngtc, &state->nnhpres, &state->nhchainlength,
2294 &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
2295 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
2297 do_cpt_state(gmx_fio_getxdr(fp), TRUE, state->flags, state, NULL);
2302 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2307 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2308 flags_enh, &state->enerhist, NULL);
2313 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2319 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2325 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2331 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE,
2332 outputfiles != NULL ? outputfiles : &files_loc,
2333 outputfiles != NULL ? nfiles : &nfiles_loc,
2334 NULL, file_version);
2335 if (files_loc != NULL)
2345 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2359 read_checkpoint_state(const char *fn, int *simulation_part,
2360 gmx_int64_t *step, double *t, t_state *state)
2364 fp = gmx_fio_open(fn, "r");
2365 read_checkpoint_data(fp, simulation_part, step, t, state, NULL, NULL);
2366 if (gmx_fio_close(fp) != 0)
2368 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2372 void read_checkpoint_trxframe(t_fileio *fp, t_trxframe *fr)
2374 /* This next line is nasty because the sub-structures of t_state
2375 * cannot be assumed to be zeroed (or even initialized in ways the
2376 * rest of the code might assume). Using snew would be better, but
2377 * this will all go away for 5.0. */
2379 int simulation_part;
2383 init_state(&state, 0, 0, 0, 0, 0);
2385 read_checkpoint_data(fp, &simulation_part, &step, &t, &state, NULL, NULL);
2387 fr->natoms = state.natoms;
2390 fr->step = gmx_int64_to_int(step,
2391 "conversion of checkpoint to trajectory");
2395 fr->lambda = state.lambda[efptFEP];
2396 fr->fep_state = state.fep_state;
2398 fr->bX = (state.flags & (1<<estX));
2404 fr->bV = (state.flags & (1<<estV));
2411 fr->bBox = (state.flags & (1<<estBOX));
2414 copy_mat(state.box, fr->box);
2419 void list_checkpoint(const char *fn, FILE *out)
2423 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2425 int eIntegrator, simulation_part, nppnodes, npme;
2430 int flags_eks, flags_enh, flags_dfh;
2434 gmx_file_position_t *outputfiles;
2437 init_state(&state, -1, -1, -1, -1, 0);
2439 fp = gmx_fio_open(fn, "r");
2440 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2441 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2442 &eIntegrator, &simulation_part, &step, &t, &nppnodes, dd_nc, &npme,
2443 &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
2444 &(state.dfhist.nlambda), &state.flags,
2445 &flags_eks, &flags_enh, &flags_dfh, &state.edsamstate.nED,
2446 &state.swapstate.eSwapCoords, out);
2447 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, state.flags, &state, out);
2452 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state.ekinstate, out);
2457 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2458 flags_enh, &state.enerhist, out);
2462 ret = do_cpt_df_hist(gmx_fio_getxdr(fp),
2463 flags_dfh, &state.dfhist, out);
2468 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state.edsamstate, out);
2473 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state.swapstate, out);
2478 do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, out, file_version);
2483 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2490 if (gmx_fio_close(fp) != 0)
2492 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2499 static gmx_bool exist_output_file(const char *fnm_cp, int nfile, const t_filenm fnm[])
2503 /* Check if the output file name stored in the checkpoint file
2504 * is one of the output file names of mdrun.
2508 !(is_output(&fnm[i]) && strcmp(fnm_cp, fnm[i].fns[0]) == 0))
2513 return (i < nfile && gmx_fexist(fnm_cp));
2516 /* This routine cannot print tons of data, since it is called before the log file is opened. */
2517 gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
2518 gmx_int64_t *cpt_step, t_commrec *cr,
2519 gmx_bool bAppendReq,
2520 int nfile, const t_filenm fnm[],
2521 const char *part_suffix, gmx_bool *bAddPart)
2524 gmx_int64_t step = 0;
2526 /* This next line is nasty because the sub-structures of t_state
2527 * cannot be assumed to be zeroed (or even initialized in ways the
2528 * rest of the code might assume). Using snew would be better, but
2529 * this will all go away for 5.0. */
2532 gmx_file_position_t *outputfiles;
2535 char *fn, suf_up[STRLEN];
2541 if (!gmx_fexist(filename) || (!(fp = gmx_fio_open(filename, "r")) ))
2543 *simulation_part = 0;
2547 init_state(&state, 0, 0, 0, 0, 0);
2549 read_checkpoint_data(fp, simulation_part, &step, &t, &state,
2550 &nfiles, &outputfiles);
2551 if (gmx_fio_close(fp) != 0)
2553 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2560 for (f = 0; f < nfiles; f++)
2562 if (exist_output_file(outputfiles[f].filename, nfile, fnm))
2567 if (nexist == nfiles)
2569 bAppend = bAppendReq;
2571 else if (nexist > 0)
2574 "Output file appending has been requested,\n"
2575 "but some output files listed in the checkpoint file %s\n"
2576 "are not present or are named differently by the current program:\n",
2578 fprintf(stderr, "output files present:");
2579 for (f = 0; f < nfiles; f++)
2581 if (exist_output_file(outputfiles[f].filename,
2584 fprintf(stderr, " %s", outputfiles[f].filename);
2587 fprintf(stderr, "\n");
2588 fprintf(stderr, "output files not present or named differently:");
2589 for (f = 0; f < nfiles; f++)
2591 if (!exist_output_file(outputfiles[f].filename,
2594 fprintf(stderr, " %s", outputfiles[f].filename);
2597 fprintf(stderr, "\n");
2599 gmx_fatal(FARGS, "File appending requested, but only %d of the %d output files are present", nexist, nfiles);
2607 gmx_fatal(FARGS, "File appending requested, but no output file information is stored in the checkpoint file");
2609 fn = outputfiles[0].filename;
2610 if (strlen(fn) < 4 ||
2611 gmx_strcasecmp(fn+strlen(fn)-4, ftp2ext(efLOG)) == 0)
2613 gmx_fatal(FARGS, "File appending requested, but the log file is not the first file listed in the checkpoint file");
2615 /* Set bAddPart to whether the suffix string '.part' is present
2616 * in the log file name.
2618 strcpy(suf_up, part_suffix);
2620 *bAddPart = (strstr(fn, part_suffix) != NULL ||
2621 strstr(fn, suf_up) != NULL);
2629 gmx_bcast(sizeof(*simulation_part), simulation_part, cr);
2631 if (*simulation_part > 0 && bAppendReq)
2633 gmx_bcast(sizeof(bAppend), &bAppend, cr);
2634 gmx_bcast(sizeof(*bAddPart), bAddPart, cr);
2637 if (NULL != cpt_step)