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. */
47 #ifdef HAVE_SYS_TIME_H
55 #ifdef GMX_NATIVE_WINDOWS
58 #include <sys/locking.h>
68 #include "gromacs/random/random.h"
69 #include "checkpoint.h"
74 #include "gromacs/fileio/filenm.h"
75 #include "gromacs/fileio/futil.h"
76 #include "gromacs/fileio/gmxfio.h"
77 #include "gromacs/fileio/xdrf.h"
78 #include "gromacs/fileio/xdr_datatype.h"
80 #include "buildinfo.h"
86 #define CPT_MAGIC1 171817
87 #define CPT_MAGIC2 171819
88 #define CPTSTRLEN 1024
91 #define GMX_CPT_BUILD_DP 1
93 #define GMX_CPT_BUILD_DP 0
96 /* cpt_version should normally only be changed
97 * when the header of footer format changes.
98 * The state data format itself is backward and forward compatible.
99 * But old code can not read a new entry that is present in the file
100 * (but can read a new format when new entries are not present).
102 static const int cpt_version = 16;
105 const char *est_names[estNR] =
108 "box", "box-rel", "box-v", "pres_prev",
109 "nosehoover-xi", "thermostat-integral",
110 "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
111 "disre_initf", "disre_rm3tav",
112 "orire_initf", "orire_Dtav",
113 "svir_prev", "nosehoover-vxi", "v_eta", "vol0", "nhpres_xi", "nhpres_vxi", "fvir_prev", "fep_state", "MC-rng", "MC-rng-i"
117 eeksEKIN_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksEKINF, eeksEKINO, eeksEKINSCALEF, eeksEKINSCALEH, eeksVSCALE, eeksEKINTOTAL, eeksNR
120 const char *eeks_names[eeksNR] =
122 "Ekin_n", "Ekinh", "dEkindlambda", "mv_cos",
123 "Ekinf", "Ekinh_old", "EkinScaleF_NHC", "EkinScaleH_NHC", "Vscale_NHC", "Ekin_Total"
127 eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
128 eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
129 eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM,
130 eenhENERGY_DELTA_H_NN,
131 eenhENERGY_DELTA_H_LIST,
132 eenhENERGY_DELTA_H_STARTTIME,
133 eenhENERGY_DELTA_H_STARTLAMBDA,
137 const char *eenh_names[eenhNR] =
139 "energy_n", "energy_aver", "energy_sum", "energy_nsum",
140 "energy_sum_sim", "energy_nsum_sim",
141 "energy_nsteps", "energy_nsteps_sim",
143 "energy_delta_h_list",
144 "energy_delta_h_start_time",
145 "energy_delta_h_start_lambda"
148 /* free energy history variables -- need to be preserved over checkpoint */
150 edfhBEQUIL, edfhNATLAMBDA, edfhWLHISTO, edfhWLDELTA, edfhSUMWEIGHTS, edfhSUMDG, edfhSUMMINVAR, edfhSUMVAR,
151 edfhACCUMP, edfhACCUMM, edfhACCUMP2, edfhACCUMM2, edfhTIJ, edfhTIJEMP, edfhNR
153 /* free energy history variable names */
154 const char *edfh_names[edfhNR] =
156 "bEquilibrated", "N_at_state", "Wang-Landau Histogram", "Wang-Landau Delta", "Weights", "Free Energies", "minvar", "variance",
157 "accumulated_plus", "accumulated_minus", "accumulated_plus_2", "accumulated_minus_2", "Tij", "Tij_empirical"
160 #ifdef GMX_NATIVE_WINDOWS
162 gmx_wintruncate(const char *filename, __int64 size)
165 /*we do this elsewhere*/
171 fp = fopen(filename, "rb+");
178 return _chsize_s( fileno(fp), size);
185 ecprREAL, ecprRVEC, ecprMATRIX
189 cptpEST, cptpEEKS, cptpEENH, cptpEDFH
191 /* enums for the different components of checkpoint variables, replacing the hard coded ones.
192 cptpEST - state variables.
193 cptpEEKS - Kinetic energy state variables.
194 cptpEENH - Energy history state variables.
195 cptpEDFH - free energy history variables.
199 static const char *st_names(int cptp, int ecpt)
203 case cptpEST: return est_names [ecpt]; break;
204 case cptpEEKS: return eeks_names[ecpt]; break;
205 case cptpEENH: return eenh_names[ecpt]; break;
206 case cptpEDFH: return edfh_names[ecpt]; break;
212 static void cp_warning(FILE *fp)
214 fprintf(fp, "\nWARNING: Checkpoint file is corrupted or truncated\n\n");
217 static void cp_error()
219 gmx_fatal(FARGS, "Checkpoint file corrupted/truncated, or maybe you are out of disk space?");
222 static void do_cpt_string_err(XDR *xd, gmx_bool bRead, const char *desc, char **s, FILE *list)
230 res = xdr_string(xd, s, CPTSTRLEN);
237 fprintf(list, "%s = %s\n", desc, *s);
242 static int do_cpt_int(XDR *xd, const char *desc, int *i, FILE *list)
246 res = xdr_int(xd, i);
253 fprintf(list, "%s = %d\n", desc, *i);
258 static int do_cpt_u_chars(XDR *xd, const char *desc, int n, unsigned char *i, FILE *list)
264 fprintf(list, "%s = ", desc);
266 for (j = 0; j < n && res; j++)
268 res &= xdr_u_char(xd, &i[j]);
271 fprintf(list, "%02x", i[j]);
286 static void do_cpt_int_err(XDR *xd, const char *desc, int *i, FILE *list)
288 if (do_cpt_int(xd, desc, i, list) < 0)
294 static void do_cpt_step_err(XDR *xd, const char *desc, gmx_int64_t *i, FILE *list)
297 char buf[STEPSTRSIZE];
299 res = xdr_int64(xd, i);
306 fprintf(list, "%s = %s\n", desc, gmx_step_str(*i, buf));
310 static void do_cpt_double_err(XDR *xd, const char *desc, double *f, FILE *list)
314 res = xdr_double(xd, f);
321 fprintf(list, "%s = %f\n", desc, *f);
325 static void do_cpt_real_err(XDR *xd, real *f)
330 res = xdr_double(xd, f);
332 res = xdr_float(xd, f);
340 static void do_cpt_n_rvecs_err(XDR *xd, const char *desc, int n, rvec f[], FILE *list)
344 for (i = 0; i < n; i++)
346 for (j = 0; j < DIM; j++)
348 do_cpt_real_err(xd, &f[i][j]);
354 pr_rvecs(list, 0, desc, f, n);
358 /* If nval >= 0, nval is used; on read this should match the passed value.
359 * If nval n<0, *nptr is used; on read the value is stored in nptr
361 static int do_cpte_reals_low(XDR *xd, int cptp, int ecpt, int sflags,
362 int nval, int *nptr, real **v,
363 FILE *list, int erealtype)
367 int dtc = xdr_datatype_float;
369 int dtc = xdr_datatype_double;
371 real *vp, *va = NULL;
386 gmx_incons("*ntpr=NULL in do_cpte_reals_low");
391 res = xdr_int(xd, &nf);
402 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), nval, nf);
411 res = xdr_int(xd, &dt);
418 fprintf(stderr, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
419 st_names(cptp, ecpt), xdr_datatype_names[dtc],
420 xdr_datatype_names[dt]);
422 if (list || !(sflags & (1<<ecpt)))
435 if (dt == xdr_datatype_float)
437 if (dtc == xdr_datatype_float)
445 res = xdr_vector(xd, (char *)vf, nf,
446 (unsigned int)sizeof(float), (xdrproc_t)xdr_float);
451 if (dtc != xdr_datatype_float)
453 for (i = 0; i < nf; i++)
462 if (dtc == xdr_datatype_double)
470 res = xdr_vector(xd, (char *)vd, nf,
471 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
476 if (dtc != xdr_datatype_double)
478 for (i = 0; i < nf; i++)
491 pr_reals(list, 0, st_names(cptp, ecpt), vp, nf);
494 pr_rvecs(list, 0, st_names(cptp, ecpt), (rvec *)vp, nf/3);
497 gmx_incons("Unknown checkpoint real type");
509 /* This function stores n along with the reals for reading,
510 * but on reading it assumes that n matches the value in the checkpoint file,
511 * a fatal error is generated when this is not the case.
513 static int do_cpte_reals(XDR *xd, int cptp, int ecpt, int sflags,
514 int n, real **v, FILE *list)
516 return do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, v, list, ecprREAL);
519 /* This function does the same as do_cpte_reals,
520 * except that on reading it ignores the passed value of *n
521 * and stored the value read from the checkpoint file in *n.
523 static int do_cpte_n_reals(XDR *xd, int cptp, int ecpt, int sflags,
524 int *n, real **v, FILE *list)
526 return do_cpte_reals_low(xd, cptp, ecpt, sflags, -1, n, v, list, ecprREAL);
529 static int do_cpte_real(XDR *xd, int cptp, int ecpt, int sflags,
534 return do_cpte_reals_low(xd, cptp, ecpt, sflags, 1, NULL, &r, list, ecprREAL);
537 static int do_cpte_ints(XDR *xd, int cptp, int ecpt, int sflags,
538 int n, int **v, FILE *list)
541 int dtc = xdr_datatype_int;
546 res = xdr_int(xd, &nf);
551 if (list == NULL && v != NULL && nf != n)
553 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
556 res = xdr_int(xd, &dt);
563 gmx_fatal(FARGS, "Type mismatch for state entry %s, code type is %s, file type is %s\n",
564 st_names(cptp, ecpt), xdr_datatype_names[dtc],
565 xdr_datatype_names[dt]);
567 if (list || !(sflags & (1<<ecpt)) || v == NULL)
580 res = xdr_vector(xd, (char *)vp, nf,
581 (unsigned int)sizeof(int), (xdrproc_t)xdr_int);
588 pr_ivec(list, 0, st_names(cptp, ecpt), vp, nf, TRUE);
598 static int do_cpte_int(XDR *xd, int cptp, int ecpt, int sflags,
601 return do_cpte_ints(xd, cptp, ecpt, sflags, 1, &i, list);
604 static int do_cpte_doubles(XDR *xd, int cptp, int ecpt, int sflags,
605 int n, double **v, FILE *list)
608 int dtc = xdr_datatype_double;
609 double *vp, *va = NULL;
613 res = xdr_int(xd, &nf);
618 if (list == NULL && nf != n)
620 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
623 res = xdr_int(xd, &dt);
630 gmx_fatal(FARGS, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
631 st_names(cptp, ecpt), xdr_datatype_names[dtc],
632 xdr_datatype_names[dt]);
634 if (list || !(sflags & (1<<ecpt)))
647 res = xdr_vector(xd, (char *)vp, nf,
648 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
655 pr_doubles(list, 0, st_names(cptp, ecpt), vp, nf);
665 static int do_cpte_double(XDR *xd, int cptp, int ecpt, int sflags,
666 double *r, FILE *list)
668 return do_cpte_doubles(xd, cptp, ecpt, sflags, 1, &r, list);
672 static int do_cpte_rvecs(XDR *xd, int cptp, int ecpt, int sflags,
673 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)
706 char name[CPTSTRLEN];
713 for (i = 0; i < n; i++)
717 reti = do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, &(v[i]), NULL, ecprREAL);
718 if (list && reti == 0)
720 sprintf(name, "%s[%d]", st_names(cptp, ecpt), i);
721 pr_reals(list, 0, name, v[i], n);
731 static int do_cpte_matrices(XDR *xd, int cptp, int ecpt, int sflags,
732 int n, matrix **v, FILE *list)
735 matrix *vp, *va = NULL;
741 res = xdr_int(xd, &nf);
746 if (list == NULL && nf != n)
748 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
750 if (list || !(sflags & (1<<ecpt)))
763 snew(vr, nf*DIM*DIM);
764 for (i = 0; i < nf; i++)
766 for (j = 0; j < DIM; j++)
768 for (k = 0; k < DIM; k++)
770 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
774 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
775 nf*DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
776 for (i = 0; i < nf; i++)
778 for (j = 0; j < DIM; j++)
780 for (k = 0; k < DIM; k++)
782 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
788 if (list && ret == 0)
790 for (i = 0; i < nf; i++)
792 pr_rvecs(list, 0, st_names(cptp, ecpt), vp[i], DIM);
803 static void do_cpt_header(XDR *xd, gmx_bool bRead, int *file_version,
804 char **version, char **btime, char **buser, char **bhost,
806 char **fprog, char **ftime,
807 int *eIntegrator, int *simulation_part,
808 gmx_int64_t *step, double *t,
809 int *nnodes, int *dd_nc, int *npme,
810 int *natoms, int *ngtc, int *nnhpres, int *nhchainlength,
811 int *nlambda, int *flags_state,
812 int *flags_eks, int *flags_enh, int *flags_dfh,
813 int *nED, int *eSwapCoords,
830 res = xdr_int(xd, &magic);
833 gmx_fatal(FARGS, "The checkpoint file is empty/corrupted, or maybe you are out of disk space?");
835 if (magic != CPT_MAGIC1)
837 gmx_fatal(FARGS, "Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
838 "The checkpoint file is corrupted or not a checkpoint file",
844 gmx_gethostname(fhost, 255);
846 do_cpt_string_err(xd, bRead, "GROMACS version", version, list);
847 do_cpt_string_err(xd, bRead, "GROMACS build time", btime, list);
848 do_cpt_string_err(xd, bRead, "GROMACS build user", buser, list);
849 do_cpt_string_err(xd, bRead, "GROMACS build host", bhost, list);
850 do_cpt_string_err(xd, bRead, "generating program", fprog, list);
851 do_cpt_string_err(xd, bRead, "generation time", ftime, list);
852 *file_version = cpt_version;
853 do_cpt_int_err(xd, "checkpoint file version", file_version, list);
854 if (*file_version > cpt_version)
856 gmx_fatal(FARGS, "Attempting to read a checkpoint file of version %d with code of version %d\n", *file_version, cpt_version);
858 if (*file_version >= 13)
860 do_cpt_int_err(xd, "GROMACS double precision", double_prec, list);
866 if (*file_version >= 12)
868 do_cpt_string_err(xd, bRead, "generating host", &fhost, list);
874 do_cpt_int_err(xd, "#atoms", natoms, list);
875 do_cpt_int_err(xd, "#T-coupling groups", ngtc, list);
876 if (*file_version >= 10)
878 do_cpt_int_err(xd, "#Nose-Hoover T-chains", nhchainlength, list);
884 if (*file_version >= 11)
886 do_cpt_int_err(xd, "#Nose-Hoover T-chains for barostat ", nnhpres, list);
892 if (*file_version >= 14)
894 do_cpt_int_err(xd, "# of total lambda states ", nlambda, list);
900 do_cpt_int_err(xd, "integrator", eIntegrator, list);
901 if (*file_version >= 3)
903 do_cpt_int_err(xd, "simulation part #", simulation_part, list);
907 *simulation_part = 1;
909 if (*file_version >= 5)
911 do_cpt_step_err(xd, "step", step, list);
915 do_cpt_int_err(xd, "step", &idum, list);
918 do_cpt_double_err(xd, "t", t, list);
919 do_cpt_int_err(xd, "#PP-nodes", nnodes, list);
921 do_cpt_int_err(xd, "dd_nc[x]", dd_nc ? &(dd_nc[0]) : &idum, list);
922 do_cpt_int_err(xd, "dd_nc[y]", dd_nc ? &(dd_nc[1]) : &idum, list);
923 do_cpt_int_err(xd, "dd_nc[z]", dd_nc ? &(dd_nc[2]) : &idum, list);
924 do_cpt_int_err(xd, "#PME-only nodes", npme, list);
925 do_cpt_int_err(xd, "state flags", flags_state, list);
926 if (*file_version >= 4)
928 do_cpt_int_err(xd, "ekin data flags", flags_eks, list);
929 do_cpt_int_err(xd, "energy history flags", flags_enh, list);
934 *flags_enh = (*flags_state >> (estORIRE_DTAV+1));
935 *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
936 (1<<(estORIRE_DTAV+2)) |
937 (1<<(estORIRE_DTAV+3))));
939 if (*file_version >= 14)
941 do_cpt_int_err(xd, "df history flags", flags_dfh, list);
948 if (*file_version >= 15)
950 do_cpt_int_err(xd, "ED data sets", nED, list);
956 if (*file_version >= 16)
958 do_cpt_int_err(xd, "swap", eSwapCoords, list);
962 static int do_cpt_footer(XDR *xd, int file_version)
967 if (file_version >= 2)
970 res = xdr_int(xd, &magic);
975 if (magic != CPT_MAGIC2)
984 static int do_cpt_state(XDR *xd, gmx_bool bRead,
985 int fflags, t_state *state,
986 gmx_bool bReadRNG, FILE *list)
989 int **rng_p, **rngi_p;
996 nnht = state->nhchainlength*state->ngtc;
997 nnhtp = state->nhchainlength*state->nnhpres;
1001 rng_p = (int **)&state->ld_rng;
1002 rngi_p = &state->ld_rngi;
1006 /* Do not read the RNG data */
1011 if (bRead) /* we need to allocate space for dfhist if we are reading */
1013 init_df_history(&state->dfhist, state->dfhist.nlambda);
1016 /* We want the MC_RNG the same across all the notes for now -- lambda MC is global */
1018 sflags = state->flags;
1019 for (i = 0; (i < estNR && ret == 0); i++)
1021 if (fflags & (1<<i))
1025 case estLAMBDA: ret = do_cpte_reals(xd, cptpEST, i, sflags, efptNR, &(state->lambda), list); break;
1026 case estFEPSTATE: ret = do_cpte_int (xd, cptpEST, i, sflags, &state->fep_state, list); break;
1027 case estBOX: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box, list); break;
1028 case estBOX_REL: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box_rel, list); break;
1029 case estBOXV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->boxv, list); break;
1030 case estPRES_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->pres_prev, list); break;
1031 case estSVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->svir_prev, list); break;
1032 case estFVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->fvir_prev, list); break;
1033 case estNH_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_xi, list); break;
1034 case estNH_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_vxi, list); break;
1035 case estNHPRES_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_xi, list); break;
1036 case estNHPRES_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_vxi, list); break;
1037 case estTC_INT: ret = do_cpte_doubles(xd, cptpEST, i, sflags, state->ngtc, &state->therm_integral, list); break;
1038 case estVETA: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->veta, list); break;
1039 case estVOL0: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->vol0, list); break;
1040 case estX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->x, list); break;
1041 case estV: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->v, list); break;
1042 case estSDX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->sd_X, list); break;
1043 case estLD_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, state->nrng, rng_p, list); break;
1044 case estLD_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, state->nrngi, rngi_p, list); break;
1045 case estMC_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, state->nmcrng, (int **)&state->mc_rng, list); break;
1046 case estMC_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 1, &state->mc_rngi, list); break;
1047 case estDISRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.disre_initf, list); break;
1048 case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.ndisrepairs, &state->hist.disre_rm3tav, list); break;
1049 case estORIRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.orire_initf, list); break;
1050 case estORIRE_DTAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.norire_Dtav, &state->hist.orire_Dtav, list); break;
1052 gmx_fatal(FARGS, "Unknown state entry %d\n"
1053 "You are probably reading a new checkpoint file with old code", i);
1061 static int do_cpt_ekinstate(XDR *xd, int fflags, ekinstate_t *ekins,
1069 for (i = 0; (i < eeksNR && ret == 0); i++)
1071 if (fflags & (1<<i))
1076 case eeksEKIN_N: ret = do_cpte_int(xd, cptpEEKS, i, fflags, &ekins->ekin_n, list); break;
1077 case eeksEKINH: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh, list); break;
1078 case eeksEKINF: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinf, list); break;
1079 case eeksEKINO: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh_old, list); break;
1080 case eeksEKINTOTAL: ret = do_cpte_matrix(xd, cptpEEKS, i, fflags, ekins->ekin_total, list); break;
1081 case eeksEKINSCALEF: ret = do_cpte_doubles(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinscalef_nhc, list); break;
1082 case eeksVSCALE: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->vscale_nhc, list); break;
1083 case eeksEKINSCALEH: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->ekinscaleh_nhc, list); break;
1084 case eeksDEKINDL: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->dekindl, list); break;
1085 case eeksMVCOS: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->mvcos, list); break;
1087 gmx_fatal(FARGS, "Unknown ekin data state entry %d\n"
1088 "You are probably reading a new checkpoint file with old code", i);
1097 static int do_cpt_swapstate(XDR *xd, gmx_bool bRead, swapstate_t *swapstate, FILE *list)
1101 int swap_cpt_version = 1;
1104 if (eswapNO == swapstate->eSwapCoords)
1109 swapstate->bFromCpt = bRead;
1111 do_cpt_int_err(xd, "swap checkpoint version", &swap_cpt_version, list);
1112 do_cpt_int_err(xd, "swap coupling steps", &swapstate->nAverage, list);
1114 /* When reading, init_swapcoords has not been called yet,
1115 * so we have to allocate memory first. */
1117 for (ic = 0; ic < eCompNR; ic++)
1119 for (ii = 0; ii < eIonNR; ii++)
1123 do_cpt_int_err(xd, "swap requested atoms", &swapstate->nat_req[ic][ii], list);
1127 do_cpt_int_err(xd, "swap requested atoms p", swapstate->nat_req_p[ic][ii], list);
1132 do_cpt_int_err(xd, "swap influx netto", &swapstate->inflow_netto[ic][ii], list);
1136 do_cpt_int_err(xd, "swap influx netto p", swapstate->inflow_netto_p[ic][ii], list);
1139 if (bRead && (NULL == swapstate->nat_past[ic][ii]) )
1141 snew(swapstate->nat_past[ic][ii], swapstate->nAverage);
1144 for (j = 0; j < swapstate->nAverage; j++)
1148 do_cpt_int_err(xd, "swap past atom counts", &swapstate->nat_past[ic][ii][j], list);
1152 do_cpt_int_err(xd, "swap past atom counts p", &swapstate->nat_past_p[ic][ii][j], list);
1158 /* Ion flux per channel */
1159 for (ic = 0; ic < eChanNR; ic++)
1161 for (ii = 0; ii < eIonNR; ii++)
1165 do_cpt_int_err(xd, "channel flux", &swapstate->fluxfromAtoB[ic][ii], list);
1169 do_cpt_int_err(xd, "channel flux p", swapstate->fluxfromAtoB_p[ic][ii], list);
1174 /* Ion flux leakage */
1177 snew(swapstate->fluxleak, 1);
1179 do_cpt_int_err(xd, "flux leakage", swapstate->fluxleak, list);
1182 do_cpt_int_err(xd, "number of ions", &swapstate->nions, list);
1186 snew(swapstate->channel_label, swapstate->nions);
1187 snew(swapstate->comp_from, swapstate->nions);
1190 do_cpt_u_chars(xd, "channel history", swapstate->nions, swapstate->channel_label, list);
1191 do_cpt_u_chars(xd, "domain history", swapstate->nions, swapstate->comp_from, list);
1193 /* Save the last known whole positions to checkpoint
1194 * file to be able to also make multimeric channels whole in PBC */
1195 do_cpt_int_err(xd, "Ch0 atoms", &swapstate->nat[eChan0], list);
1196 do_cpt_int_err(xd, "Ch1 atoms", &swapstate->nat[eChan1], list);
1199 snew(swapstate->xc_old_whole[eChan0], swapstate->nat[eChan0]);
1200 snew(swapstate->xc_old_whole[eChan1], swapstate->nat[eChan1]);
1201 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], swapstate->xc_old_whole[eChan0], list);
1202 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], swapstate->xc_old_whole[eChan1], list);
1206 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], *swapstate->xc_old_whole_p[eChan0], list);
1207 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], *swapstate->xc_old_whole_p[eChan1], list);
1214 static int do_cpt_enerhist(XDR *xd, gmx_bool bRead,
1215 int fflags, energyhistory_t *enerhist,
1226 enerhist->nsteps = 0;
1228 enerhist->nsteps_sim = 0;
1229 enerhist->nsum_sim = 0;
1230 enerhist->dht = NULL;
1232 if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
1234 snew(enerhist->dht, 1);
1235 enerhist->dht->ndh = NULL;
1236 enerhist->dht->dh = NULL;
1237 enerhist->dht->start_lambda_set = FALSE;
1241 for (i = 0; (i < eenhNR && ret == 0); i++)
1243 if (fflags & (1<<i))
1247 case eenhENERGY_N: ret = do_cpte_int(xd, cptpEENH, i, fflags, &enerhist->nener, list); break;
1248 case eenhENERGY_AVER: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_ave, list); break;
1249 case eenhENERGY_SUM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum, list); break;
1250 case eenhENERGY_NSUM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum, list); break;
1251 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum_sim, list); break;
1252 case eenhENERGY_NSUM_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum_sim, list); break;
1253 case eenhENERGY_NSTEPS: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps, list); break;
1254 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps_sim, list); break;
1255 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd, eenh_names[i], &(enerhist->dht->nndh), list);
1256 if (bRead) /* now allocate memory for it */
1258 snew(enerhist->dht->dh, enerhist->dht->nndh);
1259 snew(enerhist->dht->ndh, enerhist->dht->nndh);
1260 for (j = 0; j < enerhist->dht->nndh; j++)
1262 enerhist->dht->ndh[j] = 0;
1263 enerhist->dht->dh[j] = NULL;
1267 case eenhENERGY_DELTA_H_LIST:
1268 for (j = 0; j < enerhist->dht->nndh; j++)
1270 ret = do_cpte_n_reals(xd, cptpEENH, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list);
1273 case eenhENERGY_DELTA_H_STARTTIME:
1274 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_time), list); break;
1275 case eenhENERGY_DELTA_H_STARTLAMBDA:
1276 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_lambda), list); break;
1278 gmx_fatal(FARGS, "Unknown energy history entry %d\n"
1279 "You are probably reading a new checkpoint file with old code", i);
1284 if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
1286 /* Assume we have an old file format and copy sum to sum_sim */
1287 srenew(enerhist->ener_sum_sim, enerhist->nener);
1288 for (i = 0; i < enerhist->nener; i++)
1290 enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1292 fflags |= (1<<eenhENERGY_SUM_SIM);
1295 if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1296 !(fflags & (1<<eenhENERGY_NSTEPS)))
1298 /* Assume we have an old file format and copy nsum to nsteps */
1299 enerhist->nsteps = enerhist->nsum;
1300 fflags |= (1<<eenhENERGY_NSTEPS);
1302 if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1303 !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1305 /* Assume we have an old file format and copy nsum to nsteps */
1306 enerhist->nsteps_sim = enerhist->nsum_sim;
1307 fflags |= (1<<eenhENERGY_NSTEPS_SIM);
1313 static int do_cpt_df_hist(XDR *xd, int fflags, df_history_t *dfhist, FILE *list)
1318 nlambda = dfhist->nlambda;
1321 for (i = 0; (i < edfhNR && ret == 0); i++)
1323 if (fflags & (1<<i))
1327 case edfhBEQUIL: ret = do_cpte_int(xd, cptpEDFH, i, fflags, &dfhist->bEquil, list); break;
1328 case edfhNATLAMBDA: ret = do_cpte_ints(xd, cptpEDFH, i, fflags, nlambda, &dfhist->n_at_lam, list); break;
1329 case edfhWLHISTO: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->wl_histo, list); break;
1330 case edfhWLDELTA: ret = do_cpte_real(xd, cptpEDFH, i, fflags, &dfhist->wl_delta, list); break;
1331 case edfhSUMWEIGHTS: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_weights, list); break;
1332 case edfhSUMDG: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_dg, list); break;
1333 case edfhSUMMINVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_minvar, list); break;
1334 case edfhSUMVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_variance, list); break;
1335 case edfhACCUMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p, list); break;
1336 case edfhACCUMM: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m, list); break;
1337 case edfhACCUMP2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p2, list); break;
1338 case edfhACCUMM2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m2, list); break;
1339 case edfhTIJ: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij, list); break;
1340 case edfhTIJEMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij_empirical, list); break;
1343 gmx_fatal(FARGS, "Unknown df history entry %d\n"
1344 "You are probably reading a new checkpoint file with old code", i);
1353 /* This function stores the last whole configuration of the reference and
1354 * average structure in the .cpt file
1356 static int do_cpt_EDstate(XDR *xd, gmx_bool bRead,
1357 edsamstate_t *EDstate, FILE *list)
1364 EDstate->bFromCpt = bRead;
1366 if (EDstate->nED <= 0)
1371 /* When reading, init_edsam has not been called yet,
1372 * so we have to allocate memory first. */
1375 snew(EDstate->nref, EDstate->nED);
1376 snew(EDstate->old_sref, EDstate->nED);
1377 snew(EDstate->nav, EDstate->nED);
1378 snew(EDstate->old_sav, EDstate->nED);
1381 /* Read/write the last whole conformation of SREF and SAV for each ED dataset (usually only one) */
1382 for (i = 0; i < EDstate->nED; i++)
1384 /* Reference structure SREF */
1385 sprintf(buf, "ED%d # of atoms in reference structure", i+1);
1386 do_cpt_int_err(xd, buf, &EDstate->nref[i], list);
1387 sprintf(buf, "ED%d x_ref", i+1);
1390 snew(EDstate->old_sref[i], EDstate->nref[i]);
1391 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref[i], list);
1395 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref_p[i], list);
1398 /* Average structure SAV */
1399 sprintf(buf, "ED%d # of atoms in average structure", i+1);
1400 do_cpt_int_err(xd, buf, &EDstate->nav[i], list);
1401 sprintf(buf, "ED%d x_av", i+1);
1404 snew(EDstate->old_sav[i], EDstate->nav[i]);
1405 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav[i], list);
1409 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav_p[i], list);
1417 static int do_cpt_files(XDR *xd, gmx_bool bRead,
1418 gmx_file_position_t **p_outputfiles, int *nfiles,
1419 FILE *list, int file_version)
1423 gmx_off_t mask = 0xFFFFFFFFL;
1424 int offset_high, offset_low;
1426 gmx_file_position_t *outputfiles;
1428 if (do_cpt_int(xd, "number of output files", nfiles, list) != 0)
1435 snew(*p_outputfiles, *nfiles);
1438 outputfiles = *p_outputfiles;
1440 for (i = 0; i < *nfiles; i++)
1442 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1445 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1446 strncpy(outputfiles[i].filename, buf, CPTSTRLEN-1);
1452 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1456 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1460 outputfiles[i].offset = ( ((gmx_off_t) offset_high) << 32 ) | ( (gmx_off_t) offset_low & mask );
1464 buf = outputfiles[i].filename;
1465 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1467 offset = outputfiles[i].offset;
1475 offset_low = (int) (offset & mask);
1476 offset_high = (int) ((offset >> 32) & mask);
1478 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1482 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1487 if (file_version >= 8)
1489 if (do_cpt_int(xd, "file_checksum_size", &(outputfiles[i].chksum_size),
1494 if (do_cpt_u_chars(xd, "file_checksum", 16, outputfiles[i].chksum, list) != 0)
1501 outputfiles[i].chksum_size = -1;
1508 void write_checkpoint(const char *fn, gmx_bool bNumberAndKeep,
1509 FILE *fplog, t_commrec *cr,
1510 int eIntegrator, int simulation_part,
1511 gmx_bool bExpanded, int elamstats,
1512 gmx_int64_t step, double t, t_state *state)
1522 char *fntemp; /* the temporary checkpoint file name */
1524 char timebuf[STRLEN];
1525 int nppnodes, npmenodes, flag_64bit;
1526 char buf[1024], suffix[5+STEPSTRSIZE], sbuf[STEPSTRSIZE];
1527 gmx_file_position_t *outputfiles;
1530 int flags_eks, flags_enh, flags_dfh, i;
1535 if (DOMAINDECOMP(cr))
1537 nppnodes = cr->dd->nnodes;
1538 npmenodes = cr->npmenodes;
1542 nppnodes = cr->nnodes;
1552 #ifndef GMX_NO_RENAME
1553 /* make the new temporary filename */
1554 snew(fntemp, strlen(fn)+5+STEPSTRSIZE);
1556 fntemp[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1557 sprintf(suffix, "_%s%s", "step", gmx_step_str(step, sbuf));
1558 strcat(fntemp, suffix);
1559 strcat(fntemp, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1561 /* if we can't rename, we just overwrite the cpt file.
1562 * dangerous if interrupted.
1564 snew(fntemp, strlen(fn));
1568 gmx_ctime_r(&now, timebuf, STRLEN);
1572 fprintf(fplog, "Writing checkpoint, step %s at %s\n\n",
1573 gmx_step_str(step, buf), timebuf);
1576 /* Get offsets for open files */
1577 gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1579 fp = gmx_fio_open(fntemp, "w");
1581 if (state->ekinstate.bUpToDate)
1584 ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) |
1585 (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) |
1586 (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1594 if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1596 flags_enh |= (1<<eenhENERGY_N);
1597 if (state->enerhist.nsum > 0)
1599 flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1600 (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSUM));
1602 if (state->enerhist.nsum_sim > 0)
1604 flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSTEPS_SIM) |
1605 (1<<eenhENERGY_NSUM_SIM));
1607 if (state->enerhist.dht)
1609 flags_enh |= ( (1<< eenhENERGY_DELTA_H_NN) |
1610 (1<< eenhENERGY_DELTA_H_LIST) |
1611 (1<< eenhENERGY_DELTA_H_STARTTIME) |
1612 (1<< eenhENERGY_DELTA_H_STARTLAMBDA) );
1618 flags_dfh = ((1<<edfhBEQUIL) | (1<<edfhNATLAMBDA) | (1<<edfhSUMWEIGHTS) | (1<<edfhSUMDG) |
1619 (1<<edfhTIJ) | (1<<edfhTIJEMP));
1622 flags_dfh |= ((1<<edfhWLDELTA) | (1<<edfhWLHISTO));
1624 if ((elamstats == elamstatsMINVAR) || (elamstats == elamstatsBARKER) || (elamstats == elamstatsMETROPOLIS))
1626 flags_dfh |= ((1<<edfhACCUMP) | (1<<edfhACCUMM) | (1<<edfhACCUMP2) | (1<<edfhACCUMM2)
1627 | (1<<edfhSUMMINVAR) | (1<<edfhSUMVAR));
1635 /* We can check many more things now (CPU, acceleration, etc), but
1636 * it is highly unlikely to have two separate builds with exactly
1637 * the same version, user, time, and build host!
1640 version = gmx_strdup(VERSION);
1641 btime = gmx_strdup(BUILD_TIME);
1642 buser = gmx_strdup(BUILD_USER);
1643 bhost = gmx_strdup(BUILD_HOST);
1645 double_prec = GMX_CPT_BUILD_DP;
1646 fprog = gmx_strdup(Program());
1648 ftime = &(timebuf[0]);
1650 do_cpt_header(gmx_fio_getxdr(fp), FALSE, &file_version,
1651 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1652 &eIntegrator, &simulation_part, &step, &t, &nppnodes,
1653 DOMAINDECOMP(cr) ? cr->dd->nc : NULL, &npmenodes,
1654 &state->natoms, &state->ngtc, &state->nnhpres,
1655 &state->nhchainlength, &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
1656 &state->edsamstate.nED, &state->swapstate.eSwapCoords,
1665 if ((do_cpt_state(gmx_fio_getxdr(fp), FALSE, state->flags, state, TRUE, NULL) < 0) ||
1666 (do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL) < 0) ||
1667 (do_cpt_enerhist(gmx_fio_getxdr(fp), FALSE, flags_enh, &state->enerhist, NULL) < 0) ||
1668 (do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL) < 0) ||
1669 (do_cpt_EDstate(gmx_fio_getxdr(fp), FALSE, &state->edsamstate, NULL) < 0) ||
1670 (do_cpt_swapstate(gmx_fio_getxdr(fp), FALSE, &state->swapstate, NULL) < 0) ||
1671 (do_cpt_files(gmx_fio_getxdr(fp), FALSE, &outputfiles, &noutputfiles, NULL,
1674 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1677 do_cpt_footer(gmx_fio_getxdr(fp), file_version);
1679 /* we really, REALLY, want to make sure to physically write the checkpoint,
1680 and all the files it depends on, out to disk. Because we've
1681 opened the checkpoint with gmx_fio_open(), it's in our list
1683 ret = gmx_fio_all_output_fsync();
1689 "Cannot fsync '%s'; maybe you are out of disk space?",
1690 gmx_fio_getname(ret));
1692 if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV) == NULL)
1702 if (gmx_fio_close(fp) != 0)
1704 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1707 /* we don't move the checkpoint if the user specified they didn't want it,
1708 or if the fsyncs failed */
1709 #ifndef GMX_NO_RENAME
1710 if (!bNumberAndKeep && !ret)
1714 /* Rename the previous checkpoint file */
1716 buf[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1717 strcat(buf, "_prev");
1718 strcat(buf, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1720 /* we copy here so that if something goes wrong between now and
1721 * the rename below, there's always a state.cpt.
1722 * If renames are atomic (such as in POSIX systems),
1723 * this copying should be unneccesary.
1725 gmx_file_copy(fn, buf, FALSE);
1726 /* We don't really care if this fails:
1727 * there's already a new checkpoint.
1730 gmx_file_rename(fn, buf);
1733 if (gmx_file_rename(fntemp, fn) != 0)
1735 gmx_file("Cannot rename checkpoint file; maybe you are out of disk space?");
1738 #endif /* GMX_NO_RENAME */
1744 /*code for alternate checkpointing scheme. moved from top of loop over
1746 fcRequestCheckPoint();
1747 if (fcCheckPointParallel( cr->nodeid, NULL, 0) == 0)
1749 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", step );
1751 #endif /* end GMX_FAHCORE block */
1754 static void print_flag_mismatch(FILE *fplog, int sflags, int fflags)
1758 fprintf(fplog, "\nState entry mismatch between the simulation and the checkpoint file\n");
1759 fprintf(fplog, "Entries which are not present in the checkpoint file will not be updated\n");
1760 fprintf(fplog, " %24s %11s %11s\n", "", "simulation", "checkpoint");
1761 for (i = 0; i < estNR; i++)
1763 if ((sflags & (1<<i)) || (fflags & (1<<i)))
1765 fprintf(fplog, " %24s %11s %11s\n",
1767 (sflags & (1<<i)) ? " present " : "not present",
1768 (fflags & (1<<i)) ? " present " : "not present");
1773 static void check_int(FILE *fplog, const char *type, int p, int f, gmx_bool *mm)
1775 FILE *fp = fplog ? fplog : stderr;
1779 fprintf(fp, " %s mismatch,\n", type);
1780 fprintf(fp, " current program: %d\n", p);
1781 fprintf(fp, " checkpoint file: %d\n", f);
1787 static void check_string(FILE *fplog, const char *type, const char *p,
1788 const char *f, gmx_bool *mm)
1790 FILE *fp = fplog ? fplog : stderr;
1792 if (strcmp(p, f) != 0)
1794 fprintf(fp, " %s mismatch,\n", type);
1795 fprintf(fp, " current program: %s\n", p);
1796 fprintf(fp, " checkpoint file: %s\n", f);
1802 static void check_match(FILE *fplog,
1804 char *btime, char *buser, char *bhost, int double_prec,
1806 t_commrec *cr, gmx_bool bPartDecomp, int npp_f, int npme_f,
1807 ivec dd_nc, ivec dd_nc_f)
1814 check_string(fplog, "Version", VERSION, version, &mm);
1815 check_string(fplog, "Build time", BUILD_TIME, btime, &mm);
1816 check_string(fplog, "Build user", BUILD_USER, buser, &mm);
1817 check_string(fplog, "Build host", BUILD_HOST, bhost, &mm);
1818 check_int (fplog, "Double prec.", GMX_CPT_BUILD_DP, double_prec, &mm);
1819 check_string(fplog, "Program name", Program(), fprog, &mm);
1821 check_int (fplog, "#nodes", cr->nnodes, npp_f+npme_f, &mm);
1830 check_int (fplog, "#PME-nodes", cr->npmenodes, npme_f, &mm);
1833 if (cr->npmenodes >= 0)
1835 npp -= cr->npmenodes;
1839 check_int (fplog, "#DD-cells[x]", dd_nc[XX], dd_nc_f[XX], &mm);
1840 check_int (fplog, "#DD-cells[y]", dd_nc[YY], dd_nc_f[YY], &mm);
1841 check_int (fplog, "#DD-cells[z]", dd_nc[ZZ], dd_nc_f[ZZ], &mm);
1848 "Gromacs binary or parallel settings not identical to previous run.\n"
1849 "Continuation is exact, but is not guaranteed to be binary identical%s.\n\n",
1850 fplog ? ",\n see the log file for details" : "");
1855 "Gromacs binary or parallel settings not identical to previous run.\n"
1856 "Continuation is exact, but is not guaranteed to be binary identical.\n\n");
1861 static void read_checkpoint(const char *fn, FILE **pfplog,
1862 t_commrec *cr, gmx_bool bPartDecomp, ivec dd_nc,
1863 int eIntegrator, int *init_fep_state, gmx_int64_t *step, double *t,
1864 t_state *state, gmx_bool *bReadRNG, gmx_bool *bReadEkin,
1865 int *simulation_part,
1866 gmx_bool bAppendOutputFiles, gmx_bool bForceAppend)
1871 char *version, *btime, *buser, *bhost, *fprog, *ftime;
1873 char filename[STRLEN], buf[STEPSTRSIZE];
1874 int nppnodes, eIntegrator_f, nppnodes_f, npmenodes_f;
1876 int natoms, ngtc, nnhpres, nhchainlength, nlambda, fflags, flags_eks, flags_enh, flags_dfh;
1879 gmx_file_position_t *outputfiles;
1881 t_fileio *chksum_file;
1882 FILE * fplog = *pfplog;
1883 unsigned char digest[16];
1884 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1885 struct flock fl; /* don't initialize here: the struct order is OS
1889 const char *int_warn =
1890 "WARNING: The checkpoint file was generated with integrator %s,\n"
1891 " while the simulation uses integrator %s\n\n";
1892 const char *sd_note =
1893 "NOTE: The checkpoint file was for %d nodes doing SD or BD,\n"
1894 " while the simulation uses %d SD or BD nodes,\n"
1895 " continuation will be exact, except for the random state\n\n";
1897 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1898 fl.l_type = F_WRLCK;
1899 fl.l_whence = SEEK_SET;
1908 "read_checkpoint not (yet) supported with particle decomposition");
1911 fp = gmx_fio_open(fn, "r");
1912 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
1913 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1914 &eIntegrator_f, simulation_part, step, t,
1915 &nppnodes_f, dd_nc_f, &npmenodes_f,
1916 &natoms, &ngtc, &nnhpres, &nhchainlength, &nlambda,
1917 &fflags, &flags_eks, &flags_enh, &flags_dfh,
1918 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
1920 if (bAppendOutputFiles &&
1921 file_version >= 13 && double_prec != GMX_CPT_BUILD_DP)
1923 gmx_fatal(FARGS, "Output file appending requested, but the code and checkpoint file precision (single/double) don't match");
1926 if (cr == NULL || MASTER(cr))
1928 fprintf(stderr, "\nReading checkpoint file %s generated: %s\n\n",
1932 /* This will not be written if we do appending, since fplog is still NULL then */
1935 fprintf(fplog, "\n");
1936 fprintf(fplog, "Reading checkpoint file %s\n", fn);
1937 fprintf(fplog, " file generated by: %s\n", fprog);
1938 fprintf(fplog, " file generated at: %s\n", ftime);
1939 fprintf(fplog, " GROMACS build time: %s\n", btime);
1940 fprintf(fplog, " GROMACS build user: %s\n", buser);
1941 fprintf(fplog, " GROMACS build host: %s\n", bhost);
1942 fprintf(fplog, " GROMACS double prec.: %d\n", double_prec);
1943 fprintf(fplog, " simulation part #: %d\n", *simulation_part);
1944 fprintf(fplog, " step: %s\n", gmx_step_str(*step, buf));
1945 fprintf(fplog, " time: %f\n", *t);
1946 fprintf(fplog, "\n");
1949 if (natoms != state->natoms)
1951 gmx_fatal(FARGS, "Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms", natoms, state->natoms);
1953 if (ngtc != state->ngtc)
1955 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);
1957 if (nnhpres != state->nnhpres)
1959 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);
1962 if (nlambda != state->dfhist.nlambda)
1964 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);
1967 init_gtc_state(state, state->ngtc, state->nnhpres, nhchainlength); /* need to keep this here to keep the tpr format working */
1968 /* write over whatever was read; we use the number of Nose-Hoover chains from the checkpoint */
1970 if (eIntegrator_f != eIntegrator)
1974 fprintf(stderr, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1976 if (bAppendOutputFiles)
1979 "Output file appending requested, but input/checkpoint integrators do not match.\n"
1980 "Stopping the run to prevent you from ruining all your data...\n"
1981 "If you _really_ know what you are doing, try with the -noappend option.\n");
1985 fprintf(fplog, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1994 else if (bPartDecomp)
1996 nppnodes = cr->nnodes;
1999 else if (cr->nnodes == nppnodes_f + npmenodes_f)
2001 if (cr->npmenodes < 0)
2003 cr->npmenodes = npmenodes_f;
2005 nppnodes = cr->nnodes - cr->npmenodes;
2006 if (nppnodes == nppnodes_f)
2008 for (d = 0; d < DIM; d++)
2012 dd_nc[d] = dd_nc_f[d];
2019 /* The number of PP nodes has not been set yet */
2023 if ((EI_SD(eIntegrator) || eIntegrator == eiBD) && nppnodes > 0)
2025 /* Correct the RNG state size for the number of PP nodes.
2026 * Such assignments should all be moved to one central function.
2028 state->nrng = nppnodes*gmx_rng_n();
2029 state->nrngi = nppnodes;
2033 if (fflags != state->flags)
2038 if (bAppendOutputFiles)
2041 "Output file appending requested, but input and checkpoint states are not identical.\n"
2042 "Stopping the run to prevent you from ruining all your data...\n"
2043 "You can try with the -noappend option, and get more info in the log file.\n");
2046 if (getenv("GMX_ALLOW_CPT_MISMATCH") == NULL)
2048 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");
2053 "WARNING: The checkpoint state entries do not match the simulation,\n"
2054 " see the log file for details\n\n");
2060 print_flag_mismatch(fplog, state->flags, fflags);
2065 if ((EI_SD(eIntegrator) || eIntegrator == eiBD) &&
2066 nppnodes != nppnodes_f)
2071 fprintf(stderr, sd_note, nppnodes_f, nppnodes);
2075 fprintf(fplog, sd_note, nppnodes_f, nppnodes);
2080 check_match(fplog, version, btime, buser, bhost, double_prec, fprog,
2081 cr, bPartDecomp, nppnodes_f, npmenodes_f, dd_nc, dd_nc_f);
2084 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, fflags, state, *bReadRNG, NULL);
2085 *init_fep_state = state->fep_state; /* there should be a better way to do this than setting it here.
2086 Investigate for 5.0. */
2091 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2096 *bReadEkin = ((flags_eks & (1<<eeksEKINH)) || (flags_eks & (1<<eeksEKINF)) || (flags_eks & (1<<eeksEKINO)) ||
2097 ((flags_eks & (1<<eeksEKINSCALEF)) | (flags_eks & (1<<eeksEKINSCALEH)) | (flags_eks & (1<<eeksVSCALE))));
2099 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2100 flags_enh, &state->enerhist, NULL);
2106 if (file_version < 6)
2108 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.";
2110 fprintf(stderr, "\nWARNING: %s\n\n", warn);
2113 fprintf(fplog, "\nWARNING: %s\n\n", warn);
2115 state->enerhist.nsum = *step;
2116 state->enerhist.nsum_sim = *step;
2119 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2125 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2131 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2137 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, NULL, file_version);
2143 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2148 if (gmx_fio_close(fp) != 0)
2150 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2159 /* If the user wants to append to output files,
2160 * we use the file pointer positions of the output files stored
2161 * in the checkpoint file and truncate the files such that any frames
2162 * written after the checkpoint time are removed.
2163 * All files are md5sum checked such that we can be sure that
2164 * we do not truncate other (maybe imprortant) files.
2166 if (bAppendOutputFiles)
2168 if (fn2ftp(outputfiles[0].filename) != efLOG)
2170 /* make sure first file is log file so that it is OK to use it for
2173 gmx_fatal(FARGS, "The first output file should always be the log "
2174 "file but instead is: %s. Cannot do appending because of this condition.", outputfiles[0].filename);
2176 for (i = 0; i < nfiles; i++)
2178 if (outputfiles[i].offset < 0)
2180 gmx_fatal(FARGS, "The original run wrote a file called '%s' which "
2181 "is larger than 2 GB, but mdrun did not support large file"
2182 " offsets. Can not append. Run mdrun with -noappend",
2183 outputfiles[i].filename);
2186 chksum_file = gmx_fio_open(outputfiles[i].filename, "a");
2189 chksum_file = gmx_fio_open(outputfiles[i].filename, "r+");
2194 /* Note that there are systems where the lock operation
2195 * will succeed, but a second process can also lock the file.
2196 * We should probably try to detect this.
2198 #if defined __native_client__
2202 #elif defined GMX_NATIVE_WINDOWS
2203 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX) == -1)
2205 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl) == -1)
2208 if (errno == ENOSYS)
2212 gmx_fatal(FARGS, "File locking is not supported on this system. Use -noappend or specify -append explicitly to append anyhow.");
2216 fprintf(stderr, "\nNOTE: File locking is not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2219 fprintf(fplog, "\nNOTE: File locking not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2223 else if (errno == EACCES || errno == EAGAIN)
2225 gmx_fatal(FARGS, "Failed to lock: %s. Already running "
2226 "simulation?", outputfiles[i].filename);
2230 gmx_fatal(FARGS, "Failed to lock: %s. %s.",
2231 outputfiles[i].filename, strerror(errno));
2236 /* compute md5 chksum */
2237 if (outputfiles[i].chksum_size != -1)
2239 if (gmx_fio_get_file_md5(chksum_file, outputfiles[i].offset,
2240 digest) != outputfiles[i].chksum_size) /*at the end of the call the file position is at the end of the file*/
2242 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.",
2243 outputfiles[i].chksum_size,
2244 outputfiles[i].filename);
2247 if (i == 0) /*log file needs to be seeked in case we need to truncate (other files are truncated below)*/
2249 if (gmx_fio_seek(chksum_file, outputfiles[i].offset))
2251 gmx_fatal(FARGS, "Seek error! Failed to truncate log-file: %s.", strerror(errno));
2256 if (i == 0) /*open log file here - so that lock is never lifted
2257 after chksum is calculated */
2259 *pfplog = gmx_fio_getfp(chksum_file);
2263 gmx_fio_close(chksum_file);
2266 /* compare md5 chksum */
2267 if (outputfiles[i].chksum_size != -1 &&
2268 memcmp(digest, outputfiles[i].chksum, 16) != 0)
2272 fprintf(debug, "chksum for %s: ", outputfiles[i].filename);
2273 for (j = 0; j < 16; j++)
2275 fprintf(debug, "%02x", digest[j]);
2277 fprintf(debug, "\n");
2279 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.",
2280 outputfiles[i].filename);
2285 if (i != 0) /*log file is already seeked to correct position */
2287 #ifdef GMX_NATIVE_WINDOWS
2288 rc = gmx_wintruncate(outputfiles[i].filename, outputfiles[i].offset);
2290 rc = truncate(outputfiles[i].filename, outputfiles[i].offset);
2294 gmx_fatal(FARGS, "Truncation of file %s failed. Cannot do appending because of this failure.", outputfiles[i].filename);
2304 void load_checkpoint(const char *fn, FILE **fplog,
2305 t_commrec *cr, gmx_bool bPartDecomp, ivec dd_nc,
2306 t_inputrec *ir, t_state *state,
2307 gmx_bool *bReadRNG, gmx_bool *bReadEkin,
2308 gmx_bool bAppend, gmx_bool bForceAppend)
2315 /* Read the state from the checkpoint file */
2316 read_checkpoint(fn, fplog,
2317 cr, bPartDecomp, dd_nc,
2318 ir->eI, &(ir->fepvals->init_fep_state), &step, &t, state, bReadRNG, bReadEkin,
2319 &ir->simulation_part, bAppend, bForceAppend);
2323 gmx_bcast(sizeof(cr->npmenodes), &cr->npmenodes, cr);
2324 gmx_bcast(DIM*sizeof(dd_nc[0]), dd_nc, cr);
2325 gmx_bcast(sizeof(step), &step, cr);
2326 gmx_bcast(sizeof(*bReadRNG), bReadRNG, cr);
2327 gmx_bcast(sizeof(*bReadEkin), bReadEkin, cr);
2329 ir->bContinuation = TRUE;
2330 if (ir->nsteps >= 0)
2332 ir->nsteps += ir->init_step - step;
2334 ir->init_step = step;
2335 ir->simulation_part += 1;
2338 static void read_checkpoint_data(t_fileio *fp, int *simulation_part,
2339 gmx_int64_t *step, double *t, t_state *state,
2341 int *nfiles, gmx_file_position_t **outputfiles)
2344 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2349 int flags_eks, flags_enh, flags_dfh;
2351 gmx_file_position_t *files_loc = NULL;
2354 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2355 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2356 &eIntegrator, simulation_part, step, t, &nppnodes, dd_nc, &npme,
2357 &state->natoms, &state->ngtc, &state->nnhpres, &state->nhchainlength,
2358 &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
2359 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
2361 do_cpt_state(gmx_fio_getxdr(fp), TRUE, state->flags, state, bReadRNG, NULL);
2366 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2371 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2372 flags_enh, &state->enerhist, NULL);
2377 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2383 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2389 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2395 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE,
2396 outputfiles != NULL ? outputfiles : &files_loc,
2397 outputfiles != NULL ? nfiles : &nfiles_loc,
2398 NULL, file_version);
2399 if (files_loc != NULL)
2409 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2423 read_checkpoint_state(const char *fn, int *simulation_part,
2424 gmx_int64_t *step, double *t, t_state *state)
2428 fp = gmx_fio_open(fn, "r");
2429 read_checkpoint_data(fp, simulation_part, step, t, state, FALSE, NULL, NULL);
2430 if (gmx_fio_close(fp) != 0)
2432 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2436 void read_checkpoint_trxframe(t_fileio *fp, t_trxframe *fr)
2438 /* This next line is nasty because the sub-structures of t_state
2439 * cannot be assumed to be zeroed (or even initialized in ways the
2440 * rest of the code might assume). Using snew would be better, but
2441 * this will all go away for 5.0. */
2443 int simulation_part;
2447 init_state(&state, 0, 0, 0, 0, 0);
2449 read_checkpoint_data(fp, &simulation_part, &step, &t, &state, FALSE, NULL, NULL);
2451 fr->natoms = state.natoms;
2454 fr->step = gmx_int64_to_int(step,
2455 "conversion of checkpoint to trajectory");
2459 fr->lambda = state.lambda[efptFEP];
2460 fr->fep_state = state.fep_state;
2462 fr->bX = (state.flags & (1<<estX));
2468 fr->bV = (state.flags & (1<<estV));
2475 fr->bBox = (state.flags & (1<<estBOX));
2478 copy_mat(state.box, fr->box);
2483 void list_checkpoint(const char *fn, FILE *out)
2487 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2489 int eIntegrator, simulation_part, nppnodes, npme;
2494 int flags_eks, flags_enh, flags_dfh;
2498 gmx_file_position_t *outputfiles;
2501 init_state(&state, -1, -1, -1, -1, 0);
2503 fp = gmx_fio_open(fn, "r");
2504 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2505 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2506 &eIntegrator, &simulation_part, &step, &t, &nppnodes, dd_nc, &npme,
2507 &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
2508 &(state.dfhist.nlambda), &state.flags,
2509 &flags_eks, &flags_enh, &flags_dfh, &state.edsamstate.nED,
2510 &state.swapstate.eSwapCoords, out);
2511 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, state.flags, &state, TRUE, out);
2516 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state.ekinstate, out);
2521 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2522 flags_enh, &state.enerhist, out);
2526 ret = do_cpt_df_hist(gmx_fio_getxdr(fp),
2527 flags_dfh, &state.dfhist, out);
2532 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state.edsamstate, out);
2537 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state.swapstate, out);
2542 do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, out, file_version);
2547 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2554 if (gmx_fio_close(fp) != 0)
2556 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2563 static gmx_bool exist_output_file(const char *fnm_cp, int nfile, const t_filenm fnm[])
2567 /* Check if the output file name stored in the checkpoint file
2568 * is one of the output file names of mdrun.
2572 !(is_output(&fnm[i]) && strcmp(fnm_cp, fnm[i].fns[0]) == 0))
2577 return (i < nfile && gmx_fexist(fnm_cp));
2580 /* This routine cannot print tons of data, since it is called before the log file is opened. */
2581 gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
2582 gmx_int64_t *cpt_step, t_commrec *cr,
2583 gmx_bool bAppendReq,
2584 int nfile, const t_filenm fnm[],
2585 const char *part_suffix, gmx_bool *bAddPart)
2588 gmx_int64_t step = 0;
2590 /* This next line is nasty because the sub-structures of t_state
2591 * cannot be assumed to be zeroed (or even initialized in ways the
2592 * rest of the code might assume). Using snew would be better, but
2593 * this will all go away for 5.0. */
2596 gmx_file_position_t *outputfiles;
2599 char *fn, suf_up[STRLEN];
2605 if (!gmx_fexist(filename) || (!(fp = gmx_fio_open(filename, "r")) ))
2607 *simulation_part = 0;
2611 init_state(&state, 0, 0, 0, 0, 0);
2613 read_checkpoint_data(fp, simulation_part, &step, &t, &state, FALSE,
2614 &nfiles, &outputfiles);
2615 if (gmx_fio_close(fp) != 0)
2617 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2624 for (f = 0; f < nfiles; f++)
2626 if (exist_output_file(outputfiles[f].filename, nfile, fnm))
2631 if (nexist == nfiles)
2633 bAppend = bAppendReq;
2635 else if (nexist > 0)
2638 "Output file appending has been requested,\n"
2639 "but some output files listed in the checkpoint file %s\n"
2640 "are not present or are named differently by the current program:\n",
2642 fprintf(stderr, "output files present:");
2643 for (f = 0; f < nfiles; f++)
2645 if (exist_output_file(outputfiles[f].filename,
2648 fprintf(stderr, " %s", outputfiles[f].filename);
2651 fprintf(stderr, "\n");
2652 fprintf(stderr, "output files not present or named differently:");
2653 for (f = 0; f < nfiles; f++)
2655 if (!exist_output_file(outputfiles[f].filename,
2658 fprintf(stderr, " %s", outputfiles[f].filename);
2661 fprintf(stderr, "\n");
2663 gmx_fatal(FARGS, "File appending requested, but only %d of the %d output files are present", nexist, nfiles);
2671 gmx_fatal(FARGS, "File appending requested, but no output file information is stored in the checkpoint file");
2673 fn = outputfiles[0].filename;
2674 if (strlen(fn) < 4 ||
2675 gmx_strcasecmp(fn+strlen(fn)-4, ftp2ext(efLOG)) == 0)
2677 gmx_fatal(FARGS, "File appending requested, but the log file is not the first file listed in the checkpoint file");
2679 /* Set bAddPart to whether the suffix string '.part' is present
2680 * in the log file name.
2682 strcpy(suf_up, part_suffix);
2684 *bAddPart = (strstr(fn, part_suffix) != NULL ||
2685 strstr(fn, suf_up) != NULL);
2693 gmx_bcast(sizeof(*simulation_part), simulation_part, cr);
2695 if (*simulation_part > 0 && bAppendReq)
2697 gmx_bcast(sizeof(bAppend), &bAppend, cr);
2698 gmx_bcast(sizeof(*bAddPart), bAddPart, cr);
2701 if (NULL != cpt_step)