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>
64 #include "types/commrec.h"
65 #include "gromacs/utility/smalloc.h"
69 #include "checkpoint.h"
71 #include "gromacs/utility/cstringutil.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"
79 #include "gromacs/utility/baseversion.h"
80 #include "gmx_fatal.h"
82 #include "buildinfo.h"
88 #define CPT_MAGIC1 171817
89 #define CPT_MAGIC2 171819
90 #define CPTSTRLEN 1024
93 #define GMX_CPT_BUILD_DP 1
95 #define GMX_CPT_BUILD_DP 0
98 /* cpt_version should normally only be changed
99 * when the header of footer format changes.
100 * The state data format itself is backward and forward compatible.
101 * But old code can not read a new entry that is present in the file
102 * (but can read a new format when new entries are not present).
104 static const int cpt_version = 16;
107 const char *est_names[estNR] =
110 "box", "box-rel", "box-v", "pres_prev",
111 "nosehoover-xi", "thermostat-integral",
112 "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
113 "disre_initf", "disre_rm3tav",
114 "orire_initf", "orire_Dtav",
115 "svir_prev", "nosehoover-vxi", "v_eta", "vol0", "nhpres_xi", "nhpres_vxi", "fvir_prev", "fep_state", "MC-rng", "MC-rng-i"
119 eeksEKIN_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksEKINF, eeksEKINO, eeksEKINSCALEF, eeksEKINSCALEH, eeksVSCALE, eeksEKINTOTAL, eeksNR
122 const char *eeks_names[eeksNR] =
124 "Ekin_n", "Ekinh", "dEkindlambda", "mv_cos",
125 "Ekinf", "Ekinh_old", "EkinScaleF_NHC", "EkinScaleH_NHC", "Vscale_NHC", "Ekin_Total"
129 eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
130 eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
131 eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM,
132 eenhENERGY_DELTA_H_NN,
133 eenhENERGY_DELTA_H_LIST,
134 eenhENERGY_DELTA_H_STARTTIME,
135 eenhENERGY_DELTA_H_STARTLAMBDA,
139 const char *eenh_names[eenhNR] =
141 "energy_n", "energy_aver", "energy_sum", "energy_nsum",
142 "energy_sum_sim", "energy_nsum_sim",
143 "energy_nsteps", "energy_nsteps_sim",
145 "energy_delta_h_list",
146 "energy_delta_h_start_time",
147 "energy_delta_h_start_lambda"
150 /* free energy history variables -- need to be preserved over checkpoint */
152 edfhBEQUIL, edfhNATLAMBDA, edfhWLHISTO, edfhWLDELTA, edfhSUMWEIGHTS, edfhSUMDG, edfhSUMMINVAR, edfhSUMVAR,
153 edfhACCUMP, edfhACCUMM, edfhACCUMP2, edfhACCUMM2, edfhTIJ, edfhTIJEMP, edfhNR
155 /* free energy history variable names */
156 const char *edfh_names[edfhNR] =
158 "bEquilibrated", "N_at_state", "Wang-Landau Histogram", "Wang-Landau Delta", "Weights", "Free Energies", "minvar", "variance",
159 "accumulated_plus", "accumulated_minus", "accumulated_plus_2", "accumulated_minus_2", "Tij", "Tij_empirical"
162 #ifdef GMX_NATIVE_WINDOWS
164 gmx_wintruncate(const char *filename, __int64 size)
167 /*we do this elsewhere*/
173 fp = fopen(filename, "rb+");
181 return _chsize_s( fileno(fp), size);
183 return _chsize( fileno(fp), size);
191 ecprREAL, ecprRVEC, ecprMATRIX
195 cptpEST, cptpEEKS, cptpEENH, cptpEDFH
197 /* enums for the different components of checkpoint variables, replacing the hard coded ones.
198 cptpEST - state variables.
199 cptpEEKS - Kinetic energy state variables.
200 cptpEENH - Energy history state variables.
201 cptpEDFH - free energy history variables.
205 static const char *st_names(int cptp, int ecpt)
209 case cptpEST: return est_names [ecpt]; break;
210 case cptpEEKS: return eeks_names[ecpt]; break;
211 case cptpEENH: return eenh_names[ecpt]; break;
212 case cptpEDFH: return edfh_names[ecpt]; break;
218 static void cp_warning(FILE *fp)
220 fprintf(fp, "\nWARNING: Checkpoint file is corrupted or truncated\n\n");
223 static void cp_error()
225 gmx_fatal(FARGS, "Checkpoint file corrupted/truncated, or maybe you are out of disk space?");
228 static void do_cpt_string_err(XDR *xd, gmx_bool bRead, const char *desc, char **s, FILE *list)
236 res = xdr_string(xd, s, CPTSTRLEN);
243 fprintf(list, "%s = %s\n", desc, *s);
248 static int do_cpt_int(XDR *xd, const char *desc, int *i, FILE *list)
252 res = xdr_int(xd, i);
259 fprintf(list, "%s = %d\n", desc, *i);
264 static int do_cpt_u_chars(XDR *xd, const char *desc, int n, unsigned char *i, FILE *list)
270 fprintf(list, "%s = ", desc);
272 for (j = 0; j < n && res; j++)
274 res &= xdr_u_char(xd, &i[j]);
277 fprintf(list, "%02x", i[j]);
292 static void do_cpt_int_err(XDR *xd, const char *desc, int *i, FILE *list)
294 if (do_cpt_int(xd, desc, i, list) < 0)
300 static void do_cpt_step_err(XDR *xd, const char *desc, gmx_int64_t *i, FILE *list)
303 char buf[STEPSTRSIZE];
305 res = xdr_int64(xd, i);
312 fprintf(list, "%s = %s\n", desc, gmx_step_str(*i, buf));
316 static void do_cpt_double_err(XDR *xd, const char *desc, double *f, FILE *list)
320 res = xdr_double(xd, f);
327 fprintf(list, "%s = %f\n", desc, *f);
331 static void do_cpt_real_err(XDR *xd, real *f)
336 res = xdr_double(xd, f);
338 res = xdr_float(xd, f);
346 static void do_cpt_n_rvecs_err(XDR *xd, const char *desc, int n, rvec f[], FILE *list)
350 for (i = 0; i < n; i++)
352 for (j = 0; j < DIM; j++)
354 do_cpt_real_err(xd, &f[i][j]);
360 pr_rvecs(list, 0, desc, f, n);
364 /* If nval >= 0, nval is used; on read this should match the passed value.
365 * If nval n<0, *nptr is used; on read the value is stored in nptr
367 static int do_cpte_reals_low(XDR *xd, int cptp, int ecpt, int sflags,
368 int nval, int *nptr, real **v,
369 FILE *list, int erealtype)
373 int dtc = xdr_datatype_float;
375 int dtc = xdr_datatype_double;
377 real *vp, *va = NULL;
392 gmx_incons("*ntpr=NULL in do_cpte_reals_low");
397 res = xdr_int(xd, &nf);
408 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), nval, nf);
417 res = xdr_int(xd, &dt);
424 fprintf(stderr, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
425 st_names(cptp, ecpt), xdr_datatype_names[dtc],
426 xdr_datatype_names[dt]);
428 if (list || !(sflags & (1<<ecpt)))
441 if (dt == xdr_datatype_float)
443 if (dtc == xdr_datatype_float)
451 res = xdr_vector(xd, (char *)vf, nf,
452 (unsigned int)sizeof(float), (xdrproc_t)xdr_float);
457 if (dtc != xdr_datatype_float)
459 for (i = 0; i < nf; i++)
468 if (dtc == xdr_datatype_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++)
723 reti = do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, &(v[i]), NULL, ecprREAL);
724 if (list && reti == 0)
726 sprintf(name, "%s[%d]", st_names(cptp, ecpt), i);
727 pr_reals(list, 0, name, v[i], n);
737 static int do_cpte_matrices(XDR *xd, int cptp, int ecpt, int sflags,
738 int n, matrix **v, FILE *list)
741 matrix *vp, *va = NULL;
747 res = xdr_int(xd, &nf);
752 if (list == NULL && nf != n)
754 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
756 if (list || !(sflags & (1<<ecpt)))
769 snew(vr, nf*DIM*DIM);
770 for (i = 0; i < nf; i++)
772 for (j = 0; j < DIM; j++)
774 for (k = 0; k < DIM; k++)
776 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
780 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
781 nf*DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
782 for (i = 0; i < nf; i++)
784 for (j = 0; j < DIM; j++)
786 for (k = 0; k < DIM; k++)
788 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
794 if (list && ret == 0)
796 for (i = 0; i < nf; i++)
798 pr_rvecs(list, 0, st_names(cptp, ecpt), vp[i], DIM);
809 static void do_cpt_header(XDR *xd, gmx_bool bRead, int *file_version,
810 char **version, char **btime, char **buser, char **bhost,
812 char **fprog, char **ftime,
813 int *eIntegrator, int *simulation_part,
814 gmx_int64_t *step, double *t,
815 int *nnodes, int *dd_nc, int *npme,
816 int *natoms, int *ngtc, int *nnhpres, int *nhchainlength,
817 int *nlambda, int *flags_state,
818 int *flags_eks, int *flags_enh, int *flags_dfh,
819 int *nED, int *eSwapCoords,
836 res = xdr_int(xd, &magic);
839 gmx_fatal(FARGS, "The checkpoint file is empty/corrupted, or maybe you are out of disk space?");
841 if (magic != CPT_MAGIC1)
843 gmx_fatal(FARGS, "Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
844 "The checkpoint file is corrupted or not a checkpoint file",
850 gmx_gethostname(fhost, 255);
852 do_cpt_string_err(xd, bRead, "GROMACS version", version, list);
853 do_cpt_string_err(xd, bRead, "GROMACS build time", btime, list);
854 do_cpt_string_err(xd, bRead, "GROMACS build user", buser, list);
855 do_cpt_string_err(xd, bRead, "GROMACS build host", bhost, list);
856 do_cpt_string_err(xd, bRead, "generating program", fprog, list);
857 do_cpt_string_err(xd, bRead, "generation time", ftime, list);
858 *file_version = cpt_version;
859 do_cpt_int_err(xd, "checkpoint file version", file_version, list);
860 if (*file_version > cpt_version)
862 gmx_fatal(FARGS, "Attempting to read a checkpoint file of version %d with code of version %d\n", *file_version, cpt_version);
864 if (*file_version >= 13)
866 do_cpt_int_err(xd, "GROMACS double precision", double_prec, list);
872 if (*file_version >= 12)
874 do_cpt_string_err(xd, bRead, "generating host", &fhost, list);
880 do_cpt_int_err(xd, "#atoms", natoms, list);
881 do_cpt_int_err(xd, "#T-coupling groups", ngtc, list);
882 if (*file_version >= 10)
884 do_cpt_int_err(xd, "#Nose-Hoover T-chains", nhchainlength, list);
890 if (*file_version >= 11)
892 do_cpt_int_err(xd, "#Nose-Hoover T-chains for barostat ", nnhpres, list);
898 if (*file_version >= 14)
900 do_cpt_int_err(xd, "# of total lambda states ", nlambda, list);
906 do_cpt_int_err(xd, "integrator", eIntegrator, list);
907 if (*file_version >= 3)
909 do_cpt_int_err(xd, "simulation part #", simulation_part, list);
913 *simulation_part = 1;
915 if (*file_version >= 5)
917 do_cpt_step_err(xd, "step", step, list);
921 do_cpt_int_err(xd, "step", &idum, list);
924 do_cpt_double_err(xd, "t", t, list);
925 do_cpt_int_err(xd, "#PP-ranks", nnodes, list);
927 do_cpt_int_err(xd, "dd_nc[x]", dd_nc ? &(dd_nc[0]) : &idum, list);
928 do_cpt_int_err(xd, "dd_nc[y]", dd_nc ? &(dd_nc[1]) : &idum, list);
929 do_cpt_int_err(xd, "dd_nc[z]", dd_nc ? &(dd_nc[2]) : &idum, list);
930 do_cpt_int_err(xd, "#PME-only ranks", npme, list);
931 do_cpt_int_err(xd, "state flags", flags_state, list);
932 if (*file_version >= 4)
934 do_cpt_int_err(xd, "ekin data flags", flags_eks, list);
935 do_cpt_int_err(xd, "energy history flags", flags_enh, list);
940 *flags_enh = (*flags_state >> (estORIRE_DTAV+1));
941 *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
942 (1<<(estORIRE_DTAV+2)) |
943 (1<<(estORIRE_DTAV+3))));
945 if (*file_version >= 14)
947 do_cpt_int_err(xd, "df history flags", flags_dfh, list);
954 if (*file_version >= 15)
956 do_cpt_int_err(xd, "ED data sets", nED, list);
962 if (*file_version >= 16)
964 do_cpt_int_err(xd, "swap", eSwapCoords, list);
968 static int do_cpt_footer(XDR *xd, int file_version)
973 if (file_version >= 2)
976 res = xdr_int(xd, &magic);
981 if (magic != CPT_MAGIC2)
990 static int do_cpt_state(XDR *xd, gmx_bool bRead,
991 int fflags, t_state *state,
1001 nnht = state->nhchainlength*state->ngtc;
1002 nnhtp = state->nhchainlength*state->nnhpres;
1004 if (bRead) /* we need to allocate space for dfhist if we are reading */
1006 init_df_history(&state->dfhist, state->dfhist.nlambda);
1009 sflags = state->flags;
1010 for (i = 0; (i < estNR && ret == 0); i++)
1012 if (fflags & (1<<i))
1016 case estLAMBDA: ret = do_cpte_reals(xd, cptpEST, i, sflags, efptNR, &(state->lambda), list); break;
1017 case estFEPSTATE: ret = do_cpte_int (xd, cptpEST, i, sflags, &state->fep_state, list); break;
1018 case estBOX: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box, list); break;
1019 case estBOX_REL: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box_rel, list); break;
1020 case estBOXV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->boxv, list); break;
1021 case estPRES_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->pres_prev, list); break;
1022 case estSVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->svir_prev, list); break;
1023 case estFVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->fvir_prev, list); break;
1024 case estNH_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_xi, list); break;
1025 case estNH_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_vxi, list); break;
1026 case estNHPRES_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_xi, list); break;
1027 case estNHPRES_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_vxi, list); break;
1028 case estTC_INT: ret = do_cpte_doubles(xd, cptpEST, i, sflags, state->ngtc, &state->therm_integral, list); break;
1029 case estVETA: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->veta, list); break;
1030 case estVOL0: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->vol0, list); break;
1031 case estX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->x, list); break;
1032 case estV: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->v, list); break;
1033 case estSDX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->sd_X, list); break;
1034 /* The RNG entries are no longer written,
1035 * the next 4 lines are only for reading old files.
1037 case estLD_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1038 case estLD_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1039 case estMC_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1040 case estMC_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1041 case estDISRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.disre_initf, list); break;
1042 case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.ndisrepairs, &state->hist.disre_rm3tav, list); break;
1043 case estORIRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.orire_initf, list); break;
1044 case estORIRE_DTAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.norire_Dtav, &state->hist.orire_Dtav, list); break;
1046 gmx_fatal(FARGS, "Unknown state entry %d\n"
1047 "You are probably reading a new checkpoint file with old code", i);
1055 static int do_cpt_ekinstate(XDR *xd, int fflags, ekinstate_t *ekins,
1063 for (i = 0; (i < eeksNR && ret == 0); i++)
1065 if (fflags & (1<<i))
1070 case eeksEKIN_N: ret = do_cpte_int(xd, cptpEEKS, i, fflags, &ekins->ekin_n, list); break;
1071 case eeksEKINH: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh, list); break;
1072 case eeksEKINF: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinf, list); break;
1073 case eeksEKINO: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh_old, list); break;
1074 case eeksEKINTOTAL: ret = do_cpte_matrix(xd, cptpEEKS, i, fflags, ekins->ekin_total, list); break;
1075 case eeksEKINSCALEF: ret = do_cpte_doubles(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinscalef_nhc, list); break;
1076 case eeksVSCALE: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->vscale_nhc, list); break;
1077 case eeksEKINSCALEH: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->ekinscaleh_nhc, list); break;
1078 case eeksDEKINDL: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->dekindl, list); break;
1079 case eeksMVCOS: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->mvcos, list); break;
1081 gmx_fatal(FARGS, "Unknown ekin data state entry %d\n"
1082 "You are probably reading a new checkpoint file with old code", i);
1091 static int do_cpt_swapstate(XDR *xd, gmx_bool bRead, swapstate_t *swapstate, FILE *list)
1095 int swap_cpt_version = 1;
1098 if (eswapNO == swapstate->eSwapCoords)
1103 swapstate->bFromCpt = bRead;
1105 do_cpt_int_err(xd, "swap checkpoint version", &swap_cpt_version, list);
1106 do_cpt_int_err(xd, "swap coupling steps", &swapstate->nAverage, list);
1108 /* When reading, init_swapcoords has not been called yet,
1109 * so we have to allocate memory first. */
1111 for (ic = 0; ic < eCompNR; ic++)
1113 for (ii = 0; ii < eIonNR; ii++)
1117 do_cpt_int_err(xd, "swap requested atoms", &swapstate->nat_req[ic][ii], list);
1121 do_cpt_int_err(xd, "swap requested atoms p", swapstate->nat_req_p[ic][ii], list);
1126 do_cpt_int_err(xd, "swap influx netto", &swapstate->inflow_netto[ic][ii], list);
1130 do_cpt_int_err(xd, "swap influx netto p", swapstate->inflow_netto_p[ic][ii], list);
1133 if (bRead && (NULL == swapstate->nat_past[ic][ii]) )
1135 snew(swapstate->nat_past[ic][ii], swapstate->nAverage);
1138 for (j = 0; j < swapstate->nAverage; j++)
1142 do_cpt_int_err(xd, "swap past atom counts", &swapstate->nat_past[ic][ii][j], list);
1146 do_cpt_int_err(xd, "swap past atom counts p", &swapstate->nat_past_p[ic][ii][j], list);
1152 /* Ion flux per channel */
1153 for (ic = 0; ic < eChanNR; ic++)
1155 for (ii = 0; ii < eIonNR; ii++)
1159 do_cpt_int_err(xd, "channel flux", &swapstate->fluxfromAtoB[ic][ii], list);
1163 do_cpt_int_err(xd, "channel flux p", swapstate->fluxfromAtoB_p[ic][ii], list);
1168 /* Ion flux leakage */
1171 snew(swapstate->fluxleak, 1);
1173 do_cpt_int_err(xd, "flux leakage", swapstate->fluxleak, list);
1176 do_cpt_int_err(xd, "number of ions", &swapstate->nions, list);
1180 snew(swapstate->channel_label, swapstate->nions);
1181 snew(swapstate->comp_from, swapstate->nions);
1184 do_cpt_u_chars(xd, "channel history", swapstate->nions, swapstate->channel_label, list);
1185 do_cpt_u_chars(xd, "domain history", swapstate->nions, swapstate->comp_from, list);
1187 /* Save the last known whole positions to checkpoint
1188 * file to be able to also make multimeric channels whole in PBC */
1189 do_cpt_int_err(xd, "Ch0 atoms", &swapstate->nat[eChan0], list);
1190 do_cpt_int_err(xd, "Ch1 atoms", &swapstate->nat[eChan1], list);
1193 snew(swapstate->xc_old_whole[eChan0], swapstate->nat[eChan0]);
1194 snew(swapstate->xc_old_whole[eChan1], swapstate->nat[eChan1]);
1195 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], swapstate->xc_old_whole[eChan0], list);
1196 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], swapstate->xc_old_whole[eChan1], list);
1200 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], *swapstate->xc_old_whole_p[eChan0], list);
1201 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], *swapstate->xc_old_whole_p[eChan1], list);
1208 static int do_cpt_enerhist(XDR *xd, gmx_bool bRead,
1209 int fflags, energyhistory_t *enerhist,
1220 enerhist->nsteps = 0;
1222 enerhist->nsteps_sim = 0;
1223 enerhist->nsum_sim = 0;
1224 enerhist->dht = NULL;
1226 if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
1228 snew(enerhist->dht, 1);
1229 enerhist->dht->ndh = NULL;
1230 enerhist->dht->dh = NULL;
1231 enerhist->dht->start_lambda_set = FALSE;
1235 for (i = 0; (i < eenhNR && ret == 0); i++)
1237 if (fflags & (1<<i))
1241 case eenhENERGY_N: ret = do_cpte_int(xd, cptpEENH, i, fflags, &enerhist->nener, list); break;
1242 case eenhENERGY_AVER: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_ave, list); break;
1243 case eenhENERGY_SUM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum, list); break;
1244 case eenhENERGY_NSUM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum, list); break;
1245 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum_sim, list); break;
1246 case eenhENERGY_NSUM_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum_sim, list); break;
1247 case eenhENERGY_NSTEPS: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps, list); break;
1248 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps_sim, list); break;
1249 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd, eenh_names[i], &(enerhist->dht->nndh), list);
1250 if (bRead) /* now allocate memory for it */
1252 snew(enerhist->dht->dh, enerhist->dht->nndh);
1253 snew(enerhist->dht->ndh, enerhist->dht->nndh);
1254 for (j = 0; j < enerhist->dht->nndh; j++)
1256 enerhist->dht->ndh[j] = 0;
1257 enerhist->dht->dh[j] = NULL;
1261 case eenhENERGY_DELTA_H_LIST:
1262 for (j = 0; j < enerhist->dht->nndh; j++)
1264 ret = do_cpte_n_reals(xd, cptpEENH, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list);
1267 case eenhENERGY_DELTA_H_STARTTIME:
1268 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_time), list); break;
1269 case eenhENERGY_DELTA_H_STARTLAMBDA:
1270 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_lambda), list); break;
1272 gmx_fatal(FARGS, "Unknown energy history entry %d\n"
1273 "You are probably reading a new checkpoint file with old code", i);
1278 if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
1280 /* Assume we have an old file format and copy sum to sum_sim */
1281 srenew(enerhist->ener_sum_sim, enerhist->nener);
1282 for (i = 0; i < enerhist->nener; i++)
1284 enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1286 fflags |= (1<<eenhENERGY_SUM_SIM);
1289 if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1290 !(fflags & (1<<eenhENERGY_NSTEPS)))
1292 /* Assume we have an old file format and copy nsum to nsteps */
1293 enerhist->nsteps = enerhist->nsum;
1294 fflags |= (1<<eenhENERGY_NSTEPS);
1296 if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1297 !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1299 /* Assume we have an old file format and copy nsum to nsteps */
1300 enerhist->nsteps_sim = enerhist->nsum_sim;
1301 fflags |= (1<<eenhENERGY_NSTEPS_SIM);
1307 static int do_cpt_df_hist(XDR *xd, int fflags, df_history_t *dfhist, FILE *list)
1312 nlambda = dfhist->nlambda;
1315 for (i = 0; (i < edfhNR && ret == 0); i++)
1317 if (fflags & (1<<i))
1321 case edfhBEQUIL: ret = do_cpte_int(xd, cptpEDFH, i, fflags, &dfhist->bEquil, list); break;
1322 case edfhNATLAMBDA: ret = do_cpte_ints(xd, cptpEDFH, i, fflags, nlambda, &dfhist->n_at_lam, list); break;
1323 case edfhWLHISTO: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->wl_histo, list); break;
1324 case edfhWLDELTA: ret = do_cpte_real(xd, cptpEDFH, i, fflags, &dfhist->wl_delta, list); break;
1325 case edfhSUMWEIGHTS: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_weights, list); break;
1326 case edfhSUMDG: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_dg, list); break;
1327 case edfhSUMMINVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_minvar, list); break;
1328 case edfhSUMVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_variance, list); break;
1329 case edfhACCUMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p, list); break;
1330 case edfhACCUMM: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m, list); break;
1331 case edfhACCUMP2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p2, list); break;
1332 case edfhACCUMM2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m2, list); break;
1333 case edfhTIJ: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij, list); break;
1334 case edfhTIJEMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij_empirical, list); break;
1337 gmx_fatal(FARGS, "Unknown df history entry %d\n"
1338 "You are probably reading a new checkpoint file with old code", i);
1347 /* This function stores the last whole configuration of the reference and
1348 * average structure in the .cpt file
1350 static int do_cpt_EDstate(XDR *xd, gmx_bool bRead,
1351 edsamstate_t *EDstate, FILE *list)
1358 EDstate->bFromCpt = bRead;
1360 if (EDstate->nED <= 0)
1365 /* When reading, init_edsam has not been called yet,
1366 * so we have to allocate memory first. */
1369 snew(EDstate->nref, EDstate->nED);
1370 snew(EDstate->old_sref, EDstate->nED);
1371 snew(EDstate->nav, EDstate->nED);
1372 snew(EDstate->old_sav, EDstate->nED);
1375 /* Read/write the last whole conformation of SREF and SAV for each ED dataset (usually only one) */
1376 for (i = 0; i < EDstate->nED; i++)
1378 /* Reference structure SREF */
1379 sprintf(buf, "ED%d # of atoms in reference structure", i+1);
1380 do_cpt_int_err(xd, buf, &EDstate->nref[i], list);
1381 sprintf(buf, "ED%d x_ref", i+1);
1384 snew(EDstate->old_sref[i], EDstate->nref[i]);
1385 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref[i], list);
1389 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref_p[i], list);
1392 /* Average structure SAV */
1393 sprintf(buf, "ED%d # of atoms in average structure", i+1);
1394 do_cpt_int_err(xd, buf, &EDstate->nav[i], list);
1395 sprintf(buf, "ED%d x_av", i+1);
1398 snew(EDstate->old_sav[i], EDstate->nav[i]);
1399 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav[i], list);
1403 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav_p[i], list);
1411 static int do_cpt_files(XDR *xd, gmx_bool bRead,
1412 gmx_file_position_t **p_outputfiles, int *nfiles,
1413 FILE *list, int file_version)
1417 gmx_off_t mask = 0xFFFFFFFFL;
1418 int offset_high, offset_low;
1420 gmx_file_position_t *outputfiles;
1422 if (do_cpt_int(xd, "number of output files", nfiles, list) != 0)
1429 snew(*p_outputfiles, *nfiles);
1432 outputfiles = *p_outputfiles;
1434 for (i = 0; i < *nfiles; i++)
1436 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1439 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1440 strncpy(outputfiles[i].filename, buf, CPTSTRLEN-1);
1446 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1450 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1454 outputfiles[i].offset = ( ((gmx_off_t) offset_high) << 32 ) | ( (gmx_off_t) offset_low & mask );
1458 buf = outputfiles[i].filename;
1459 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1461 offset = outputfiles[i].offset;
1469 offset_low = (int) (offset & mask);
1470 offset_high = (int) ((offset >> 32) & mask);
1472 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1476 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1481 if (file_version >= 8)
1483 if (do_cpt_int(xd, "file_checksum_size", &(outputfiles[i].chksum_size),
1488 if (do_cpt_u_chars(xd, "file_checksum", 16, outputfiles[i].chksum, list) != 0)
1495 outputfiles[i].chksum_size = -1;
1502 void write_checkpoint(const char *fn, gmx_bool bNumberAndKeep,
1503 FILE *fplog, t_commrec *cr,
1504 int eIntegrator, int simulation_part,
1505 gmx_bool bExpanded, int elamstats,
1506 gmx_int64_t step, double t, t_state *state)
1516 char *fntemp; /* the temporary checkpoint file name */
1518 char timebuf[STRLEN];
1519 int nppnodes, npmenodes, flag_64bit;
1520 char buf[1024], suffix[5+STEPSTRSIZE], sbuf[STEPSTRSIZE];
1521 gmx_file_position_t *outputfiles;
1524 int flags_eks, flags_enh, flags_dfh, i;
1527 if (DOMAINDECOMP(cr))
1529 nppnodes = cr->dd->nnodes;
1530 npmenodes = cr->npmenodes;
1538 #ifndef GMX_NO_RENAME
1539 /* make the new temporary filename */
1540 snew(fntemp, strlen(fn)+5+STEPSTRSIZE);
1542 fntemp[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1543 sprintf(suffix, "_%s%s", "step", gmx_step_str(step, sbuf));
1544 strcat(fntemp, suffix);
1545 strcat(fntemp, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1547 /* if we can't rename, we just overwrite the cpt file.
1548 * dangerous if interrupted.
1550 snew(fntemp, strlen(fn));
1554 gmx_ctime_r(&now, timebuf, STRLEN);
1558 fprintf(fplog, "Writing checkpoint, step %s at %s\n\n",
1559 gmx_step_str(step, buf), timebuf);
1562 /* Get offsets for open files */
1563 gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1565 fp = gmx_fio_open(fntemp, "w");
1567 if (state->ekinstate.bUpToDate)
1570 ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) |
1571 (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) |
1572 (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1580 if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1582 flags_enh |= (1<<eenhENERGY_N);
1583 if (state->enerhist.nsum > 0)
1585 flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1586 (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSUM));
1588 if (state->enerhist.nsum_sim > 0)
1590 flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSTEPS_SIM) |
1591 (1<<eenhENERGY_NSUM_SIM));
1593 if (state->enerhist.dht)
1595 flags_enh |= ( (1<< eenhENERGY_DELTA_H_NN) |
1596 (1<< eenhENERGY_DELTA_H_LIST) |
1597 (1<< eenhENERGY_DELTA_H_STARTTIME) |
1598 (1<< eenhENERGY_DELTA_H_STARTLAMBDA) );
1604 flags_dfh = ((1<<edfhBEQUIL) | (1<<edfhNATLAMBDA) | (1<<edfhSUMWEIGHTS) | (1<<edfhSUMDG) |
1605 (1<<edfhTIJ) | (1<<edfhTIJEMP));
1608 flags_dfh |= ((1<<edfhWLDELTA) | (1<<edfhWLHISTO));
1610 if ((elamstats == elamstatsMINVAR) || (elamstats == elamstatsBARKER) || (elamstats == elamstatsMETROPOLIS))
1612 flags_dfh |= ((1<<edfhACCUMP) | (1<<edfhACCUMM) | (1<<edfhACCUMP2) | (1<<edfhACCUMM2)
1613 | (1<<edfhSUMMINVAR) | (1<<edfhSUMVAR));
1621 /* We can check many more things now (CPU, acceleration, etc), but
1622 * it is highly unlikely to have two separate builds with exactly
1623 * the same version, user, time, and build host!
1626 version = gmx_strdup(gmx_version());
1627 btime = gmx_strdup(BUILD_TIME);
1628 buser = gmx_strdup(BUILD_USER);
1629 bhost = gmx_strdup(BUILD_HOST);
1631 double_prec = GMX_CPT_BUILD_DP;
1632 fprog = gmx_strdup(Program());
1634 ftime = &(timebuf[0]);
1636 do_cpt_header(gmx_fio_getxdr(fp), FALSE, &file_version,
1637 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1638 &eIntegrator, &simulation_part, &step, &t, &nppnodes,
1639 DOMAINDECOMP(cr) ? cr->dd->nc : NULL, &npmenodes,
1640 &state->natoms, &state->ngtc, &state->nnhpres,
1641 &state->nhchainlength, &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
1642 &state->edsamstate.nED, &state->swapstate.eSwapCoords,
1651 if ((do_cpt_state(gmx_fio_getxdr(fp), FALSE, state->flags, state, NULL) < 0) ||
1652 (do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL) < 0) ||
1653 (do_cpt_enerhist(gmx_fio_getxdr(fp), FALSE, flags_enh, &state->enerhist, NULL) < 0) ||
1654 (do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL) < 0) ||
1655 (do_cpt_EDstate(gmx_fio_getxdr(fp), FALSE, &state->edsamstate, NULL) < 0) ||
1656 (do_cpt_swapstate(gmx_fio_getxdr(fp), FALSE, &state->swapstate, NULL) < 0) ||
1657 (do_cpt_files(gmx_fio_getxdr(fp), FALSE, &outputfiles, &noutputfiles, NULL,
1660 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1663 do_cpt_footer(gmx_fio_getxdr(fp), file_version);
1665 /* we really, REALLY, want to make sure to physically write the checkpoint,
1666 and all the files it depends on, out to disk. Because we've
1667 opened the checkpoint with gmx_fio_open(), it's in our list
1669 ret = gmx_fio_all_output_fsync();
1675 "Cannot fsync '%s'; maybe you are out of disk space?",
1676 gmx_fio_getname(ret));
1678 if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV) == NULL)
1688 if (gmx_fio_close(fp) != 0)
1690 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1693 /* we don't move the checkpoint if the user specified they didn't want it,
1694 or if the fsyncs failed */
1695 #ifndef GMX_NO_RENAME
1696 if (!bNumberAndKeep && !ret)
1700 /* Rename the previous checkpoint file */
1702 buf[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1703 strcat(buf, "_prev");
1704 strcat(buf, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1706 /* we copy here so that if something goes wrong between now and
1707 * the rename below, there's always a state.cpt.
1708 * If renames are atomic (such as in POSIX systems),
1709 * this copying should be unneccesary.
1711 gmx_file_copy(fn, buf, FALSE);
1712 /* We don't really care if this fails:
1713 * there's already a new checkpoint.
1716 gmx_file_rename(fn, buf);
1719 if (gmx_file_rename(fntemp, fn) != 0)
1721 gmx_file("Cannot rename checkpoint file; maybe you are out of disk space?");
1724 #endif /* GMX_NO_RENAME */
1730 /*code for alternate checkpointing scheme. moved from top of loop over
1732 fcRequestCheckPoint();
1733 if (fcCheckPointParallel( cr->nodeid, NULL, 0) == 0)
1735 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", step );
1737 #endif /* end GMX_FAHCORE block */
1740 static void print_flag_mismatch(FILE *fplog, int sflags, int fflags)
1744 fprintf(fplog, "\nState entry mismatch between the simulation and the checkpoint file\n");
1745 fprintf(fplog, "Entries which are not present in the checkpoint file will not be updated\n");
1746 fprintf(fplog, " %24s %11s %11s\n", "", "simulation", "checkpoint");
1747 for (i = 0; i < estNR; i++)
1749 if ((sflags & (1<<i)) || (fflags & (1<<i)))
1751 fprintf(fplog, " %24s %11s %11s\n",
1753 (sflags & (1<<i)) ? " present " : "not present",
1754 (fflags & (1<<i)) ? " present " : "not present");
1759 static void check_int(FILE *fplog, const char *type, int p, int f, gmx_bool *mm)
1761 FILE *fp = fplog ? fplog : stderr;
1765 fprintf(fp, " %s mismatch,\n", type);
1766 fprintf(fp, " current program: %d\n", p);
1767 fprintf(fp, " checkpoint file: %d\n", f);
1773 static void check_string(FILE *fplog, const char *type, const char *p,
1774 const char *f, gmx_bool *mm)
1776 FILE *fp = fplog ? fplog : stderr;
1778 if (strcmp(p, f) != 0)
1780 fprintf(fp, " %s mismatch,\n", type);
1781 fprintf(fp, " current program: %s\n", p);
1782 fprintf(fp, " checkpoint file: %s\n", f);
1788 static void check_match(FILE *fplog,
1790 char *btime, char *buser, char *bhost, int double_prec,
1792 t_commrec *cr, int npp_f, int npme_f,
1793 ivec dd_nc, ivec dd_nc_f)
1796 gmx_bool mm = FALSE;
1797 gmx_bool patchlevel_differs = FALSE;
1798 gmx_bool version_differs = FALSE;
1800 check_string(fplog, "Version", gmx_version(), version, &mm);
1801 patchlevel_differs = mm;
1803 if (patchlevel_differs)
1805 /* Gromacs should be able to continue from checkpoints between
1806 * different patch level versions, but we do not guarantee
1807 * compatibility between different major/minor versions - check this.
1809 int gmx_major, gmx_minor;
1810 int cpt_major, cpt_minor;
1811 sscanf(gmx_version(), "VERSION %d.%d", &gmx_major, &gmx_minor);
1812 sscanf(version, "VERSION %d.%d", &cpt_major, &cpt_minor);
1813 version_differs = (gmx_major != cpt_major || gmx_minor != cpt_minor);
1816 check_string(fplog, "Build time", BUILD_TIME, btime, &mm);
1817 check_string(fplog, "Build user", BUILD_USER, buser, &mm);
1818 check_string(fplog, "Build host", BUILD_HOST, bhost, &mm);
1819 check_int (fplog, "Double prec.", GMX_CPT_BUILD_DP, double_prec, &mm);
1820 check_string(fplog, "Program name", Program(), fprog, &mm);
1822 check_int (fplog, "#ranks", cr->nnodes, npp_f+npme_f, &mm);
1825 check_int (fplog, "#PME-ranks", cr->npmenodes, npme_f, &mm);
1828 if (cr->npmenodes >= 0)
1830 npp -= cr->npmenodes;
1834 check_int (fplog, "#DD-cells[x]", dd_nc[XX], dd_nc_f[XX], &mm);
1835 check_int (fplog, "#DD-cells[y]", dd_nc[YY], dd_nc_f[YY], &mm);
1836 check_int (fplog, "#DD-cells[z]", dd_nc[ZZ], dd_nc_f[ZZ], &mm);
1842 const char msg_version_difference[] =
1843 "The current Gromacs major & minor version are not identical to those that\n"
1844 "generated the checkpoint file. In principle Gromacs does not support\n"
1845 "continuation from checkpoints between different versions, so we advise\n"
1846 "against this. If you still want to try your luck we recommend that you use\n"
1847 "the -noappend flag to keep your output files from the two versions separate.\n"
1848 "This might also work around errors where the output fields in the energy\n"
1849 "file have changed between the different major & minor versions.\n";
1851 const char msg_mismatch_notice[] =
1852 "Gromacs patchlevel, binary or parallel settings differ from previous run.\n"
1853 "Continuation is exact, but not guaranteed to be binary identical.\n";
1855 const char msg_logdetails[] =
1856 "See the log file for details.\n";
1858 if (version_differs)
1860 fprintf(stderr, "%s%s\n", msg_version_difference, fplog ? msg_logdetails : "");
1864 fprintf(fplog, "%s\n", msg_version_difference);
1869 /* Major & minor versions match at least, but something is different. */
1870 fprintf(stderr, "%s%s\n", msg_mismatch_notice, fplog ? msg_logdetails : "");
1873 fprintf(fplog, "%s\n", msg_mismatch_notice);
1879 static void read_checkpoint(const char *fn, FILE **pfplog,
1880 t_commrec *cr, ivec dd_nc,
1881 int eIntegrator, int *init_fep_state, gmx_int64_t *step, double *t,
1882 t_state *state, gmx_bool *bReadEkin,
1883 int *simulation_part,
1884 gmx_bool bAppendOutputFiles, gmx_bool bForceAppend)
1889 char *version, *btime, *buser, *bhost, *fprog, *ftime;
1891 char filename[STRLEN], buf[STEPSTRSIZE];
1892 int nppnodes, eIntegrator_f, nppnodes_f, npmenodes_f;
1894 int natoms, ngtc, nnhpres, nhchainlength, nlambda, fflags, flags_eks, flags_enh, flags_dfh;
1897 gmx_file_position_t *outputfiles;
1899 t_fileio *chksum_file;
1900 FILE * fplog = *pfplog;
1901 unsigned char digest[16];
1902 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1903 struct flock fl; /* don't initialize here: the struct order is OS
1907 const char *int_warn =
1908 "WARNING: The checkpoint file was generated with integrator %s,\n"
1909 " while the simulation uses integrator %s\n\n";
1911 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1912 fl.l_type = F_WRLCK;
1913 fl.l_whence = SEEK_SET;
1919 fp = gmx_fio_open(fn, "r");
1920 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
1921 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1922 &eIntegrator_f, simulation_part, step, t,
1923 &nppnodes_f, dd_nc_f, &npmenodes_f,
1924 &natoms, &ngtc, &nnhpres, &nhchainlength, &nlambda,
1925 &fflags, &flags_eks, &flags_enh, &flags_dfh,
1926 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
1928 if (bAppendOutputFiles &&
1929 file_version >= 13 && double_prec != GMX_CPT_BUILD_DP)
1931 gmx_fatal(FARGS, "Output file appending requested, but the code and checkpoint file precision (single/double) don't match");
1934 if (cr == NULL || MASTER(cr))
1936 fprintf(stderr, "\nReading checkpoint file %s generated: %s\n\n",
1940 /* This will not be written if we do appending, since fplog is still NULL then */
1943 fprintf(fplog, "\n");
1944 fprintf(fplog, "Reading checkpoint file %s\n", fn);
1945 fprintf(fplog, " file generated by: %s\n", fprog);
1946 fprintf(fplog, " file generated at: %s\n", ftime);
1947 fprintf(fplog, " GROMACS build time: %s\n", btime);
1948 fprintf(fplog, " GROMACS build user: %s\n", buser);
1949 fprintf(fplog, " GROMACS build host: %s\n", bhost);
1950 fprintf(fplog, " GROMACS double prec.: %d\n", double_prec);
1951 fprintf(fplog, " simulation part #: %d\n", *simulation_part);
1952 fprintf(fplog, " step: %s\n", gmx_step_str(*step, buf));
1953 fprintf(fplog, " time: %f\n", *t);
1954 fprintf(fplog, "\n");
1957 if (natoms != state->natoms)
1959 gmx_fatal(FARGS, "Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms", natoms, state->natoms);
1961 if (ngtc != state->ngtc)
1963 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);
1965 if (nnhpres != state->nnhpres)
1967 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);
1970 if (nlambda != state->dfhist.nlambda)
1972 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);
1975 init_gtc_state(state, state->ngtc, state->nnhpres, nhchainlength); /* need to keep this here to keep the tpr format working */
1976 /* write over whatever was read; we use the number of Nose-Hoover chains from the checkpoint */
1978 if (eIntegrator_f != eIntegrator)
1982 fprintf(stderr, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1984 if (bAppendOutputFiles)
1987 "Output file appending requested, but input/checkpoint integrators do not match.\n"
1988 "Stopping the run to prevent you from ruining all your data...\n"
1989 "If you _really_ know what you are doing, try with the -noappend option.\n");
1993 fprintf(fplog, int_warn, EI(eIntegrator_f), EI(eIntegrator));
2002 else if (cr->nnodes == nppnodes_f + npmenodes_f)
2004 if (cr->npmenodes < 0)
2006 cr->npmenodes = npmenodes_f;
2008 nppnodes = cr->nnodes - cr->npmenodes;
2009 if (nppnodes == nppnodes_f)
2011 for (d = 0; d < DIM; d++)
2015 dd_nc[d] = dd_nc_f[d];
2022 /* The number of PP nodes has not been set yet */
2026 if (fflags != state->flags)
2031 if (bAppendOutputFiles)
2034 "Output file appending requested, but input and checkpoint states are not identical.\n"
2035 "Stopping the run to prevent you from ruining all your data...\n"
2036 "You can try with the -noappend option, and get more info in the log file.\n");
2039 if (getenv("GMX_ALLOW_CPT_MISMATCH") == NULL)
2041 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");
2046 "WARNING: The checkpoint state entries do not match the simulation,\n"
2047 " see the log file for details\n\n");
2053 print_flag_mismatch(fplog, state->flags, fflags);
2060 check_match(fplog, version, btime, buser, bhost, double_prec, fprog,
2061 cr, nppnodes_f, npmenodes_f, dd_nc, dd_nc_f);
2064 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, fflags, state, NULL);
2065 *init_fep_state = state->fep_state; /* there should be a better way to do this than setting it here.
2066 Investigate for 5.0. */
2071 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2076 *bReadEkin = ((flags_eks & (1<<eeksEKINH)) || (flags_eks & (1<<eeksEKINF)) || (flags_eks & (1<<eeksEKINO)) ||
2077 ((flags_eks & (1<<eeksEKINSCALEF)) | (flags_eks & (1<<eeksEKINSCALEH)) | (flags_eks & (1<<eeksVSCALE))));
2079 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2080 flags_enh, &state->enerhist, NULL);
2086 if (file_version < 6)
2088 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.";
2090 fprintf(stderr, "\nWARNING: %s\n\n", warn);
2093 fprintf(fplog, "\nWARNING: %s\n\n", warn);
2095 state->enerhist.nsum = *step;
2096 state->enerhist.nsum_sim = *step;
2099 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2105 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2111 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2117 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, NULL, file_version);
2123 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2128 if (gmx_fio_close(fp) != 0)
2130 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2139 /* If the user wants to append to output files,
2140 * we use the file pointer positions of the output files stored
2141 * in the checkpoint file and truncate the files such that any frames
2142 * written after the checkpoint time are removed.
2143 * All files are md5sum checked such that we can be sure that
2144 * we do not truncate other (maybe imprortant) files.
2146 if (bAppendOutputFiles)
2148 if (fn2ftp(outputfiles[0].filename) != efLOG)
2150 /* make sure first file is log file so that it is OK to use it for
2153 gmx_fatal(FARGS, "The first output file should always be the log "
2154 "file but instead is: %s. Cannot do appending because of this condition.", outputfiles[0].filename);
2156 for (i = 0; i < nfiles; i++)
2158 if (outputfiles[i].offset < 0)
2160 gmx_fatal(FARGS, "The original run wrote a file called '%s' which "
2161 "is larger than 2 GB, but mdrun did not support large file"
2162 " offsets. Can not append. Run mdrun with -noappend",
2163 outputfiles[i].filename);
2166 chksum_file = gmx_fio_open(outputfiles[i].filename, "a");
2169 chksum_file = gmx_fio_open(outputfiles[i].filename, "r+");
2174 /* Note that there are systems where the lock operation
2175 * will succeed, but a second process can also lock the file.
2176 * We should probably try to detect this.
2178 #if defined __native_client__
2182 #elif defined GMX_NATIVE_WINDOWS
2183 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX) == -1)
2185 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl) == -1)
2188 if (errno == ENOSYS)
2192 gmx_fatal(FARGS, "File locking is not supported on this system. Use -noappend or specify -append explicitly to append anyhow.");
2196 fprintf(stderr, "\nNOTE: File locking is not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2199 fprintf(fplog, "\nNOTE: File locking not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2203 else if (errno == EACCES || errno == EAGAIN)
2205 gmx_fatal(FARGS, "Failed to lock: %s. Already running "
2206 "simulation?", outputfiles[i].filename);
2210 gmx_fatal(FARGS, "Failed to lock: %s. %s.",
2211 outputfiles[i].filename, strerror(errno));
2216 /* compute md5 chksum */
2217 if (outputfiles[i].chksum_size != -1)
2219 if (gmx_fio_get_file_md5(chksum_file, outputfiles[i].offset,
2220 digest) != outputfiles[i].chksum_size) /*at the end of the call the file position is at the end of the file*/
2222 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.",
2223 outputfiles[i].chksum_size,
2224 outputfiles[i].filename);
2227 if (i == 0) /*log file needs to be seeked in case we need to truncate (other files are truncated below)*/
2229 if (gmx_fio_seek(chksum_file, outputfiles[i].offset))
2231 gmx_fatal(FARGS, "Seek error! Failed to truncate log-file: %s.", strerror(errno));
2236 if (i == 0) /*open log file here - so that lock is never lifted
2237 after chksum is calculated */
2239 *pfplog = gmx_fio_getfp(chksum_file);
2243 gmx_fio_close(chksum_file);
2246 /* compare md5 chksum */
2247 if (outputfiles[i].chksum_size != -1 &&
2248 memcmp(digest, outputfiles[i].chksum, 16) != 0)
2252 fprintf(debug, "chksum for %s: ", outputfiles[i].filename);
2253 for (j = 0; j < 16; j++)
2255 fprintf(debug, "%02x", digest[j]);
2257 fprintf(debug, "\n");
2259 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.",
2260 outputfiles[i].filename);
2265 if (i != 0) /*log file is already seeked to correct position */
2267 #ifdef GMX_NATIVE_WINDOWS
2268 rc = gmx_wintruncate(outputfiles[i].filename, outputfiles[i].offset);
2270 rc = truncate(outputfiles[i].filename, outputfiles[i].offset);
2274 gmx_fatal(FARGS, "Truncation of file %s failed. Cannot do appending because of this failure.", outputfiles[i].filename);
2284 void load_checkpoint(const char *fn, FILE **fplog,
2285 t_commrec *cr, ivec dd_nc,
2286 t_inputrec *ir, t_state *state,
2287 gmx_bool *bReadEkin,
2288 gmx_bool bAppend, gmx_bool bForceAppend)
2295 /* Read the state from the checkpoint file */
2296 read_checkpoint(fn, fplog,
2298 ir->eI, &(ir->fepvals->init_fep_state), &step, &t, state, bReadEkin,
2299 &ir->simulation_part, bAppend, bForceAppend);
2303 gmx_bcast(sizeof(cr->npmenodes), &cr->npmenodes, cr);
2304 gmx_bcast(DIM*sizeof(dd_nc[0]), dd_nc, cr);
2305 gmx_bcast(sizeof(step), &step, cr);
2306 gmx_bcast(sizeof(*bReadEkin), bReadEkin, cr);
2308 ir->bContinuation = TRUE;
2309 if (ir->nsteps >= 0)
2311 ir->nsteps += ir->init_step - step;
2313 ir->init_step = step;
2314 ir->simulation_part += 1;
2317 static void read_checkpoint_data(t_fileio *fp, int *simulation_part,
2318 gmx_int64_t *step, double *t, t_state *state,
2319 int *nfiles, gmx_file_position_t **outputfiles)
2322 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2327 int flags_eks, flags_enh, flags_dfh;
2329 gmx_file_position_t *files_loc = NULL;
2332 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2333 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2334 &eIntegrator, simulation_part, step, t, &nppnodes, dd_nc, &npme,
2335 &state->natoms, &state->ngtc, &state->nnhpres, &state->nhchainlength,
2336 &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
2337 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
2339 do_cpt_state(gmx_fio_getxdr(fp), TRUE, state->flags, state, NULL);
2344 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2349 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2350 flags_enh, &state->enerhist, NULL);
2355 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2361 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2367 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2373 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE,
2374 outputfiles != NULL ? outputfiles : &files_loc,
2375 outputfiles != NULL ? nfiles : &nfiles_loc,
2376 NULL, file_version);
2377 if (files_loc != NULL)
2387 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2401 read_checkpoint_state(const char *fn, int *simulation_part,
2402 gmx_int64_t *step, double *t, t_state *state)
2406 fp = gmx_fio_open(fn, "r");
2407 read_checkpoint_data(fp, simulation_part, step, t, state, NULL, NULL);
2408 if (gmx_fio_close(fp) != 0)
2410 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2414 void read_checkpoint_trxframe(t_fileio *fp, t_trxframe *fr)
2416 /* This next line is nasty because the sub-structures of t_state
2417 * cannot be assumed to be zeroed (or even initialized in ways the
2418 * rest of the code might assume). Using snew would be better, but
2419 * this will all go away for 5.0. */
2421 int simulation_part;
2425 init_state(&state, 0, 0, 0, 0, 0);
2427 read_checkpoint_data(fp, &simulation_part, &step, &t, &state, NULL, NULL);
2429 fr->natoms = state.natoms;
2432 fr->step = gmx_int64_to_int(step,
2433 "conversion of checkpoint to trajectory");
2437 fr->lambda = state.lambda[efptFEP];
2438 fr->fep_state = state.fep_state;
2440 fr->bX = (state.flags & (1<<estX));
2446 fr->bV = (state.flags & (1<<estV));
2453 fr->bBox = (state.flags & (1<<estBOX));
2456 copy_mat(state.box, fr->box);
2461 void list_checkpoint(const char *fn, FILE *out)
2465 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2467 int eIntegrator, simulation_part, nppnodes, npme;
2472 int flags_eks, flags_enh, flags_dfh;
2476 gmx_file_position_t *outputfiles;
2479 init_state(&state, -1, -1, -1, -1, 0);
2481 fp = gmx_fio_open(fn, "r");
2482 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2483 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2484 &eIntegrator, &simulation_part, &step, &t, &nppnodes, dd_nc, &npme,
2485 &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
2486 &(state.dfhist.nlambda), &state.flags,
2487 &flags_eks, &flags_enh, &flags_dfh, &state.edsamstate.nED,
2488 &state.swapstate.eSwapCoords, out);
2489 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, state.flags, &state, out);
2494 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state.ekinstate, out);
2499 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2500 flags_enh, &state.enerhist, out);
2504 ret = do_cpt_df_hist(gmx_fio_getxdr(fp),
2505 flags_dfh, &state.dfhist, out);
2510 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state.edsamstate, out);
2515 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state.swapstate, out);
2520 do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, out, file_version);
2525 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2532 if (gmx_fio_close(fp) != 0)
2534 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2541 static gmx_bool exist_output_file(const char *fnm_cp, int nfile, const t_filenm fnm[])
2545 /* Check if the output file name stored in the checkpoint file
2546 * is one of the output file names of mdrun.
2550 !(is_output(&fnm[i]) && strcmp(fnm_cp, fnm[i].fns[0]) == 0))
2555 return (i < nfile && gmx_fexist(fnm_cp));
2558 /* This routine cannot print tons of data, since it is called before the log file is opened. */
2559 gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
2560 gmx_int64_t *cpt_step, t_commrec *cr,
2561 gmx_bool bAppendReq,
2562 int nfile, const t_filenm fnm[],
2563 const char *part_suffix, gmx_bool *bAddPart)
2566 gmx_int64_t step = 0;
2568 /* This next line is nasty because the sub-structures of t_state
2569 * cannot be assumed to be zeroed (or even initialized in ways the
2570 * rest of the code might assume). Using snew would be better, but
2571 * this will all go away for 5.0. */
2574 gmx_file_position_t *outputfiles;
2577 char *fn, suf_up[STRLEN];
2583 if (!gmx_fexist(filename) || (!(fp = gmx_fio_open(filename, "r")) ))
2585 *simulation_part = 0;
2589 init_state(&state, 0, 0, 0, 0, 0);
2591 read_checkpoint_data(fp, simulation_part, &step, &t, &state,
2592 &nfiles, &outputfiles);
2593 if (gmx_fio_close(fp) != 0)
2595 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2602 for (f = 0; f < nfiles; f++)
2604 if (exist_output_file(outputfiles[f].filename, nfile, fnm))
2609 if (nexist == nfiles)
2611 bAppend = bAppendReq;
2613 else if (nexist > 0)
2616 "Output file appending has been requested,\n"
2617 "but some output files listed in the checkpoint file %s\n"
2618 "are not present or are named differently by the current program:\n",
2620 fprintf(stderr, "output files present:");
2621 for (f = 0; f < nfiles; f++)
2623 if (exist_output_file(outputfiles[f].filename,
2626 fprintf(stderr, " %s", outputfiles[f].filename);
2629 fprintf(stderr, "\n");
2630 fprintf(stderr, "output files not present or named differently:");
2631 for (f = 0; f < nfiles; f++)
2633 if (!exist_output_file(outputfiles[f].filename,
2636 fprintf(stderr, " %s", outputfiles[f].filename);
2639 fprintf(stderr, "\n");
2641 gmx_fatal(FARGS, "File appending requested, but %d of the %d output files are not present or are named differently", nfiles-nexist, nfiles);
2649 gmx_fatal(FARGS, "File appending requested, but no output file information is stored in the checkpoint file");
2651 fn = outputfiles[0].filename;
2652 if (strlen(fn) < 4 ||
2653 gmx_strcasecmp(fn+strlen(fn)-4, ftp2ext(efLOG)) == 0)
2655 gmx_fatal(FARGS, "File appending requested, but the log file is not the first file listed in the checkpoint file");
2657 /* Set bAddPart to whether the suffix string '.part' is present
2658 * in the log file name.
2660 strcpy(suf_up, part_suffix);
2662 *bAddPart = (strstr(fn, part_suffix) != NULL ||
2663 strstr(fn, suf_up) != NULL);
2671 gmx_bcast(sizeof(*simulation_part), simulation_part, cr);
2673 if (*simulation_part > 0 && bAppendReq)
2675 gmx_bcast(sizeof(bAppend), &bAppend, cr);
2676 gmx_bcast(sizeof(*bAddPart), bAddPart, cr);
2679 if (NULL != cpt_step)