2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2008,2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 /* The source code in this file should be thread-safe.
37 Please keep it that way. */
50 #ifdef HAVE_SYS_TIME_H
58 #ifdef GMX_NATIVE_WINDOWS
61 #include <sys/locking.h>
67 #include "types/commrec.h"
68 #include "gromacs/utility/smalloc.h"
70 #include "gromacs/math/vec.h"
72 #include "checkpoint.h"
73 #include "gromacs/utility/cstringutil.h"
75 #include "gromacs/fileio/filenm.h"
76 #include "gromacs/utility/futil.h"
77 #include "gromacs/fileio/gmxfio.h"
78 #include "gromacs/fileio/xdrf.h"
79 #include "gromacs/fileio/xdr_datatype.h"
80 #include "gromacs/utility/basenetwork.h"
81 #include "gromacs/utility/baseversion.h"
82 #include "gromacs/utility/fatalerror.h"
84 #include "buildinfo.h"
90 #define CPT_MAGIC1 171817
91 #define CPT_MAGIC2 171819
92 #define CPTSTRLEN 1024
95 #define GMX_CPT_BUILD_DP 1
97 #define GMX_CPT_BUILD_DP 0
100 /* cpt_version should normally only be changed
101 * when the header of footer format changes.
102 * The state data format itself is backward and forward compatible.
103 * But old code can not read a new entry that is present in the file
104 * (but can read a new format when new entries are not present).
106 static const int cpt_version = 16;
109 const char *est_names[estNR] =
112 "box", "box-rel", "box-v", "pres_prev",
113 "nosehoover-xi", "thermostat-integral",
114 "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
115 "disre_initf", "disre_rm3tav",
116 "orire_initf", "orire_Dtav",
117 "svir_prev", "nosehoover-vxi", "v_eta", "vol0", "nhpres_xi", "nhpres_vxi", "fvir_prev", "fep_state", "MC-rng", "MC-rng-i"
121 eeksEKIN_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksEKINF, eeksEKINO, eeksEKINSCALEF, eeksEKINSCALEH, eeksVSCALE, eeksEKINTOTAL, eeksNR
124 const char *eeks_names[eeksNR] =
126 "Ekin_n", "Ekinh", "dEkindlambda", "mv_cos",
127 "Ekinf", "Ekinh_old", "EkinScaleF_NHC", "EkinScaleH_NHC", "Vscale_NHC", "Ekin_Total"
131 eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
132 eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
133 eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM,
134 eenhENERGY_DELTA_H_NN,
135 eenhENERGY_DELTA_H_LIST,
136 eenhENERGY_DELTA_H_STARTTIME,
137 eenhENERGY_DELTA_H_STARTLAMBDA,
141 const char *eenh_names[eenhNR] =
143 "energy_n", "energy_aver", "energy_sum", "energy_nsum",
144 "energy_sum_sim", "energy_nsum_sim",
145 "energy_nsteps", "energy_nsteps_sim",
147 "energy_delta_h_list",
148 "energy_delta_h_start_time",
149 "energy_delta_h_start_lambda"
152 /* free energy history variables -- need to be preserved over checkpoint */
154 edfhBEQUIL, edfhNATLAMBDA, edfhWLHISTO, edfhWLDELTA, edfhSUMWEIGHTS, edfhSUMDG, edfhSUMMINVAR, edfhSUMVAR,
155 edfhACCUMP, edfhACCUMM, edfhACCUMP2, edfhACCUMM2, edfhTIJ, edfhTIJEMP, edfhNR
157 /* free energy history variable names */
158 const char *edfh_names[edfhNR] =
160 "bEquilibrated", "N_at_state", "Wang-Landau Histogram", "Wang-Landau Delta", "Weights", "Free Energies", "minvar", "variance",
161 "accumulated_plus", "accumulated_minus", "accumulated_plus_2", "accumulated_minus_2", "Tij", "Tij_empirical"
164 #ifdef GMX_NATIVE_WINDOWS
166 gmx_wintruncate(const char *filename, __int64 size)
169 /*we do this elsewhere*/
175 fp = fopen(filename, "rb+");
182 return _chsize_s( fileno(fp), size);
189 ecprREAL, ecprRVEC, ecprMATRIX
193 cptpEST, cptpEEKS, cptpEENH, cptpEDFH
195 /* enums for the different components of checkpoint variables, replacing the hard coded ones.
196 cptpEST - state variables.
197 cptpEEKS - Kinetic energy state variables.
198 cptpEENH - Energy history state variables.
199 cptpEDFH - free energy history variables.
203 static const char *st_names(int cptp, int ecpt)
207 case cptpEST: return est_names [ecpt]; break;
208 case cptpEEKS: return eeks_names[ecpt]; break;
209 case cptpEENH: return eenh_names[ecpt]; break;
210 case cptpEDFH: return edfh_names[ecpt]; break;
216 static void cp_warning(FILE *fp)
218 fprintf(fp, "\nWARNING: Checkpoint file is corrupted or truncated\n\n");
221 static void cp_error()
223 gmx_fatal(FARGS, "Checkpoint file corrupted/truncated, or maybe you are out of disk space?");
226 static void do_cpt_string_err(XDR *xd, gmx_bool bRead, const char *desc, char **s, FILE *list)
234 res = xdr_string(xd, s, CPTSTRLEN);
241 fprintf(list, "%s = %s\n", desc, *s);
246 static int do_cpt_int(XDR *xd, const char *desc, int *i, FILE *list)
250 res = xdr_int(xd, i);
257 fprintf(list, "%s = %d\n", desc, *i);
262 static int do_cpt_u_chars(XDR *xd, const char *desc, int n, unsigned char *i, FILE *list)
268 fprintf(list, "%s = ", desc);
270 for (j = 0; j < n && res; j++)
272 res &= xdr_u_char(xd, &i[j]);
275 fprintf(list, "%02x", i[j]);
290 static void do_cpt_int_err(XDR *xd, const char *desc, int *i, FILE *list)
292 if (do_cpt_int(xd, desc, i, list) < 0)
298 static void do_cpt_step_err(XDR *xd, const char *desc, gmx_int64_t *i, FILE *list)
301 char buf[STEPSTRSIZE];
303 res = xdr_int64(xd, i);
310 fprintf(list, "%s = %s\n", desc, gmx_step_str(*i, buf));
314 static void do_cpt_double_err(XDR *xd, const char *desc, double *f, FILE *list)
318 res = xdr_double(xd, f);
325 fprintf(list, "%s = %f\n", desc, *f);
329 static void do_cpt_real_err(XDR *xd, real *f)
334 res = xdr_double(xd, f);
336 res = xdr_float(xd, f);
344 static void do_cpt_n_rvecs_err(XDR *xd, const char *desc, int n, rvec f[], FILE *list)
348 for (i = 0; i < n; i++)
350 for (j = 0; j < DIM; j++)
352 do_cpt_real_err(xd, &f[i][j]);
358 pr_rvecs(list, 0, desc, f, n);
362 /* If nval >= 0, nval is used; on read this should match the passed value.
363 * If nval n<0, *nptr is used; on read the value is stored in nptr
365 static int do_cpte_reals_low(XDR *xd, int cptp, int ecpt, int sflags,
366 int nval, int *nptr, real **v,
367 FILE *list, int erealtype)
371 int dtc = xdr_datatype_float;
373 int dtc = xdr_datatype_double;
375 real *vp, *va = NULL;
390 gmx_incons("*ntpr=NULL in do_cpte_reals_low");
395 res = xdr_int(xd, &nf);
406 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), nval, nf);
415 res = xdr_int(xd, &dt);
422 fprintf(stderr, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
423 st_names(cptp, ecpt), xdr_datatype_names[dtc],
424 xdr_datatype_names[dt]);
426 if (list || !(sflags & (1<<ecpt)))
439 if (dt == xdr_datatype_float)
441 if (dtc == xdr_datatype_float)
449 res = xdr_vector(xd, (char *)vf, nf,
450 (unsigned int)sizeof(float), (xdrproc_t)xdr_float);
455 if (dtc != xdr_datatype_float)
457 for (i = 0; i < nf; i++)
466 if (dtc == xdr_datatype_double)
474 res = xdr_vector(xd, (char *)vd, nf,
475 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
480 if (dtc != xdr_datatype_double)
482 for (i = 0; i < nf; i++)
495 pr_reals(list, 0, st_names(cptp, ecpt), vp, nf);
498 pr_rvecs(list, 0, st_names(cptp, ecpt), (rvec *)vp, nf/3);
501 gmx_incons("Unknown checkpoint real type");
513 /* This function stores n along with the reals for reading,
514 * but on reading it assumes that n matches the value in the checkpoint file,
515 * a fatal error is generated when this is not the case.
517 static int do_cpte_reals(XDR *xd, int cptp, int ecpt, int sflags,
518 int n, real **v, FILE *list)
520 return do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, v, list, ecprREAL);
523 /* This function does the same as do_cpte_reals,
524 * except that on reading it ignores the passed value of *n
525 * and stored the value read from the checkpoint file in *n.
527 static int do_cpte_n_reals(XDR *xd, int cptp, int ecpt, int sflags,
528 int *n, real **v, FILE *list)
530 return do_cpte_reals_low(xd, cptp, ecpt, sflags, -1, n, v, list, ecprREAL);
533 static int do_cpte_real(XDR *xd, int cptp, int ecpt, int sflags,
538 return do_cpte_reals_low(xd, cptp, ecpt, sflags, 1, NULL, &r, list, ecprREAL);
541 static int do_cpte_ints(XDR *xd, int cptp, int ecpt, int sflags,
542 int n, int **v, FILE *list)
545 int dtc = xdr_datatype_int;
550 res = xdr_int(xd, &nf);
555 if (list == NULL && v != NULL && nf != n)
557 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
560 res = xdr_int(xd, &dt);
567 gmx_fatal(FARGS, "Type mismatch for state entry %s, code type is %s, file type is %s\n",
568 st_names(cptp, ecpt), xdr_datatype_names[dtc],
569 xdr_datatype_names[dt]);
571 if (list || !(sflags & (1<<ecpt)) || v == NULL)
584 res = xdr_vector(xd, (char *)vp, nf,
585 (unsigned int)sizeof(int), (xdrproc_t)xdr_int);
592 pr_ivec(list, 0, st_names(cptp, ecpt), vp, nf, TRUE);
602 static int do_cpte_int(XDR *xd, int cptp, int ecpt, int sflags,
605 return do_cpte_ints(xd, cptp, ecpt, sflags, 1, &i, list);
608 static int do_cpte_doubles(XDR *xd, int cptp, int ecpt, int sflags,
609 int n, double **v, FILE *list)
612 int dtc = xdr_datatype_double;
613 double *vp, *va = NULL;
617 res = xdr_int(xd, &nf);
622 if (list == NULL && nf != n)
624 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
627 res = xdr_int(xd, &dt);
634 gmx_fatal(FARGS, "Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
635 st_names(cptp, ecpt), xdr_datatype_names[dtc],
636 xdr_datatype_names[dt]);
638 if (list || !(sflags & (1<<ecpt)))
651 res = xdr_vector(xd, (char *)vp, nf,
652 (unsigned int)sizeof(double), (xdrproc_t)xdr_double);
659 pr_doubles(list, 0, st_names(cptp, ecpt), vp, nf);
669 static int do_cpte_double(XDR *xd, int cptp, int ecpt, int sflags,
670 double *r, FILE *list)
672 return do_cpte_doubles(xd, cptp, ecpt, sflags, 1, &r, list);
676 static int do_cpte_rvecs(XDR *xd, int cptp, int ecpt, int sflags,
677 int n, rvec **v, FILE *list)
681 return do_cpte_reals_low(xd, cptp, ecpt, sflags,
682 n*DIM, NULL, (real **)v, list, ecprRVEC);
685 static int do_cpte_matrix(XDR *xd, int cptp, int ecpt, int sflags,
686 matrix v, FILE *list)
691 vr = (real *)&(v[0][0]);
692 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
693 DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
695 if (list && ret == 0)
697 pr_rvecs(list, 0, st_names(cptp, ecpt), v, DIM);
704 static int do_cpte_nmatrix(XDR *xd, int cptp, int ecpt, int sflags,
705 int n, real **v, FILE *list)
710 char name[CPTSTRLEN];
717 for (i = 0; i < n; i++)
721 reti = do_cpte_reals_low(xd, cptp, ecpt, sflags, n, NULL, &(v[i]), NULL, ecprREAL);
722 if (list && reti == 0)
724 sprintf(name, "%s[%d]", st_names(cptp, ecpt), i);
725 pr_reals(list, 0, name, v[i], n);
735 static int do_cpte_matrices(XDR *xd, int cptp, int ecpt, int sflags,
736 int n, matrix **v, FILE *list)
739 matrix *vp, *va = NULL;
745 res = xdr_int(xd, &nf);
750 if (list == NULL && nf != n)
752 gmx_fatal(FARGS, "Count mismatch for state entry %s, code count is %d, file count is %d\n", st_names(cptp, ecpt), n, nf);
754 if (list || !(sflags & (1<<ecpt)))
767 snew(vr, nf*DIM*DIM);
768 for (i = 0; i < nf; i++)
770 for (j = 0; j < DIM; j++)
772 for (k = 0; k < DIM; k++)
774 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
778 ret = do_cpte_reals_low(xd, cptp, ecpt, sflags,
779 nf*DIM*DIM, NULL, &vr, NULL, ecprMATRIX);
780 for (i = 0; i < nf; i++)
782 for (j = 0; j < DIM; j++)
784 for (k = 0; k < DIM; k++)
786 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
792 if (list && ret == 0)
794 for (i = 0; i < nf; i++)
796 pr_rvecs(list, 0, st_names(cptp, ecpt), vp[i], DIM);
807 static void do_cpt_header(XDR *xd, gmx_bool bRead, int *file_version,
808 char **version, char **btime, char **buser, char **bhost,
810 char **fprog, char **ftime,
811 int *eIntegrator, int *simulation_part,
812 gmx_int64_t *step, double *t,
813 int *nnodes, int *dd_nc, int *npme,
814 int *natoms, int *ngtc, int *nnhpres, int *nhchainlength,
815 int *nlambda, int *flags_state,
816 int *flags_eks, int *flags_enh, int *flags_dfh,
817 int *nED, int *eSwapCoords,
834 res = xdr_int(xd, &magic);
837 gmx_fatal(FARGS, "The checkpoint file is empty/corrupted, or maybe you are out of disk space?");
839 if (magic != CPT_MAGIC1)
841 gmx_fatal(FARGS, "Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
842 "The checkpoint file is corrupted or not a checkpoint file",
848 gmx_gethostname(fhost, 255);
850 do_cpt_string_err(xd, bRead, "GROMACS version", version, list);
851 do_cpt_string_err(xd, bRead, "GROMACS build time", btime, list);
852 do_cpt_string_err(xd, bRead, "GROMACS build user", buser, list);
853 do_cpt_string_err(xd, bRead, "GROMACS build host", bhost, list);
854 do_cpt_string_err(xd, bRead, "generating program", fprog, list);
855 do_cpt_string_err(xd, bRead, "generation time", ftime, list);
856 *file_version = cpt_version;
857 do_cpt_int_err(xd, "checkpoint file version", file_version, list);
858 if (*file_version > cpt_version)
860 gmx_fatal(FARGS, "Attempting to read a checkpoint file of version %d with code of version %d\n", *file_version, cpt_version);
862 if (*file_version >= 13)
864 do_cpt_int_err(xd, "GROMACS double precision", double_prec, list);
870 if (*file_version >= 12)
872 do_cpt_string_err(xd, bRead, "generating host", &fhost, list);
878 do_cpt_int_err(xd, "#atoms", natoms, list);
879 do_cpt_int_err(xd, "#T-coupling groups", ngtc, list);
880 if (*file_version >= 10)
882 do_cpt_int_err(xd, "#Nose-Hoover T-chains", nhchainlength, list);
888 if (*file_version >= 11)
890 do_cpt_int_err(xd, "#Nose-Hoover T-chains for barostat ", nnhpres, list);
896 if (*file_version >= 14)
898 do_cpt_int_err(xd, "# of total lambda states ", nlambda, list);
904 do_cpt_int_err(xd, "integrator", eIntegrator, list);
905 if (*file_version >= 3)
907 do_cpt_int_err(xd, "simulation part #", simulation_part, list);
911 *simulation_part = 1;
913 if (*file_version >= 5)
915 do_cpt_step_err(xd, "step", step, list);
919 do_cpt_int_err(xd, "step", &idum, list);
922 do_cpt_double_err(xd, "t", t, list);
923 do_cpt_int_err(xd, "#PP-nodes", nnodes, list);
925 do_cpt_int_err(xd, "dd_nc[x]", dd_nc ? &(dd_nc[0]) : &idum, list);
926 do_cpt_int_err(xd, "dd_nc[y]", dd_nc ? &(dd_nc[1]) : &idum, list);
927 do_cpt_int_err(xd, "dd_nc[z]", dd_nc ? &(dd_nc[2]) : &idum, list);
928 do_cpt_int_err(xd, "#PME-only nodes", npme, list);
929 do_cpt_int_err(xd, "state flags", flags_state, list);
930 if (*file_version >= 4)
932 do_cpt_int_err(xd, "ekin data flags", flags_eks, list);
933 do_cpt_int_err(xd, "energy history flags", flags_enh, list);
938 *flags_enh = (*flags_state >> (estORIRE_DTAV+1));
939 *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
940 (1<<(estORIRE_DTAV+2)) |
941 (1<<(estORIRE_DTAV+3))));
943 if (*file_version >= 14)
945 do_cpt_int_err(xd, "df history flags", flags_dfh, list);
952 if (*file_version >= 15)
954 do_cpt_int_err(xd, "ED data sets", nED, list);
960 if (*file_version >= 16)
962 do_cpt_int_err(xd, "swap", eSwapCoords, list);
966 static int do_cpt_footer(XDR *xd, int file_version)
971 if (file_version >= 2)
974 res = xdr_int(xd, &magic);
979 if (magic != CPT_MAGIC2)
988 static int do_cpt_state(XDR *xd, gmx_bool bRead,
989 int fflags, t_state *state,
999 nnht = state->nhchainlength*state->ngtc;
1000 nnhtp = state->nhchainlength*state->nnhpres;
1002 if (bRead) /* we need to allocate space for dfhist if we are reading */
1004 init_df_history(&state->dfhist, state->dfhist.nlambda);
1007 sflags = state->flags;
1008 for (i = 0; (i < estNR && ret == 0); i++)
1010 if (fflags & (1<<i))
1014 case estLAMBDA: ret = do_cpte_reals(xd, cptpEST, i, sflags, efptNR, &(state->lambda), list); break;
1015 case estFEPSTATE: ret = do_cpte_int (xd, cptpEST, i, sflags, &state->fep_state, list); break;
1016 case estBOX: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box, list); break;
1017 case estBOX_REL: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->box_rel, list); break;
1018 case estBOXV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->boxv, list); break;
1019 case estPRES_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->pres_prev, list); break;
1020 case estSVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->svir_prev, list); break;
1021 case estFVIR_PREV: ret = do_cpte_matrix(xd, cptpEST, i, sflags, state->fvir_prev, list); break;
1022 case estNH_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_xi, list); break;
1023 case estNH_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnht, &state->nosehoover_vxi, list); break;
1024 case estNHPRES_XI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_xi, list); break;
1025 case estNHPRES_VXI: ret = do_cpte_doubles(xd, cptpEST, i, sflags, nnhtp, &state->nhpres_vxi, list); break;
1026 case estTC_INT: ret = do_cpte_doubles(xd, cptpEST, i, sflags, state->ngtc, &state->therm_integral, list); break;
1027 case estVETA: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->veta, list); break;
1028 case estVOL0: ret = do_cpte_real(xd, cptpEST, i, sflags, &state->vol0, list); break;
1029 case estX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->x, list); break;
1030 case estV: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->v, list); break;
1031 case estSDX: ret = do_cpte_rvecs(xd, cptpEST, i, sflags, state->natoms, &state->sd_X, list); break;
1032 /* The RNG entries are no longer written,
1033 * the next 4 lines are only for reading old files.
1035 case estLD_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1036 case estLD_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1037 case estMC_RNG: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1038 case estMC_RNGI: ret = do_cpte_ints(xd, cptpEST, i, sflags, 0, NULL, list); break;
1039 case estDISRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.disre_initf, list); break;
1040 case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.ndisrepairs, &state->hist.disre_rm3tav, list); break;
1041 case estORIRE_INITF: ret = do_cpte_real (xd, cptpEST, i, sflags, &state->hist.orire_initf, list); break;
1042 case estORIRE_DTAV: ret = do_cpte_n_reals(xd, cptpEST, i, sflags, &state->hist.norire_Dtav, &state->hist.orire_Dtav, list); break;
1044 gmx_fatal(FARGS, "Unknown state entry %d\n"
1045 "You are probably reading a new checkpoint file with old code", i);
1053 static int do_cpt_ekinstate(XDR *xd, int fflags, ekinstate_t *ekins,
1061 for (i = 0; (i < eeksNR && ret == 0); i++)
1063 if (fflags & (1<<i))
1068 case eeksEKIN_N: ret = do_cpte_int(xd, cptpEEKS, i, fflags, &ekins->ekin_n, list); break;
1069 case eeksEKINH: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh, list); break;
1070 case eeksEKINF: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinf, list); break;
1071 case eeksEKINO: ret = do_cpte_matrices(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinh_old, list); break;
1072 case eeksEKINTOTAL: ret = do_cpte_matrix(xd, cptpEEKS, i, fflags, ekins->ekin_total, list); break;
1073 case eeksEKINSCALEF: ret = do_cpte_doubles(xd, cptpEEKS, i, fflags, ekins->ekin_n, &ekins->ekinscalef_nhc, list); break;
1074 case eeksVSCALE: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->vscale_nhc, list); break;
1075 case eeksEKINSCALEH: ret = do_cpte_doubles(xd, 1, cptpEEKS, fflags, ekins->ekin_n, &ekins->ekinscaleh_nhc, list); break;
1076 case eeksDEKINDL: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->dekindl, list); break;
1077 case eeksMVCOS: ret = do_cpte_real(xd, 1, cptpEEKS, fflags, &ekins->mvcos, list); break;
1079 gmx_fatal(FARGS, "Unknown ekin data state entry %d\n"
1080 "You are probably reading a new checkpoint file with old code", i);
1089 static int do_cpt_swapstate(XDR *xd, gmx_bool bRead, swapstate_t *swapstate, FILE *list)
1093 int swap_cpt_version = 1;
1096 if (eswapNO == swapstate->eSwapCoords)
1101 swapstate->bFromCpt = bRead;
1103 do_cpt_int_err(xd, "swap checkpoint version", &swap_cpt_version, list);
1104 do_cpt_int_err(xd, "swap coupling steps", &swapstate->nAverage, list);
1106 /* When reading, init_swapcoords has not been called yet,
1107 * so we have to allocate memory first. */
1109 for (ic = 0; ic < eCompNR; ic++)
1111 for (ii = 0; ii < eIonNR; ii++)
1115 do_cpt_int_err(xd, "swap requested atoms", &swapstate->nat_req[ic][ii], list);
1119 do_cpt_int_err(xd, "swap requested atoms p", swapstate->nat_req_p[ic][ii], list);
1124 do_cpt_int_err(xd, "swap influx netto", &swapstate->inflow_netto[ic][ii], list);
1128 do_cpt_int_err(xd, "swap influx netto p", swapstate->inflow_netto_p[ic][ii], list);
1131 if (bRead && (NULL == swapstate->nat_past[ic][ii]) )
1133 snew(swapstate->nat_past[ic][ii], swapstate->nAverage);
1136 for (j = 0; j < swapstate->nAverage; j++)
1140 do_cpt_int_err(xd, "swap past atom counts", &swapstate->nat_past[ic][ii][j], list);
1144 do_cpt_int_err(xd, "swap past atom counts p", &swapstate->nat_past_p[ic][ii][j], list);
1150 /* Ion flux per channel */
1151 for (ic = 0; ic < eChanNR; ic++)
1153 for (ii = 0; ii < eIonNR; ii++)
1157 do_cpt_int_err(xd, "channel flux", &swapstate->fluxfromAtoB[ic][ii], list);
1161 do_cpt_int_err(xd, "channel flux p", swapstate->fluxfromAtoB_p[ic][ii], list);
1166 /* Ion flux leakage */
1169 snew(swapstate->fluxleak, 1);
1171 do_cpt_int_err(xd, "flux leakage", swapstate->fluxleak, list);
1174 do_cpt_int_err(xd, "number of ions", &swapstate->nions, list);
1178 snew(swapstate->channel_label, swapstate->nions);
1179 snew(swapstate->comp_from, swapstate->nions);
1182 do_cpt_u_chars(xd, "channel history", swapstate->nions, swapstate->channel_label, list);
1183 do_cpt_u_chars(xd, "domain history", swapstate->nions, swapstate->comp_from, list);
1185 /* Save the last known whole positions to checkpoint
1186 * file to be able to also make multimeric channels whole in PBC */
1187 do_cpt_int_err(xd, "Ch0 atoms", &swapstate->nat[eChan0], list);
1188 do_cpt_int_err(xd, "Ch1 atoms", &swapstate->nat[eChan1], list);
1191 snew(swapstate->xc_old_whole[eChan0], swapstate->nat[eChan0]);
1192 snew(swapstate->xc_old_whole[eChan1], swapstate->nat[eChan1]);
1193 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], swapstate->xc_old_whole[eChan0], list);
1194 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], swapstate->xc_old_whole[eChan1], list);
1198 do_cpt_n_rvecs_err(xd, "Ch0 whole x", swapstate->nat[eChan0], *swapstate->xc_old_whole_p[eChan0], list);
1199 do_cpt_n_rvecs_err(xd, "Ch1 whole x", swapstate->nat[eChan1], *swapstate->xc_old_whole_p[eChan1], list);
1206 static int do_cpt_enerhist(XDR *xd, gmx_bool bRead,
1207 int fflags, energyhistory_t *enerhist,
1218 enerhist->nsteps = 0;
1220 enerhist->nsteps_sim = 0;
1221 enerhist->nsum_sim = 0;
1222 enerhist->dht = NULL;
1224 if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
1226 snew(enerhist->dht, 1);
1227 enerhist->dht->ndh = NULL;
1228 enerhist->dht->dh = NULL;
1229 enerhist->dht->start_lambda_set = FALSE;
1233 for (i = 0; (i < eenhNR && ret == 0); i++)
1235 if (fflags & (1<<i))
1239 case eenhENERGY_N: ret = do_cpte_int(xd, cptpEENH, i, fflags, &enerhist->nener, list); break;
1240 case eenhENERGY_AVER: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_ave, list); break;
1241 case eenhENERGY_SUM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum, list); break;
1242 case eenhENERGY_NSUM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum, list); break;
1243 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd, cptpEENH, i, fflags, enerhist->nener, &enerhist->ener_sum_sim, list); break;
1244 case eenhENERGY_NSUM_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsum_sim, list); break;
1245 case eenhENERGY_NSTEPS: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps, list); break;
1246 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd, eenh_names[i], &enerhist->nsteps_sim, list); break;
1247 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd, eenh_names[i], &(enerhist->dht->nndh), list);
1248 if (bRead) /* now allocate memory for it */
1250 snew(enerhist->dht->dh, enerhist->dht->nndh);
1251 snew(enerhist->dht->ndh, enerhist->dht->nndh);
1252 for (j = 0; j < enerhist->dht->nndh; j++)
1254 enerhist->dht->ndh[j] = 0;
1255 enerhist->dht->dh[j] = NULL;
1259 case eenhENERGY_DELTA_H_LIST:
1260 for (j = 0; j < enerhist->dht->nndh; j++)
1262 ret = do_cpte_n_reals(xd, cptpEENH, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list);
1265 case eenhENERGY_DELTA_H_STARTTIME:
1266 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_time), list); break;
1267 case eenhENERGY_DELTA_H_STARTLAMBDA:
1268 ret = do_cpte_double(xd, cptpEENH, i, fflags, &(enerhist->dht->start_lambda), list); break;
1270 gmx_fatal(FARGS, "Unknown energy history entry %d\n"
1271 "You are probably reading a new checkpoint file with old code", i);
1276 if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
1278 /* Assume we have an old file format and copy sum to sum_sim */
1279 srenew(enerhist->ener_sum_sim, enerhist->nener);
1280 for (i = 0; i < enerhist->nener; i++)
1282 enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1284 fflags |= (1<<eenhENERGY_SUM_SIM);
1287 if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1288 !(fflags & (1<<eenhENERGY_NSTEPS)))
1290 /* Assume we have an old file format and copy nsum to nsteps */
1291 enerhist->nsteps = enerhist->nsum;
1292 fflags |= (1<<eenhENERGY_NSTEPS);
1294 if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1295 !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1297 /* Assume we have an old file format and copy nsum to nsteps */
1298 enerhist->nsteps_sim = enerhist->nsum_sim;
1299 fflags |= (1<<eenhENERGY_NSTEPS_SIM);
1305 static int do_cpt_df_hist(XDR *xd, int fflags, df_history_t *dfhist, FILE *list)
1310 nlambda = dfhist->nlambda;
1313 for (i = 0; (i < edfhNR && ret == 0); i++)
1315 if (fflags & (1<<i))
1319 case edfhBEQUIL: ret = do_cpte_int(xd, cptpEDFH, i, fflags, &dfhist->bEquil, list); break;
1320 case edfhNATLAMBDA: ret = do_cpte_ints(xd, cptpEDFH, i, fflags, nlambda, &dfhist->n_at_lam, list); break;
1321 case edfhWLHISTO: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->wl_histo, list); break;
1322 case edfhWLDELTA: ret = do_cpte_real(xd, cptpEDFH, i, fflags, &dfhist->wl_delta, list); break;
1323 case edfhSUMWEIGHTS: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_weights, list); break;
1324 case edfhSUMDG: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_dg, list); break;
1325 case edfhSUMMINVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_minvar, list); break;
1326 case edfhSUMVAR: ret = do_cpte_reals(xd, cptpEDFH, i, fflags, nlambda, &dfhist->sum_variance, list); break;
1327 case edfhACCUMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p, list); break;
1328 case edfhACCUMM: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m, list); break;
1329 case edfhACCUMP2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_p2, list); break;
1330 case edfhACCUMM2: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->accum_m2, list); break;
1331 case edfhTIJ: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij, list); break;
1332 case edfhTIJEMP: ret = do_cpte_nmatrix(xd, cptpEDFH, i, fflags, nlambda, dfhist->Tij_empirical, list); break;
1335 gmx_fatal(FARGS, "Unknown df history entry %d\n"
1336 "You are probably reading a new checkpoint file with old code", i);
1345 /* This function stores the last whole configuration of the reference and
1346 * average structure in the .cpt file
1348 static int do_cpt_EDstate(XDR *xd, gmx_bool bRead,
1349 edsamstate_t *EDstate, FILE *list)
1356 EDstate->bFromCpt = bRead;
1358 if (EDstate->nED <= 0)
1363 /* When reading, init_edsam has not been called yet,
1364 * so we have to allocate memory first. */
1367 snew(EDstate->nref, EDstate->nED);
1368 snew(EDstate->old_sref, EDstate->nED);
1369 snew(EDstate->nav, EDstate->nED);
1370 snew(EDstate->old_sav, EDstate->nED);
1373 /* Read/write the last whole conformation of SREF and SAV for each ED dataset (usually only one) */
1374 for (i = 0; i < EDstate->nED; i++)
1376 /* Reference structure SREF */
1377 sprintf(buf, "ED%d # of atoms in reference structure", i+1);
1378 do_cpt_int_err(xd, buf, &EDstate->nref[i], list);
1379 sprintf(buf, "ED%d x_ref", i+1);
1382 snew(EDstate->old_sref[i], EDstate->nref[i]);
1383 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref[i], list);
1387 do_cpt_n_rvecs_err(xd, buf, EDstate->nref[i], EDstate->old_sref_p[i], list);
1390 /* Average structure SAV */
1391 sprintf(buf, "ED%d # of atoms in average structure", i+1);
1392 do_cpt_int_err(xd, buf, &EDstate->nav[i], list);
1393 sprintf(buf, "ED%d x_av", i+1);
1396 snew(EDstate->old_sav[i], EDstate->nav[i]);
1397 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav[i], list);
1401 do_cpt_n_rvecs_err(xd, buf, EDstate->nav[i], EDstate->old_sav_p[i], list);
1409 static int do_cpt_files(XDR *xd, gmx_bool bRead,
1410 gmx_file_position_t **p_outputfiles, int *nfiles,
1411 FILE *list, int file_version)
1415 gmx_off_t mask = 0xFFFFFFFFL;
1416 int offset_high, offset_low;
1418 gmx_file_position_t *outputfiles;
1420 if (do_cpt_int(xd, "number of output files", nfiles, list) != 0)
1427 snew(*p_outputfiles, *nfiles);
1430 outputfiles = *p_outputfiles;
1432 for (i = 0; i < *nfiles; i++)
1434 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1437 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1438 strncpy(outputfiles[i].filename, buf, CPTSTRLEN-1);
1444 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1448 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1452 outputfiles[i].offset = ( ((gmx_off_t) offset_high) << 32 ) | ( (gmx_off_t) offset_low & mask );
1456 buf = outputfiles[i].filename;
1457 do_cpt_string_err(xd, bRead, "output filename", &buf, list);
1459 offset = outputfiles[i].offset;
1467 offset_low = (int) (offset & mask);
1468 offset_high = (int) ((offset >> 32) & mask);
1470 if (do_cpt_int(xd, "file_offset_high", &offset_high, list) != 0)
1474 if (do_cpt_int(xd, "file_offset_low", &offset_low, list) != 0)
1479 if (file_version >= 8)
1481 if (do_cpt_int(xd, "file_checksum_size", &(outputfiles[i].chksum_size),
1486 if (do_cpt_u_chars(xd, "file_checksum", 16, outputfiles[i].chksum, list) != 0)
1493 outputfiles[i].chksum_size = -1;
1500 void write_checkpoint(const char *fn, gmx_bool bNumberAndKeep,
1501 FILE *fplog, t_commrec *cr,
1502 int eIntegrator, int simulation_part,
1503 gmx_bool bExpanded, int elamstats,
1504 gmx_int64_t step, double t, t_state *state)
1514 char *fntemp; /* the temporary checkpoint file name */
1516 char timebuf[STRLEN];
1517 int nppnodes, npmenodes, flag_64bit;
1518 char buf[1024], suffix[5+STEPSTRSIZE], sbuf[STEPSTRSIZE];
1519 gmx_file_position_t *outputfiles;
1522 int flags_eks, flags_enh, flags_dfh, i;
1525 if (DOMAINDECOMP(cr))
1527 nppnodes = cr->dd->nnodes;
1528 npmenodes = cr->npmenodes;
1536 #ifndef GMX_NO_RENAME
1537 /* make the new temporary filename */
1538 snew(fntemp, strlen(fn)+5+STEPSTRSIZE);
1540 fntemp[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1541 sprintf(suffix, "_%s%s", "step", gmx_step_str(step, sbuf));
1542 strcat(fntemp, suffix);
1543 strcat(fntemp, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1545 /* if we can't rename, we just overwrite the cpt file.
1546 * dangerous if interrupted.
1548 snew(fntemp, strlen(fn));
1552 gmx_ctime_r(&now, timebuf, STRLEN);
1556 fprintf(fplog, "Writing checkpoint, step %s at %s\n\n",
1557 gmx_step_str(step, buf), timebuf);
1560 /* Get offsets for open files */
1561 gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1563 fp = gmx_fio_open(fntemp, "w");
1565 if (state->ekinstate.bUpToDate)
1568 ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) |
1569 (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) |
1570 (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1578 if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1580 flags_enh |= (1<<eenhENERGY_N);
1581 if (state->enerhist.nsum > 0)
1583 flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1584 (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSUM));
1586 if (state->enerhist.nsum_sim > 0)
1588 flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSTEPS_SIM) |
1589 (1<<eenhENERGY_NSUM_SIM));
1591 if (state->enerhist.dht)
1593 flags_enh |= ( (1<< eenhENERGY_DELTA_H_NN) |
1594 (1<< eenhENERGY_DELTA_H_LIST) |
1595 (1<< eenhENERGY_DELTA_H_STARTTIME) |
1596 (1<< eenhENERGY_DELTA_H_STARTLAMBDA) );
1602 flags_dfh = ((1<<edfhBEQUIL) | (1<<edfhNATLAMBDA) | (1<<edfhSUMWEIGHTS) | (1<<edfhSUMDG) |
1603 (1<<edfhTIJ) | (1<<edfhTIJEMP));
1606 flags_dfh |= ((1<<edfhWLDELTA) | (1<<edfhWLHISTO));
1608 if ((elamstats == elamstatsMINVAR) || (elamstats == elamstatsBARKER) || (elamstats == elamstatsMETROPOLIS))
1610 flags_dfh |= ((1<<edfhACCUMP) | (1<<edfhACCUMM) | (1<<edfhACCUMP2) | (1<<edfhACCUMM2)
1611 | (1<<edfhSUMMINVAR) | (1<<edfhSUMVAR));
1619 /* We can check many more things now (CPU, acceleration, etc), but
1620 * it is highly unlikely to have two separate builds with exactly
1621 * the same version, user, time, and build host!
1624 version = gmx_strdup(gmx_version());
1625 btime = gmx_strdup(BUILD_TIME);
1626 buser = gmx_strdup(BUILD_USER);
1627 bhost = gmx_strdup(BUILD_HOST);
1629 double_prec = GMX_CPT_BUILD_DP;
1630 fprog = gmx_strdup(Program());
1632 ftime = &(timebuf[0]);
1634 do_cpt_header(gmx_fio_getxdr(fp), FALSE, &file_version,
1635 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1636 &eIntegrator, &simulation_part, &step, &t, &nppnodes,
1637 DOMAINDECOMP(cr) ? cr->dd->nc : NULL, &npmenodes,
1638 &state->natoms, &state->ngtc, &state->nnhpres,
1639 &state->nhchainlength, &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
1640 &state->edsamstate.nED, &state->swapstate.eSwapCoords,
1649 if ((do_cpt_state(gmx_fio_getxdr(fp), FALSE, state->flags, state, NULL) < 0) ||
1650 (do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL) < 0) ||
1651 (do_cpt_enerhist(gmx_fio_getxdr(fp), FALSE, flags_enh, &state->enerhist, NULL) < 0) ||
1652 (do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL) < 0) ||
1653 (do_cpt_EDstate(gmx_fio_getxdr(fp), FALSE, &state->edsamstate, NULL) < 0) ||
1654 (do_cpt_swapstate(gmx_fio_getxdr(fp), FALSE, &state->swapstate, NULL) < 0) ||
1655 (do_cpt_files(gmx_fio_getxdr(fp), FALSE, &outputfiles, &noutputfiles, NULL,
1658 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1661 do_cpt_footer(gmx_fio_getxdr(fp), file_version);
1663 /* we really, REALLY, want to make sure to physically write the checkpoint,
1664 and all the files it depends on, out to disk. Because we've
1665 opened the checkpoint with gmx_fio_open(), it's in our list
1667 ret = gmx_fio_all_output_fsync();
1673 "Cannot fsync '%s'; maybe you are out of disk space?",
1674 gmx_fio_getname(ret));
1676 if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV) == NULL)
1686 if (gmx_fio_close(fp) != 0)
1688 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1691 /* we don't move the checkpoint if the user specified they didn't want it,
1692 or if the fsyncs failed */
1693 #ifndef GMX_NO_RENAME
1694 if (!bNumberAndKeep && !ret)
1698 /* Rename the previous checkpoint file */
1700 buf[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1701 strcat(buf, "_prev");
1702 strcat(buf, fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1704 /* we copy here so that if something goes wrong between now and
1705 * the rename below, there's always a state.cpt.
1706 * If renames are atomic (such as in POSIX systems),
1707 * this copying should be unneccesary.
1709 gmx_file_copy(fn, buf, FALSE);
1710 /* We don't really care if this fails:
1711 * there's already a new checkpoint.
1714 gmx_file_rename(fn, buf);
1717 if (gmx_file_rename(fntemp, fn) != 0)
1719 gmx_file("Cannot rename checkpoint file; maybe you are out of disk space?");
1722 #endif /* GMX_NO_RENAME */
1728 /*code for alternate checkpointing scheme. moved from top of loop over
1730 fcRequestCheckPoint();
1731 if (fcCheckPointParallel( cr->nodeid, NULL, 0) == 0)
1733 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", step );
1735 #endif /* end GMX_FAHCORE block */
1738 static void print_flag_mismatch(FILE *fplog, int sflags, int fflags)
1742 fprintf(fplog, "\nState entry mismatch between the simulation and the checkpoint file\n");
1743 fprintf(fplog, "Entries which are not present in the checkpoint file will not be updated\n");
1744 fprintf(fplog, " %24s %11s %11s\n", "", "simulation", "checkpoint");
1745 for (i = 0; i < estNR; i++)
1747 if ((sflags & (1<<i)) || (fflags & (1<<i)))
1749 fprintf(fplog, " %24s %11s %11s\n",
1751 (sflags & (1<<i)) ? " present " : "not present",
1752 (fflags & (1<<i)) ? " present " : "not present");
1757 static void check_int(FILE *fplog, const char *type, int p, int f, gmx_bool *mm)
1759 FILE *fp = fplog ? fplog : stderr;
1763 fprintf(fp, " %s mismatch,\n", type);
1764 fprintf(fp, " current program: %d\n", p);
1765 fprintf(fp, " checkpoint file: %d\n", f);
1771 static void check_string(FILE *fplog, const char *type, const char *p,
1772 const char *f, gmx_bool *mm)
1774 FILE *fp = fplog ? fplog : stderr;
1776 if (strcmp(p, f) != 0)
1778 fprintf(fp, " %s mismatch,\n", type);
1779 fprintf(fp, " current program: %s\n", p);
1780 fprintf(fp, " checkpoint file: %s\n", f);
1786 static void check_match(FILE *fplog,
1788 char *btime, char *buser, char *bhost, int double_prec,
1790 t_commrec *cr, int npp_f, int npme_f,
1791 ivec dd_nc, ivec dd_nc_f)
1798 check_string(fplog, "Version", gmx_version(), version, &mm);
1799 check_string(fplog, "Build time", BUILD_TIME, btime, &mm);
1800 check_string(fplog, "Build user", BUILD_USER, buser, &mm);
1801 check_string(fplog, "Build host", BUILD_HOST, bhost, &mm);
1802 check_int (fplog, "Double prec.", GMX_CPT_BUILD_DP, double_prec, &mm);
1803 check_string(fplog, "Program name", Program(), fprog, &mm);
1805 check_int (fplog, "#nodes", cr->nnodes, npp_f+npme_f, &mm);
1808 check_int (fplog, "#PME-nodes", cr->npmenodes, npme_f, &mm);
1811 if (cr->npmenodes >= 0)
1813 npp -= cr->npmenodes;
1817 check_int (fplog, "#DD-cells[x]", dd_nc[XX], dd_nc_f[XX], &mm);
1818 check_int (fplog, "#DD-cells[y]", dd_nc[YY], dd_nc_f[YY], &mm);
1819 check_int (fplog, "#DD-cells[z]", dd_nc[ZZ], dd_nc_f[ZZ], &mm);
1826 "Gromacs binary or parallel settings not identical to previous run.\n"
1827 "Continuation is exact, but is not guaranteed to be binary identical%s.\n\n",
1828 fplog ? ",\n see the log file for details" : "");
1833 "Gromacs binary or parallel settings not identical to previous run.\n"
1834 "Continuation is exact, but is not guaranteed to be binary identical.\n\n");
1839 static void read_checkpoint(const char *fn, FILE **pfplog,
1840 t_commrec *cr, ivec dd_nc,
1841 int eIntegrator, int *init_fep_state, gmx_int64_t *step, double *t,
1842 t_state *state, gmx_bool *bReadEkin,
1843 int *simulation_part,
1844 gmx_bool bAppendOutputFiles, gmx_bool bForceAppend)
1849 char *version, *btime, *buser, *bhost, *fprog, *ftime;
1851 char filename[STRLEN], buf[STEPSTRSIZE];
1852 int nppnodes, eIntegrator_f, nppnodes_f, npmenodes_f;
1854 int natoms, ngtc, nnhpres, nhchainlength, nlambda, fflags, flags_eks, flags_enh, flags_dfh;
1857 gmx_file_position_t *outputfiles;
1859 t_fileio *chksum_file;
1860 FILE * fplog = *pfplog;
1861 unsigned char digest[16];
1862 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1863 struct flock fl; /* don't initialize here: the struct order is OS
1867 const char *int_warn =
1868 "WARNING: The checkpoint file was generated with integrator %s,\n"
1869 " while the simulation uses integrator %s\n\n";
1871 #if !defined __native_client__ && !defined GMX_NATIVE_WINDOWS
1872 fl.l_type = F_WRLCK;
1873 fl.l_whence = SEEK_SET;
1879 fp = gmx_fio_open(fn, "r");
1880 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
1881 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
1882 &eIntegrator_f, simulation_part, step, t,
1883 &nppnodes_f, dd_nc_f, &npmenodes_f,
1884 &natoms, &ngtc, &nnhpres, &nhchainlength, &nlambda,
1885 &fflags, &flags_eks, &flags_enh, &flags_dfh,
1886 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
1888 if (bAppendOutputFiles &&
1889 file_version >= 13 && double_prec != GMX_CPT_BUILD_DP)
1891 gmx_fatal(FARGS, "Output file appending requested, but the code and checkpoint file precision (single/double) don't match");
1894 if (cr == NULL || MASTER(cr))
1896 fprintf(stderr, "\nReading checkpoint file %s generated: %s\n\n",
1900 /* This will not be written if we do appending, since fplog is still NULL then */
1903 fprintf(fplog, "\n");
1904 fprintf(fplog, "Reading checkpoint file %s\n", fn);
1905 fprintf(fplog, " file generated by: %s\n", fprog);
1906 fprintf(fplog, " file generated at: %s\n", ftime);
1907 fprintf(fplog, " GROMACS build time: %s\n", btime);
1908 fprintf(fplog, " GROMACS build user: %s\n", buser);
1909 fprintf(fplog, " GROMACS build host: %s\n", bhost);
1910 fprintf(fplog, " GROMACS double prec.: %d\n", double_prec);
1911 fprintf(fplog, " simulation part #: %d\n", *simulation_part);
1912 fprintf(fplog, " step: %s\n", gmx_step_str(*step, buf));
1913 fprintf(fplog, " time: %f\n", *t);
1914 fprintf(fplog, "\n");
1917 if (natoms != state->natoms)
1919 gmx_fatal(FARGS, "Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms", natoms, state->natoms);
1921 if (ngtc != state->ngtc)
1923 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);
1925 if (nnhpres != state->nnhpres)
1927 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);
1930 if (nlambda != state->dfhist.nlambda)
1932 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);
1935 init_gtc_state(state, state->ngtc, state->nnhpres, nhchainlength); /* need to keep this here to keep the tpr format working */
1936 /* write over whatever was read; we use the number of Nose-Hoover chains from the checkpoint */
1938 if (eIntegrator_f != eIntegrator)
1942 fprintf(stderr, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1944 if (bAppendOutputFiles)
1947 "Output file appending requested, but input/checkpoint integrators do not match.\n"
1948 "Stopping the run to prevent you from ruining all your data...\n"
1949 "If you _really_ know what you are doing, try with the -noappend option.\n");
1953 fprintf(fplog, int_warn, EI(eIntegrator_f), EI(eIntegrator));
1962 else if (cr->nnodes == nppnodes_f + npmenodes_f)
1964 if (cr->npmenodes < 0)
1966 cr->npmenodes = npmenodes_f;
1968 nppnodes = cr->nnodes - cr->npmenodes;
1969 if (nppnodes == nppnodes_f)
1971 for (d = 0; d < DIM; d++)
1975 dd_nc[d] = dd_nc_f[d];
1982 /* The number of PP nodes has not been set yet */
1986 if (fflags != state->flags)
1991 if (bAppendOutputFiles)
1994 "Output file appending requested, but input and checkpoint states are not identical.\n"
1995 "Stopping the run to prevent you from ruining all your data...\n"
1996 "You can try with the -noappend option, and get more info in the log file.\n");
1999 if (getenv("GMX_ALLOW_CPT_MISMATCH") == NULL)
2001 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");
2006 "WARNING: The checkpoint state entries do not match the simulation,\n"
2007 " see the log file for details\n\n");
2013 print_flag_mismatch(fplog, state->flags, fflags);
2020 check_match(fplog, version, btime, buser, bhost, double_prec, fprog,
2021 cr, nppnodes_f, npmenodes_f, dd_nc, dd_nc_f);
2024 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, fflags, state, NULL);
2025 *init_fep_state = state->fep_state; /* there should be a better way to do this than setting it here.
2026 Investigate for 5.0. */
2031 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2036 *bReadEkin = ((flags_eks & (1<<eeksEKINH)) || (flags_eks & (1<<eeksEKINF)) || (flags_eks & (1<<eeksEKINO)) ||
2037 ((flags_eks & (1<<eeksEKINSCALEF)) | (flags_eks & (1<<eeksEKINSCALEH)) | (flags_eks & (1<<eeksVSCALE))));
2039 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2040 flags_enh, &state->enerhist, NULL);
2046 if (file_version < 6)
2048 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.";
2050 fprintf(stderr, "\nWARNING: %s\n\n", warn);
2053 fprintf(fplog, "\nWARNING: %s\n\n", warn);
2055 state->enerhist.nsum = *step;
2056 state->enerhist.nsum_sim = *step;
2059 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2065 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2071 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2077 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, NULL, file_version);
2083 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2088 if (gmx_fio_close(fp) != 0)
2090 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2099 /* If the user wants to append to output files,
2100 * we use the file pointer positions of the output files stored
2101 * in the checkpoint file and truncate the files such that any frames
2102 * written after the checkpoint time are removed.
2103 * All files are md5sum checked such that we can be sure that
2104 * we do not truncate other (maybe imprortant) files.
2106 if (bAppendOutputFiles)
2108 if (fn2ftp(outputfiles[0].filename) != efLOG)
2110 /* make sure first file is log file so that it is OK to use it for
2113 gmx_fatal(FARGS, "The first output file should always be the log "
2114 "file but instead is: %s. Cannot do appending because of this condition.", outputfiles[0].filename);
2116 for (i = 0; i < nfiles; i++)
2118 if (outputfiles[i].offset < 0)
2120 gmx_fatal(FARGS, "The original run wrote a file called '%s' which "
2121 "is larger than 2 GB, but mdrun did not support large file"
2122 " offsets. Can not append. Run mdrun with -noappend",
2123 outputfiles[i].filename);
2126 chksum_file = gmx_fio_open(outputfiles[i].filename, "a");
2129 chksum_file = gmx_fio_open(outputfiles[i].filename, "r+");
2134 /* Note that there are systems where the lock operation
2135 * will succeed, but a second process can also lock the file.
2136 * We should probably try to detect this.
2138 #if defined __native_client__
2142 #elif defined GMX_NATIVE_WINDOWS
2143 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX) == -1)
2145 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl) == -1)
2148 if (errno == ENOSYS)
2152 gmx_fatal(FARGS, "File locking is not supported on this system. Use -noappend or specify -append explicitly to append anyhow.");
2156 fprintf(stderr, "\nNOTE: File locking is not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2159 fprintf(fplog, "\nNOTE: File locking not supported on this system, will not lock %s\n\n", outputfiles[i].filename);
2163 else if (errno == EACCES || errno == EAGAIN)
2165 gmx_fatal(FARGS, "Failed to lock: %s. Already running "
2166 "simulation?", outputfiles[i].filename);
2170 gmx_fatal(FARGS, "Failed to lock: %s. %s.",
2171 outputfiles[i].filename, strerror(errno));
2176 /* compute md5 chksum */
2177 if (outputfiles[i].chksum_size != -1)
2179 if (gmx_fio_get_file_md5(chksum_file, outputfiles[i].offset,
2180 digest) != outputfiles[i].chksum_size) /*at the end of the call the file position is at the end of the file*/
2182 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.",
2183 outputfiles[i].chksum_size,
2184 outputfiles[i].filename);
2187 if (i == 0) /*log file needs to be seeked in case we need to truncate (other files are truncated below)*/
2189 if (gmx_fio_seek(chksum_file, outputfiles[i].offset))
2191 gmx_fatal(FARGS, "Seek error! Failed to truncate log-file: %s.", strerror(errno));
2196 if (i == 0) /*open log file here - so that lock is never lifted
2197 after chksum is calculated */
2199 *pfplog = gmx_fio_getfp(chksum_file);
2203 gmx_fio_close(chksum_file);
2206 /* compare md5 chksum */
2207 if (outputfiles[i].chksum_size != -1 &&
2208 memcmp(digest, outputfiles[i].chksum, 16) != 0)
2212 fprintf(debug, "chksum for %s: ", outputfiles[i].filename);
2213 for (j = 0; j < 16; j++)
2215 fprintf(debug, "%02x", digest[j]);
2217 fprintf(debug, "\n");
2219 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.",
2220 outputfiles[i].filename);
2225 if (i != 0) /*log file is already seeked to correct position */
2227 #ifdef GMX_NATIVE_WINDOWS
2228 rc = gmx_wintruncate(outputfiles[i].filename, outputfiles[i].offset);
2230 rc = truncate(outputfiles[i].filename, outputfiles[i].offset);
2234 gmx_fatal(FARGS, "Truncation of file %s failed. Cannot do appending because of this failure.", outputfiles[i].filename);
2244 void load_checkpoint(const char *fn, FILE **fplog,
2245 t_commrec *cr, ivec dd_nc,
2246 t_inputrec *ir, t_state *state,
2247 gmx_bool *bReadEkin,
2248 gmx_bool bAppend, gmx_bool bForceAppend)
2255 /* Read the state from the checkpoint file */
2256 read_checkpoint(fn, fplog,
2258 ir->eI, &(ir->fepvals->init_fep_state), &step, &t, state, bReadEkin,
2259 &ir->simulation_part, bAppend, bForceAppend);
2263 gmx_bcast(sizeof(cr->npmenodes), &cr->npmenodes, cr);
2264 gmx_bcast(DIM*sizeof(dd_nc[0]), dd_nc, cr);
2265 gmx_bcast(sizeof(step), &step, cr);
2266 gmx_bcast(sizeof(*bReadEkin), bReadEkin, cr);
2268 ir->bContinuation = TRUE;
2269 if (ir->nsteps >= 0)
2271 ir->nsteps += ir->init_step - step;
2273 ir->init_step = step;
2274 ir->simulation_part += 1;
2277 static void read_checkpoint_data(t_fileio *fp, int *simulation_part,
2278 gmx_int64_t *step, double *t, t_state *state,
2279 int *nfiles, gmx_file_position_t **outputfiles)
2282 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2287 int flags_eks, flags_enh, flags_dfh;
2289 gmx_file_position_t *files_loc = NULL;
2292 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2293 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2294 &eIntegrator, simulation_part, step, t, &nppnodes, dd_nc, &npme,
2295 &state->natoms, &state->ngtc, &state->nnhpres, &state->nhchainlength,
2296 &(state->dfhist.nlambda), &state->flags, &flags_eks, &flags_enh, &flags_dfh,
2297 &state->edsamstate.nED, &state->swapstate.eSwapCoords, NULL);
2299 do_cpt_state(gmx_fio_getxdr(fp), TRUE, state->flags, state, NULL);
2304 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state->ekinstate, NULL);
2309 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2310 flags_enh, &state->enerhist, NULL);
2315 ret = do_cpt_df_hist(gmx_fio_getxdr(fp), flags_dfh, &state->dfhist, NULL);
2321 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state->edsamstate, NULL);
2327 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state->swapstate, NULL);
2333 ret = do_cpt_files(gmx_fio_getxdr(fp), TRUE,
2334 outputfiles != NULL ? outputfiles : &files_loc,
2335 outputfiles != NULL ? nfiles : &nfiles_loc,
2336 NULL, file_version);
2337 if (files_loc != NULL)
2347 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2361 read_checkpoint_state(const char *fn, int *simulation_part,
2362 gmx_int64_t *step, double *t, t_state *state)
2366 fp = gmx_fio_open(fn, "r");
2367 read_checkpoint_data(fp, simulation_part, step, t, state, NULL, NULL);
2368 if (gmx_fio_close(fp) != 0)
2370 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2374 void read_checkpoint_trxframe(t_fileio *fp, t_trxframe *fr)
2376 /* This next line is nasty because the sub-structures of t_state
2377 * cannot be assumed to be zeroed (or even initialized in ways the
2378 * rest of the code might assume). Using snew would be better, but
2379 * this will all go away for 5.0. */
2381 int simulation_part;
2385 init_state(&state, 0, 0, 0, 0, 0);
2387 read_checkpoint_data(fp, &simulation_part, &step, &t, &state, NULL, NULL);
2389 fr->natoms = state.natoms;
2392 fr->step = gmx_int64_to_int(step,
2393 "conversion of checkpoint to trajectory");
2397 fr->lambda = state.lambda[efptFEP];
2398 fr->fep_state = state.fep_state;
2400 fr->bX = (state.flags & (1<<estX));
2406 fr->bV = (state.flags & (1<<estV));
2413 fr->bBox = (state.flags & (1<<estBOX));
2416 copy_mat(state.box, fr->box);
2421 void list_checkpoint(const char *fn, FILE *out)
2425 char *version, *btime, *buser, *bhost, *fprog, *ftime;
2427 int eIntegrator, simulation_part, nppnodes, npme;
2432 int flags_eks, flags_enh, flags_dfh;
2436 gmx_file_position_t *outputfiles;
2439 init_state(&state, -1, -1, -1, -1, 0);
2441 fp = gmx_fio_open(fn, "r");
2442 do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
2443 &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
2444 &eIntegrator, &simulation_part, &step, &t, &nppnodes, dd_nc, &npme,
2445 &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
2446 &(state.dfhist.nlambda), &state.flags,
2447 &flags_eks, &flags_enh, &flags_dfh, &state.edsamstate.nED,
2448 &state.swapstate.eSwapCoords, out);
2449 ret = do_cpt_state(gmx_fio_getxdr(fp), TRUE, state.flags, &state, out);
2454 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp), flags_eks, &state.ekinstate, out);
2459 ret = do_cpt_enerhist(gmx_fio_getxdr(fp), TRUE,
2460 flags_enh, &state.enerhist, out);
2464 ret = do_cpt_df_hist(gmx_fio_getxdr(fp),
2465 flags_dfh, &state.dfhist, out);
2470 ret = do_cpt_EDstate(gmx_fio_getxdr(fp), TRUE, &state.edsamstate, out);
2475 ret = do_cpt_swapstate(gmx_fio_getxdr(fp), TRUE, &state.swapstate, out);
2480 do_cpt_files(gmx_fio_getxdr(fp), TRUE, &outputfiles, &nfiles, out, file_version);
2485 ret = do_cpt_footer(gmx_fio_getxdr(fp), file_version);
2492 if (gmx_fio_close(fp) != 0)
2494 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2501 static gmx_bool exist_output_file(const char *fnm_cp, int nfile, const t_filenm fnm[])
2505 /* Check if the output file name stored in the checkpoint file
2506 * is one of the output file names of mdrun.
2510 !(is_output(&fnm[i]) && strcmp(fnm_cp, fnm[i].fns[0]) == 0))
2515 return (i < nfile && gmx_fexist(fnm_cp));
2518 /* This routine cannot print tons of data, since it is called before the log file is opened. */
2519 gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
2520 gmx_int64_t *cpt_step, t_commrec *cr,
2521 gmx_bool bAppendReq,
2522 int nfile, const t_filenm fnm[],
2523 const char *part_suffix, gmx_bool *bAddPart)
2526 gmx_int64_t step = 0;
2528 /* This next line is nasty because the sub-structures of t_state
2529 * cannot be assumed to be zeroed (or even initialized in ways the
2530 * rest of the code might assume). Using snew would be better, but
2531 * this will all go away for 5.0. */
2534 gmx_file_position_t *outputfiles;
2537 char *fn, suf_up[STRLEN];
2543 if (!gmx_fexist(filename) || (!(fp = gmx_fio_open(filename, "r")) ))
2545 *simulation_part = 0;
2549 init_state(&state, 0, 0, 0, 0, 0);
2551 read_checkpoint_data(fp, simulation_part, &step, &t, &state,
2552 &nfiles, &outputfiles);
2553 if (gmx_fio_close(fp) != 0)
2555 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2562 for (f = 0; f < nfiles; f++)
2564 if (exist_output_file(outputfiles[f].filename, nfile, fnm))
2569 if (nexist == nfiles)
2571 bAppend = bAppendReq;
2573 else if (nexist > 0)
2576 "Output file appending has been requested,\n"
2577 "but some output files listed in the checkpoint file %s\n"
2578 "are not present or are named differently by the current program:\n",
2580 fprintf(stderr, "output files present:");
2581 for (f = 0; f < nfiles; f++)
2583 if (exist_output_file(outputfiles[f].filename,
2586 fprintf(stderr, " %s", outputfiles[f].filename);
2589 fprintf(stderr, "\n");
2590 fprintf(stderr, "output files not present or named differently:");
2591 for (f = 0; f < nfiles; f++)
2593 if (!exist_output_file(outputfiles[f].filename,
2596 fprintf(stderr, " %s", outputfiles[f].filename);
2599 fprintf(stderr, "\n");
2601 gmx_fatal(FARGS, "File appending requested, but only %d of the %d output files are present", nexist, nfiles);
2609 gmx_fatal(FARGS, "File appending requested, but no output file information is stored in the checkpoint file");
2611 fn = outputfiles[0].filename;
2612 if (strlen(fn) < 4 ||
2613 gmx_strcasecmp(fn+strlen(fn)-4, ftp2ext(efLOG)) == 0)
2615 gmx_fatal(FARGS, "File appending requested, but the log file is not the first file listed in the checkpoint file");
2617 /* Set bAddPart to whether the suffix string '.part' is present
2618 * in the log file name.
2620 strcpy(suf_up, part_suffix);
2622 *bAddPart = (strstr(fn, part_suffix) != NULL ||
2623 strstr(fn, suf_up) != NULL);
2631 gmx_bcast(sizeof(*simulation_part), simulation_part, cr);
2633 if (*simulation_part > 0 && bAppendReq)
2635 gmx_bcast(sizeof(bAppend), &bAppend, cr);
2636 gmx_bcast(sizeof(*bAddPart), bAddPart, cr);
2639 if (NULL != cpt_step)