1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This file is part of Gromacs Copyright (c) 1991-2008
5 * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
7 * This program is free software; you can redistribute it and/or
8 * modify it under the terms of the GNU General Public License
9 * as published by the Free Software Foundation; either version 2
10 * of the License, or (at your option) any later version.
12 * To help us fund GROMACS development, we humbly ask that you cite
13 * the research papers on the package. Check out http://www.gromacs.org
16 * Gnomes, ROck Monsters And Chili Sauce
26 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
29 #include <sys/locking.h>
43 #include "gmx_random.h"
44 #include "checkpoint.h"
54 #define CPT_MAGIC1 171817
55 #define CPT_MAGIC2 171819
57 /* The source code in this file should be thread-safe.
58 Please keep it that way. */
60 /* cpt_version should normally only be changed
61 * when the header of footer format changes.
62 * The state data format itself is backward and forward compatible.
63 * But old code can not read a new entry that is present in the file
64 * (but can read a new format when new entries are not present).
66 static const int cpt_version
= 8;
68 enum { ecpdtINT
, ecpdtFLOAT
, ecpdtDOUBLE
, ecpdtNR
};
70 const char *ecpdt_names
[ecpdtNR
] = { "int", "float", "double" };
72 const char *est_names
[estNR
]=
75 "box", "box-rel", "box-v", "pres_prev",
76 "nosehoover-xi", "thermostat-integral",
77 "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
78 "disre_initf", "disre_rm3tav",
79 "orire_initf", "orire_Dtav",
82 enum { eeksEKINH_N
, eeksEKINH
, eeksDEKINDL
, eeksMVCOS
, eeksNR
};
84 const char *eeks_names
[eeksNR
]=
86 "Ekinh_n", "Ekinh", "dEkindlambda", "mv_cos"
89 enum { eenhENERGY_N
, eenhENERGY_AVER
, eenhENERGY_SUM
, eenhENERGY_NSUM
,
90 eenhENERGY_SUM_SIM
, eenhENERGY_NSUM_SIM
,
91 eenhENERGY_NSTEPS
, eenhENERGY_NSTEPS_SIM
, eenhNR
};
93 const char *eenh_names
[eenhNR
]=
95 "energy_n", "energy_aver", "energy_sum", "energy_nsum",
96 "energy_sum_sim", "energy_nsum_sim",
97 "energy_nsteps", "energy_nsteps_sim"
102 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
104 gmx_wintruncate(const char *filename
, __int64 size
)
109 fp
=fopen(filename
,"r+");
116 return _chsize_s( fileno(fp
), size
);
122 enum { ecprREAL
, ecprRVEC
, ecprMATRIX
};
124 static const char *st_names(int cptp
,int ecpt
)
128 case 0: return est_names
[ecpt
]; break;
129 case 1: return eeks_names
[ecpt
]; break;
130 case 2: return eenh_names
[ecpt
]; break;
136 static void cp_warning(FILE *fp
)
138 fprintf(fp
,"\nWARNING: Checkpoint file is corrupted or truncated\n\n");
141 static void cp_error()
143 gmx_fatal(FARGS
,"Checkpoint file corrupted/truncated, or maybe you are out of quota?");
146 static void do_cpt_string_err(XDR
*xd
,bool bRead
,const char *desc
,char **s
,FILE *list
)
148 #define CPTSTRLEN 1024
155 res
= xdr_string(xd
,s
,CPTSTRLEN
);
162 fprintf(list
,"%s = %s\n",desc
,*s
);
167 static int do_cpt_int(XDR
*xd
,const char *desc
,int *i
,FILE *list
)
178 fprintf(list
,"%s = %d\n",desc
,*i
);
183 static int do_cpt_u_chars(XDR
*xd
,const char *desc
,int n
,unsigned char *i
,FILE *list
)
189 fprintf(list
,"%s = ",desc
);
191 for (j
=0; j
<n
&& res
; j
++)
193 res
&= xdr_u_char(xd
,&i
[j
]);
196 fprintf(list
,"%02x",i
[j
]);
211 static void do_cpt_int_err(XDR
*xd
,const char *desc
,int *i
,FILE *list
)
213 if (do_cpt_int(xd
,desc
,i
,list
) < 0)
219 static void do_cpt_step_err(XDR
*xd
,const char *desc
,gmx_large_int_t
*i
,FILE *list
)
224 res
= xdr_gmx_large_int(xd
,i
,"reading checkpoint file");
231 fprintf(list
,"%s = %s\n",desc
,gmx_step_str(*i
,buf
));
235 static void do_cpt_double_err(XDR
*xd
,const char *desc
,double *f
,FILE *list
)
239 res
= xdr_double(xd
,f
);
246 fprintf(list
,"%s = %f\n",desc
,*f
);
251 static int do_cpte_reals_low(XDR
*xd
,int cptp
,int ecpt
,int sflags
,
253 FILE *list
,int erealtype
)
267 res
= xdr_int(xd
,&nf
);
272 if (list
== NULL
&& nf
!= n
)
274 gmx_fatal(FARGS
,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp
,ecpt
),n
,nf
);
277 res
= xdr_int(xd
,&dt
);
284 fprintf(stderr
,"Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
285 st_names(cptp
,ecpt
),ecpdt_names
[dtc
],ecpdt_names
[dt
]);
287 if (list
|| !(sflags
& (1<<ecpt
)))
300 if (dt
== ecpdtFLOAT
)
302 if (dtc
== ecpdtFLOAT
)
310 res
= xdr_vector(xd
,(char *)vf
,nf
,
311 (unsigned int)sizeof(float),(xdrproc_t
)xdr_float
);
316 if (dtc
!= ecpdtFLOAT
)
327 if (dtc
== ecpdtDOUBLE
)
335 res
= xdr_vector(xd
,(char *)vd
,nf
,
336 (unsigned int)sizeof(double),(xdrproc_t
)xdr_double
);
341 if (dtc
!= ecpdtDOUBLE
)
356 pr_reals(list
,0,st_names(cptp
,ecpt
),vp
,nf
);
359 pr_rvecs(list
,0,st_names(cptp
,ecpt
),(rvec
*)vp
,nf
/3);
362 gmx_incons("Unknown checkpoint real type");
375 static int do_cpte_reals(XDR
*xd
,int cptp
,int ecpt
,int sflags
,
376 int n
,real
**v
,FILE *list
)
378 return do_cpte_reals_low(xd
,cptp
,ecpt
,sflags
,n
,v
,list
,ecprREAL
);
381 static int do_cpte_real(XDR
*xd
,int cptp
,int ecpt
,int sflags
,
384 return do_cpte_reals_low(xd
,cptp
,ecpt
,sflags
,1,&r
,list
,ecprREAL
);
387 static int do_cpte_ints(XDR
*xd
,int cptp
,int ecpt
,int sflags
,
388 int n
,int **v
,FILE *list
)
396 res
= xdr_int(xd
,&nf
);
401 if (list
== NULL
&& v
!= NULL
&& nf
!= n
)
403 gmx_fatal(FARGS
,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp
,ecpt
),n
,nf
);
406 res
= xdr_int(xd
,&dt
);
413 gmx_fatal(FARGS
,"Type mismatch for state entry %s, code type is %s, file type is %s\n",
414 st_names(cptp
,ecpt
),ecpdt_names
[dtc
],ecpdt_names
[dt
]);
416 if (list
|| !(sflags
& (1<<ecpt
)) || v
== NULL
)
429 res
= xdr_vector(xd
,(char *)vp
,nf
,
430 (unsigned int)sizeof(int),(xdrproc_t
)xdr_int
);
437 pr_ivec(list
,0,st_names(cptp
,ecpt
),vp
,nf
,TRUE
);
447 static int do_cpte_int(XDR
*xd
,int cptp
,int ecpt
,int sflags
,
450 return do_cpte_ints(xd
,cptp
,ecpt
,sflags
,1,&i
,list
);
453 static int do_cpte_doubles(XDR
*xd
,int cptp
,int ecpt
,int sflags
,
454 int n
,double **v
,FILE *list
)
462 res
= xdr_int(xd
,&nf
);
467 if (list
== NULL
&& nf
!= n
)
469 gmx_fatal(FARGS
,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp
,ecpt
),n
,nf
);
472 res
= xdr_int(xd
,&dt
);
479 gmx_fatal(FARGS
,"Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
480 st_names(cptp
,ecpt
),ecpdt_names
[dtc
],ecpdt_names
[dt
]);
482 if (list
|| !(sflags
& (1<<ecpt
)))
495 res
= xdr_vector(xd
,(char *)vp
,nf
,
496 (unsigned int)sizeof(double),(xdrproc_t
)xdr_double
);
503 pr_doubles(list
,0,st_names(cptp
,ecpt
),vp
,nf
);
514 static int do_cpte_rvecs(XDR
*xd
,int cptp
,int ecpt
,int sflags
,
515 int n
,rvec
**v
,FILE *list
)
517 return do_cpte_reals_low(xd
,cptp
,ecpt
,sflags
,
518 n
*DIM
,(real
**)v
,list
,ecprRVEC
);
521 static int do_cpte_matrix(XDR
*xd
,int cptp
,int ecpt
,int sflags
,
527 vr
= (real
*)&(v
[0][0]);
528 ret
= do_cpte_reals_low(xd
,cptp
,ecpt
,sflags
,DIM
*DIM
,&vr
,NULL
,ecprMATRIX
);
530 if (list
&& ret
== 0)
532 pr_rvecs(list
,0,st_names(cptp
,ecpt
),v
,DIM
);
538 static int do_cpte_matrices(XDR
*xd
,int cptp
,int ecpt
,int sflags
,
539 int n
,matrix
**v
,FILE *list
)
548 res
= xdr_int(xd
,&nf
);
553 if (list
== NULL
&& nf
!= n
)
555 gmx_fatal(FARGS
,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp
,ecpt
),n
,nf
);
557 if (list
|| !(sflags
& (1<<ecpt
)))
577 vr
[(i
*DIM
+j
)*DIM
+k
] = vp
[i
][j
][k
];
581 ret
= do_cpte_reals_low(xd
,cptp
,ecpt
,sflags
,
582 nf
*DIM
*DIM
,&vr
,NULL
,ecprMATRIX
);
589 vp
[i
][j
][k
] = vr
[(i
*DIM
+j
)*DIM
+k
];
595 if (list
&& ret
== 0)
599 pr_rvecs(list
,0,st_names(cptp
,ecpt
),vp
[i
],DIM
);
610 static void do_cpt_header(XDR
*xd
,bool bRead
,int *file_version
,
611 char **version
,char **btime
,char **buser
,char **bmach
,
612 char **fprog
,char **ftime
,
613 int *eIntegrator
,int *simulation_part
,
614 gmx_large_int_t
*step
,double *t
,
615 int *nnodes
,int *dd_nc
,int *npme
,
616 int *natoms
,int *ngtc
,
617 int *flags_state
,int *flags_eks
,int *flags_enh
,
633 res
= xdr_int(xd
,&magic
);
636 gmx_fatal(FARGS
,"The checkpoint file is empty/corrupted, or maybe you are out of quota?");
638 if (magic
!= CPT_MAGIC1
)
640 gmx_fatal(FARGS
,"Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
641 "The checkpoint file is corrupted or not a checkpoint file",
644 do_cpt_string_err(xd
,bRead
,"GROMACS version" ,version
,list
);
645 do_cpt_string_err(xd
,bRead
,"GROMACS build time" ,btime
,list
);
646 do_cpt_string_err(xd
,bRead
,"GROMACS build user" ,buser
,list
);
647 do_cpt_string_err(xd
,bRead
,"GROMACS build machine" ,bmach
,list
);
648 do_cpt_string_err(xd
,bRead
,"generating program" ,fprog
,list
);
649 do_cpt_string_err(xd
,bRead
,"generation time" ,ftime
,list
);
650 *file_version
= cpt_version
;
651 do_cpt_int_err(xd
,"checkpoint file version",file_version
,list
);
652 if (*file_version
> cpt_version
)
654 gmx_fatal(FARGS
,"Attempting to read a checkpoint file of version %d with code of version %d\n",*file_version
,cpt_version
);
656 do_cpt_int_err(xd
,"#atoms" ,natoms
,list
);
657 do_cpt_int_err(xd
,"#T-coupling groups",ngtc
,list
);
658 do_cpt_int_err(xd
,"integrator" ,eIntegrator
,list
);
659 if (*file_version
>= 3)
661 do_cpt_int_err(xd
,"simulation part #", simulation_part
,list
);
665 *simulation_part
= 1;
667 if (*file_version
>= 5)
669 do_cpt_step_err(xd
,"step" ,step
,list
);
673 do_cpt_int_err(xd
,"step" ,&idum
,list
);
676 do_cpt_double_err(xd
,"t" ,t
,list
);
677 do_cpt_int_err(xd
,"#PP-nodes" ,nnodes
,list
);
679 do_cpt_int_err(xd
,"dd_nc[x]",dd_nc
? &(dd_nc
[0]) : &idum
,list
);
680 do_cpt_int_err(xd
,"dd_nc[y]",dd_nc
? &(dd_nc
[1]) : &idum
,list
);
681 do_cpt_int_err(xd
,"dd_nc[z]",dd_nc
? &(dd_nc
[2]) : &idum
,list
);
682 do_cpt_int_err(xd
,"#PME-only nodes",npme
,list
);
683 do_cpt_int_err(xd
,"state flags",flags_state
,list
);
684 if (*file_version
>= 4)
686 do_cpt_int_err(xd
,"ekin data flags",flags_eks
,list
);
687 do_cpt_int_err(xd
,"energy history flags",flags_enh
,list
);
692 *flags_enh
= (*flags_state
>> (estORIRE_DTAV
+1));
693 *flags_state
= (*flags_state
& ~((1<<(estORIRE_DTAV
+1)) |
694 (1<<(estORIRE_DTAV
+2)) |
695 (1<<(estORIRE_DTAV
+3))));
699 static int do_cpt_footer(XDR
*xd
,bool bRead
,int file_version
)
704 if (file_version
>= 2)
707 res
= xdr_int(xd
,&magic
);
712 if (magic
!= CPT_MAGIC2
)
721 static int do_cpt_state(XDR
*xd
,bool bRead
,
722 int fflags
,t_state
*state
,
723 bool bReadRNG
,FILE *list
)
726 int **rng_p
,**rngi_p
;
734 rng_p
= (int **)&state
->ld_rng
;
735 rngi_p
= &state
->ld_rngi
;
739 /* Do not read the RNG data */
744 sflags
= state
->flags
;
745 for(i
=0; (i
<estNR
&& ret
== 0); i
++)
751 case estLAMBDA
: ret
= do_cpte_real (xd
,0,i
,sflags
,&state
->lambda
,list
); break;
752 case estBOX
: ret
= do_cpte_matrix(xd
,0,i
,sflags
,state
->box
,list
); break;
753 case estBOX_REL
: ret
= do_cpte_matrix(xd
,0,i
,sflags
,state
->box_rel
,list
); break;
754 case estBOXV
: ret
= do_cpte_matrix(xd
,0,i
,sflags
,state
->boxv
,list
); break;
755 case estPRES_PREV
: ret
= do_cpte_matrix(xd
,0,i
,sflags
,state
->pres_prev
,list
); break;
756 case estNH_XI
: ret
= do_cpte_reals (xd
,0,i
,sflags
,state
->ngtc
,&state
->nosehoover_xi
,list
); break;
757 case estTC_INT
: ret
= do_cpte_doubles(xd
,0,i
,sflags
,state
->ngtc
,&state
->therm_integral
,list
); break;
758 case estX
: ret
= do_cpte_rvecs (xd
,0,i
,sflags
,state
->natoms
,&state
->x
,list
); break;
759 case estV
: ret
= do_cpte_rvecs (xd
,0,i
,sflags
,state
->natoms
,&state
->v
,list
); break;
760 case estSDX
: ret
= do_cpte_rvecs (xd
,0,i
,sflags
,state
->natoms
,&state
->sd_X
,list
); break;
761 case estLD_RNG
: ret
= do_cpte_ints (xd
,0,i
,sflags
,state
->nrng
,rng_p
,list
); break;
762 case estLD_RNGI
: ret
= do_cpte_ints (xd
,0,i
,sflags
,state
->nrngi
,rngi_p
,list
); break;
763 case estDISRE_INITF
: ret
= do_cpte_real (xd
,0,i
,sflags
,&state
->hist
.disre_initf
,list
); break;
764 case estDISRE_RM3TAV
: ret
= do_cpte_reals(xd
,0,i
,sflags
,state
->hist
.ndisrepairs
,&state
->hist
.disre_rm3tav
,list
); break;
765 case estORIRE_INITF
: ret
= do_cpte_real (xd
,0,i
,sflags
,&state
->hist
.orire_initf
,list
); break;
766 case estORIRE_DTAV
: ret
= do_cpte_reals(xd
,0,i
,sflags
,state
->hist
.norire_Dtav
,&state
->hist
.orire_Dtav
,list
); break;
768 gmx_fatal(FARGS
,"Unknown state entry %d\n"
769 "You are probably reading a new checkpoint file with old code",i
);
777 static int do_cpt_ekinstate(XDR
*xd
,bool bRead
,
778 int fflags
,ekinstate_t
*ekins
,
786 for(i
=0; (i
<eeksNR
&& ret
== 0); i
++)
793 case eeksEKINH_N
: ret
= do_cpte_int(xd
,1,i
,fflags
,&ekins
->ekinh_n
,list
); break;
794 case eeksEKINH
: ret
= do_cpte_matrices(xd
,1,i
,fflags
,ekins
->ekinh_n
,&ekins
->ekinh
,list
); break;
795 case eeksDEKINDL
: ret
= do_cpte_real(xd
,1,i
,fflags
,&ekins
->dekindl
,list
); break;
796 case eeksMVCOS
: ret
= do_cpte_real(xd
,1,i
,fflags
,&ekins
->mvcos
,list
); break;
798 gmx_fatal(FARGS
,"Unknown ekin data state entry %d\n"
799 "You are probably reading a new checkpoint file with old code",i
);
807 static int do_cpt_enerhist(XDR
*xd
,bool bRead
,
808 int fflags
,energyhistory_t
*enerhist
,
818 enerhist
->nsteps
= 0;
820 enerhist
->nsteps_sim
= 0;
821 enerhist
->nsum_sim
= 0;
824 for(i
=0; (i
<eenhNR
&& ret
== 0); i
++)
830 case eenhENERGY_N
: ret
= do_cpte_int(xd
,2,i
,fflags
,&enerhist
->nener
,list
); break;
831 case eenhENERGY_AVER
: ret
= do_cpte_doubles(xd
,2,i
,fflags
,enerhist
->nener
,&enerhist
->ener_ave
,list
); break;
832 case eenhENERGY_SUM
: ret
= do_cpte_doubles(xd
,2,i
,fflags
,enerhist
->nener
,&enerhist
->ener_sum
,list
); break;
833 case eenhENERGY_NSUM
: do_cpt_step_err(xd
,eenh_names
[i
],&enerhist
->nsum
,list
); break;
834 case eenhENERGY_SUM_SIM
: ret
= do_cpte_doubles(xd
,2,i
,fflags
,enerhist
->nener
,&enerhist
->ener_sum_sim
,list
); break;
835 case eenhENERGY_NSUM_SIM
: do_cpt_step_err(xd
,eenh_names
[i
],&enerhist
->nsum_sim
,list
); break;
836 case eenhENERGY_NSTEPS
: do_cpt_step_err(xd
,eenh_names
[i
],&enerhist
->nsteps
,list
); break;
837 case eenhENERGY_NSTEPS_SIM
: do_cpt_step_err(xd
,eenh_names
[i
],&enerhist
->nsteps_sim
,list
); break;
839 gmx_fatal(FARGS
,"Unknown energy history entry %d\n"
840 "You are probably reading a new checkpoint file with old code",i
);
845 if ((fflags
& (1<<eenhENERGY_SUM
)) && !(fflags
& (1<<eenhENERGY_SUM_SIM
)))
847 /* Assume we have an old file format and copy sum to sum_sim */
848 srenew(enerhist
->ener_sum_sim
,enerhist
->nener
);
849 for(i
=0; i
<enerhist
->nener
; i
++)
851 enerhist
->ener_sum_sim
[i
] = enerhist
->ener_sum
[i
];
853 fflags
|= (1<<eenhENERGY_SUM_SIM
);
856 if ( (fflags
& (1<<eenhENERGY_NSUM
)) &&
857 !(fflags
& (1<<eenhENERGY_NSTEPS
)))
859 /* Assume we have an old file format and copy nsum to nsteps */
860 enerhist
->nsteps
= enerhist
->nsum
;
861 fflags
|= (1<<eenhENERGY_NSTEPS
);
863 if ( (fflags
& (1<<eenhENERGY_NSUM_SIM
)) &&
864 !(fflags
& (1<<eenhENERGY_NSTEPS_SIM
)))
866 /* Assume we have an old file format and copy nsum to nsteps */
867 enerhist
->nsteps_sim
= enerhist
->nsum_sim
;
868 fflags
|= (1<<eenhENERGY_NSTEPS_SIM
);
874 static int do_cpt_files(XDR
*xd
, bool bRead
,
875 gmx_file_position_t
**p_outputfiles
, int *nfiles
,
876 FILE *list
, int file_version
)
880 off_t mask
= 0xFFFFFFFFL
;
881 int offset_high
,offset_low
;
883 gmx_file_position_t
*outputfiles
;
885 if (do_cpt_int(xd
,"number of output files",nfiles
,list
) != 0)
892 snew(*p_outputfiles
,*nfiles
);
895 outputfiles
= *p_outputfiles
;
897 for(i
=0;i
<*nfiles
;i
++)
899 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
902 do_cpt_string_err(xd
,bRead
,"output filename",&buf
,list
);
903 strncpy(outputfiles
[i
].filename
,buf
,CPTSTRLEN
-1);
909 if (do_cpt_int(xd
,"file_offset_high",&offset_high
,list
) != 0)
913 if (do_cpt_int(xd
,"file_offset_low",&offset_low
,list
) != 0)
917 #if (SIZEOF_OFF_T > 4)
918 outputfiles
[i
].offset
= ( ((off_t
) offset_high
) << 32 ) | ( (off_t
) offset_low
& mask
);
920 outputfiles
[i
].offset
= offset_low
;
925 buf
= outputfiles
[i
].filename
;
926 do_cpt_string_err(xd
,bRead
,"output filename",&buf
,list
);
928 offset
= outputfiles
[i
].offset
;
936 #if (SIZEOF_OFF_T > 4)
937 offset_low
= (int) (offset
& mask
);
938 offset_high
= (int) ((offset
>> 32) & mask
);
944 if (do_cpt_int(xd
,"file_offset_high",&offset_high
,list
) != 0)
948 if (do_cpt_int(xd
,"file_offset_low",&offset_low
,list
) != 0)
953 if (file_version
>= 8)
955 if (do_cpt_int(xd
,"file_checksum_size",&(outputfiles
[i
].chksum_size
),
960 if (do_cpt_u_chars(xd
,"file_checksum",16,outputfiles
[i
].chksum
,list
) != 0)
967 outputfiles
[i
].chksum_size
= -1;
974 void write_checkpoint(const char *fn
,FILE *fplog
,t_commrec
*cr
,
975 int eIntegrator
,int simulation_part
,
976 gmx_large_int_t step
,double t
,t_state
*state
)
987 int nppnodes
,npmenodes
,flag_64bit
;
989 gmx_file_position_t
*outputfiles
;
991 int flags_eks
,flags_enh
,i
;
995 if (DOMAINDECOMP(cr
))
997 nppnodes
= cr
->dd
->nnodes
;
998 npmenodes
= cr
->npmenodes
;
1002 nppnodes
= cr
->nnodes
;
1014 /* Rename the previous checkpoint file */
1016 buf
[strlen(fn
) - strlen(ftp2ext(fn2ftp(fn
))) - 1] = '\0';
1017 strcat(buf
,"_prev");
1018 strcat(buf
,fn
+strlen(fn
) - strlen(ftp2ext(fn2ftp(fn
))) - 1);
1019 (void)remove(buf
); /* Unix will overwrite buf if it exists, but for windows we need to remove first */
1020 if(rename(fn
,buf
) != 0)
1022 gmx_file("Cannot rename checkpoint file; maybe you are out of quota?");
1027 ftime
= strdup(ctime(&now
));
1028 ftime
[strlen(ftime
)-1] = '\0';
1030 /* No need to pollute stderr every time we write a checkpoint file */
1031 /* fprintf(stderr,"\nWriting checkpoint, step %d at %s\n",step,ftime); */
1034 fprintf(fplog
,"Writing checkpoint, step %s at %s\n\n",
1035 gmx_step_str(step
,buf
),ftime
);
1038 /* Get offsets for open files */
1039 gmx_fio_get_output_file_positions(&outputfiles
, &noutputfiles
);
1041 fp
= gmx_fio_open(fn
,"w");
1043 if (state
->ekinstate
.bUpToDate
)
1046 ((1<<eeksEKINH_N
) | (1<<eeksEKINH
) | (1<<eeksDEKINDL
) |
1055 if (state
->enerhist
.nsum
> 0 || state
->enerhist
.nsum_sim
> 0)
1057 flags_enh
|= (1<<eenhENERGY_N
);
1058 if (state
->enerhist
.nsum
> 0)
1060 flags_enh
|= ((1<<eenhENERGY_AVER
) | (1<<eenhENERGY_SUM
) |
1061 (1<<eenhENERGY_NSTEPS
) | (1<<eenhENERGY_NSUM
));
1063 if (state
->enerhist
.nsum_sim
> 0)
1065 flags_enh
|= ((1<<eenhENERGY_SUM_SIM
) | (1<<eenhENERGY_NSTEPS_SIM
) |
1066 (1<<eenhENERGY_NSUM_SIM
));
1071 version
= strdup(VERSION
);
1072 btime
= strdup(BUILD_TIME
);
1073 buser
= strdup(BUILD_USER
);
1074 bmach
= strdup(BUILD_MACHINE
);
1075 fprog
= strdup(Program());
1077 do_cpt_header(gmx_fio_getxdr(fp
),FALSE
,&file_version
,
1078 &version
,&btime
,&buser
,&bmach
,&fprog
,&ftime
,
1079 &eIntegrator
,&simulation_part
,&step
,&t
,&nppnodes
,
1080 DOMAINDECOMP(cr
) ? cr
->dd
->nc
: NULL
,&npmenodes
,
1081 &state
->natoms
,&state
->ngtc
,
1082 &state
->flags
,&flags_eks
,&flags_enh
,NULL
);
1090 if( (do_cpt_state(gmx_fio_getxdr(fp
),FALSE
,state
->flags
,state
,TRUE
,NULL
) < 0) ||
1091 (do_cpt_ekinstate(gmx_fio_getxdr(fp
),FALSE
,flags_eks
,&state
->ekinstate
,NULL
) < 0) ||
1092 (do_cpt_enerhist(gmx_fio_getxdr(fp
),FALSE
,flags_enh
,&state
->enerhist
,NULL
) < 0) ||
1093 (do_cpt_files(gmx_fio_getxdr(fp
),FALSE
,&outputfiles
,&noutputfiles
,NULL
,file_version
) < 0))
1095 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1098 do_cpt_footer(gmx_fio_getxdr(fp
),FALSE
,file_version
);
1100 if( gmx_fio_close(fp
) != 0)
1102 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1108 /*code for alternate checkpointing scheme. moved from top of loop over steps */
1109 fcRequestCheckPoint();
1110 if ( fcCheckPointParallel( cr
->nodeid
, NULL
,0) == 0 ) {
1111 gmx_fatal( 3,__FILE__
,__LINE__
, "Checkpoint error on step %d\n", step
);
1113 #endif /* end FAHCORE block */
1117 static void print_flag_mismatch(FILE *fplog
,int sflags
,int fflags
)
1121 fprintf(fplog
,"\nState entry mismatch between the simulation and the checkpoint file\n");
1122 fprintf(fplog
,"Entries which are not present in the checkpoint file will not be updated\n");
1123 fprintf(fplog
," %24s %11s %11s\n","","simulation","checkpoint");
1124 for(i
=0; i
<estNR
; i
++)
1126 if ((sflags
& (1<<i
)) || (fflags
& (1<<i
)))
1128 fprintf(fplog
," %24s %11s %11s\n",
1130 (sflags
& (1<<i
)) ? " present " : "not present",
1131 (fflags
& (1<<i
)) ? " present " : "not present");
1136 static void check_int(FILE *fplog
,const char *type
,int p
,int f
,bool *mm
)
1138 FILE *fp
= fplog
? fplog
: stderr
;
1142 fprintf(fp
," %s mismatch,\n",type
);
1143 fprintf(fp
," current program: %d\n",p
);
1144 fprintf(fp
," checkpoint file: %d\n",f
);
1150 static void check_string(FILE *fplog
,const char *type
,const char *p
,
1151 const char *f
,bool *mm
)
1153 FILE *fp
= fplog
? fplog
: stderr
;
1155 if (strcmp(p
,f
) != 0)
1157 fprintf(fp
," %s mismatch,\n",type
);
1158 fprintf(fp
," current program: %s\n",p
);
1159 fprintf(fp
," checkpoint file: %s\n",f
);
1165 static void check_match(FILE *fplog
,
1167 char *btime
,char *buser
,char *bmach
,char *fprog
,
1168 t_commrec
*cr
,bool bPartDecomp
,int npp_f
,int npme_f
,
1169 ivec dd_nc
,ivec dd_nc_f
)
1176 check_string(fplog
,"Version" ,VERSION
,version
,&mm
);
1177 check_string(fplog
,"Build time" ,BUILD_TIME
,btime
,&mm
);
1178 check_string(fplog
,"Build user" ,BUILD_USER
,buser
,&mm
);
1179 check_string(fplog
,"Build machine",BUILD_MACHINE
,bmach
,&mm
);
1180 check_string(fplog
,"Program name" ,Program() ,fprog
,&mm
);
1182 npp
= cr
->nnodes
- cr
->npmenodes
;
1183 check_int (fplog
,"#nodes" ,cr
->nnodes
,npp_f
+npme_f
,&mm
);
1192 check_int (fplog
,"#PME-nodes" ,cr
->npmenodes
,npme_f
,&mm
);
1195 check_int (fplog
,"#DD-cells[x]",dd_nc
[XX
] ,dd_nc_f
[XX
],&mm
);
1196 check_int (fplog
,"#DD-cells[y]",dd_nc
[YY
] ,dd_nc_f
[YY
],&mm
);
1197 check_int (fplog
,"#DD-cells[z]",dd_nc
[ZZ
] ,dd_nc_f
[ZZ
],&mm
);
1204 "Gromacs binary or parallel settings not identical to previous run.\n"
1205 "Continuation is exact, but is not guaranteed to be binary identical%s.\n\n",
1206 fplog
? ",\n see the log file for details" : "");
1211 "Gromacs binary or parallel settings not identical to previous run.\n"
1212 "Continuation is exact, but is not guaranteed to be binary identical.\n\n");
1217 static void read_checkpoint(const char *fn
,FILE **pfplog
,
1218 t_commrec
*cr
,bool bPartDecomp
,ivec dd_nc
,
1219 int eIntegrator
,gmx_large_int_t
*step
,double *t
,
1220 t_state
*state
,bool *bReadRNG
,bool *bReadEkin
,
1221 int *simulation_part
,bool bAppendOutputFiles
)
1225 char *version
,*btime
,*buser
,*bmach
,*fprog
,*ftime
;
1226 char filename
[STRLEN
],buf
[22];
1227 int nppnodes
,eIntegrator_f
,nppnodes_f
,npmenodes_f
;
1229 int natoms
,ngtc
,fflags
,flags_eks
,flags_enh
;
1232 gmx_file_position_t
*outputfiles
;
1235 FILE* fplog
= *pfplog
;
1236 unsigned char digest
[16];
1237 #if !((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
1238 struct flock fl
= { F_WRLCK
, SEEK_SET
, 0, 0, 0 };
1241 const char *int_warn
=
1242 "WARNING: The checkpoint file was generator with integrator %s,\n"
1243 " while the simulation uses integrator %s\n\n";
1244 const char *sd_note
=
1245 "NOTE: The checkpoint file was for %d nodes doing SD or BD,\n"
1246 " while the simulation uses %d SD or BD nodes,\n"
1247 " continuation will be exact, except for the random state\n\n";
1252 "read_checkpoint not (yet) supported with particle decomposition");
1255 fp
= gmx_fio_open(fn
,"r");
1256 do_cpt_header(gmx_fio_getxdr(fp
),TRUE
,&file_version
,
1257 &version
,&btime
,&buser
,&bmach
,&fprog
,&ftime
,
1258 &eIntegrator_f
,simulation_part
,step
,t
,
1259 &nppnodes_f
,dd_nc_f
,&npmenodes_f
,
1261 &fflags
,&flags_eks
,&flags_enh
,NULL
);
1263 if (cr
== NULL
|| MASTER(cr
))
1265 fprintf(stderr
,"\nReading checkpoint file %s generated: %s\n\n",
1269 /* This will not be written if we do appending, since fplog is still NULL then */
1272 fprintf(fplog
,"\n");
1273 fprintf(fplog
,"Reading checkpoint file %s\n",fn
);
1274 fprintf(fplog
," file generated by: %s\n",fprog
);
1275 fprintf(fplog
," file generated at: %s\n",ftime
);
1276 fprintf(fplog
," GROMACS build time: %s\n",btime
);
1277 fprintf(fplog
," GROMACS build user: %s\n",buser
);
1278 fprintf(fplog
," GROMACS build machine: %s\n",bmach
);
1279 fprintf(fplog
," simulation part #: %d\n",*simulation_part
);
1280 fprintf(fplog
," step: %s\n",gmx_step_str(*step
,buf
));
1281 fprintf(fplog
," time: %f\n",*t
);
1282 fprintf(fplog
,"\n");
1285 if (natoms
!= state
->natoms
)
1287 gmx_fatal(FARGS
,"Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms",natoms
,state
->natoms
);
1289 if (ngtc
!= state
->ngtc
)
1291 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
);
1293 if (eIntegrator_f
!= eIntegrator
)
1297 fprintf(stderr
,int_warn
,EI(eIntegrator_f
),EI(eIntegrator
));
1299 if(bAppendOutputFiles
)
1302 "Output file appending requested, but input/checkpoint integrators do not match.\n"
1303 "Stopping the run to prevent you from ruining all your data...\n"
1304 "If you _really_ know what you are doing, try without the -append option.\n");
1308 fprintf(fplog
,int_warn
,EI(eIntegrator_f
),EI(eIntegrator
));
1316 else if (bPartDecomp
)
1318 nppnodes
= cr
->nnodes
;
1320 else if (cr
->nnodes
== nppnodes_f
+ npmenodes_f
)
1322 if (cr
->npmenodes
< 0)
1324 cr
->npmenodes
= npmenodes_f
;
1326 nppnodes
= cr
->nnodes
- cr
->npmenodes
;
1327 if (nppnodes
== nppnodes_f
)
1329 for(d
=0; d
<DIM
; d
++)
1333 dd_nc
[d
] = dd_nc_f
[d
];
1340 /* The number of PP nodes has not been set yet */
1344 if ((EI_SD(eIntegrator
) || eIntegrator
== eiBD
) && nppnodes
> 0)
1346 /* Correct the RNG state size for the number of PP nodes.
1347 * Such assignments should all be moved to one central function.
1349 state
->nrng
= nppnodes
*gmx_rng_n();
1350 state
->nrngi
= nppnodes
;
1354 if (fflags
!= state
->flags
)
1359 if(bAppendOutputFiles
)
1362 "Output file appending requested, but input and checkpoint states are not identical.\n"
1363 "Stopping the run to prevent you from ruining all your data...\n"
1364 "You can try without the -append option, and get more info in the log file.\n");
1368 "WARNING: The checkpoint state entries do not match the simulation,\n"
1369 " see the log file for details\n\n");
1374 print_flag_mismatch(fplog
,state
->flags
,fflags
);
1379 if ((EI_SD(eIntegrator
) || eIntegrator
== eiBD
) &&
1380 nppnodes
!= nppnodes_f
)
1385 fprintf(stderr
,sd_note
,nppnodes_f
,nppnodes
);
1389 fprintf(fplog
,sd_note
,nppnodes_f
,nppnodes
);
1394 check_match(fplog
,version
,btime
,buser
,bmach
,fprog
,
1395 cr
,bPartDecomp
,nppnodes_f
,npmenodes_f
,dd_nc
,dd_nc_f
);
1398 ret
= do_cpt_state(gmx_fio_getxdr(fp
),TRUE
,fflags
,state
,*bReadRNG
,NULL
);
1403 ret
= do_cpt_ekinstate(gmx_fio_getxdr(fp
),TRUE
,
1404 flags_eks
,&state
->ekinstate
,NULL
);
1409 *bReadEkin
= (flags_eks
& (1<<eeksEKINH
));
1411 ret
= do_cpt_enerhist(gmx_fio_getxdr(fp
),TRUE
,
1412 flags_enh
,&state
->enerhist
,NULL
);
1418 if (file_version
< 6)
1420 fprintf(stderr
,"\nWARNING: 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 energy averages will be incorrect.\n\n");
1423 fprintf(fplog
,"\nWARNING: 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 energy averages will be incorrect.\n\n");
1425 state
->enerhist
.nsum
= *step
;
1426 state
->enerhist
.nsum_sim
= *step
;
1429 ret
= do_cpt_files(gmx_fio_getxdr(fp
),TRUE
,&outputfiles
,&nfiles
,NULL
,file_version
);
1435 ret
= do_cpt_footer(gmx_fio_getxdr(fp
),TRUE
,file_version
);
1440 if( gmx_fio_close(fp
) != 0)
1442 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1451 /* If the user wants to append to output files (danger), we use the file pointer
1452 * positions of the output files stored in the checkpoint file and truncate the
1453 * files such that any frames written after the checkpoint time are removed.
1455 * You will get REALLY fun problems if you use the -append option by provide
1456 * mdrun with other input files (in-frame truncation in the wrong places). Suit yourself!
1458 if (bAppendOutputFiles
)
1460 if (fn2ftp(outputfiles
[0].filename
)!=efLOG
)
1462 /* make sure first file is log file so that it is OK to use it for
1465 gmx_fatal(FARGS
,"The first output file should always be the log "
1466 "file but instead is: %s", outputfiles
[0].filename
);
1468 for(i
=0;i
<nfiles
;i
++)
1470 if (outputfiles
[i
].filename
,outputfiles
[i
].offset
< 0)
1472 gmx_fatal(FARGS
,"The original run wrote a file called '%s' which "
1473 "is larger than 2 GB, but mdrun did not support large file"
1474 " offsets. Can not append. Run mdrun without -append",
1475 outputfiles
[i
].filename
);
1478 chksum_file
=gmx_fio_open(outputfiles
[i
].filename
,"a");
1481 chksum_file
=gmx_fio_open(outputfiles
[i
].filename
,"r+");
1486 #if !((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
1487 if (fcntl(fileno(gmx_fio_getfp(chksum_file
)), F_SETLK
, &fl
)
1490 if (_locking(fileno(gmx_fio_getfp(chksum_file
)), _LK_NBLCK
, LONG_MAX
)==-1)
1493 gmx_fatal(FARGS
,"Failed to lock: %s. Already running "
1494 "simulation?", outputfiles
[i
].filename
);
1498 /* compute md5 chksum */
1499 if (outputfiles
[i
].chksum_size
!= -1)
1501 if (gmx_fio_get_file_md5(chksum_file
,outputfiles
[i
].offset
,
1502 digest
) != outputfiles
[i
].chksum_size
)
1504 gmx_fatal(FARGS
,"Can't read %d bytes of '%s' to compute"
1505 " checksum.", outputfiles
[i
].chksum_size
,
1506 outputfiles
[i
].filename
);
1509 else if (i
==0) /*log file need to be seeked even when not reading md5*/
1511 gmx_fio_seek(chksum_file
,outputfiles
[i
].offset
);
1515 if (i
==0) /*open log file here - so that lock is never lifted
1516 after chksum is calculated */
1518 *pfplog
= gmx_fio_getfp(chksum_file
);
1522 gmx_fio_close(chksum_file
);
1525 /* compare md5 chksum */
1526 if (outputfiles
[i
].chksum_size
!= -1 &&
1527 memcmp(digest
,outputfiles
[i
].chksum
,16)!=0)
1531 fprintf(debug
,"chksum for %s: ",outputfiles
[i
].filename
);
1532 for (j
=0; j
<16; j
++)
1534 fprintf(debug
,"%02x",digest
[j
]);
1536 fprintf(debug
,"\n");
1538 gmx_fatal(FARGS
,"Checksum wrong for '%s'.",
1539 outputfiles
[i
].filename
);
1544 if (i
!=0) /*log file is already seeked to correct position */
1546 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
1547 rc
= gmx_wintruncate(outputfiles
[i
].filename
,outputfiles
[i
].offset
);
1549 rc
= truncate(outputfiles
[i
].filename
,outputfiles
[i
].offset
);
1553 gmx_fatal(FARGS
,"Truncation of file %s failed.",outputfiles
[i
].filename
);
1563 void load_checkpoint(const char *fn
,FILE **fplog
,
1564 t_commrec
*cr
,bool bPartDecomp
,ivec dd_nc
,
1565 t_inputrec
*ir
,t_state
*state
,
1566 bool *bReadRNG
,bool *bReadEkin
,bool bAppend
)
1568 gmx_large_int_t step
;
1571 if (SIMMASTER(cr
)) {
1572 /* Read the state from the checkpoint file */
1573 read_checkpoint(fn
,fplog
,
1574 cr
,bPartDecomp
,dd_nc
,
1575 ir
->eI
,&step
,&t
,state
,bReadRNG
,bReadEkin
,
1576 &ir
->simulation_part
,bAppend
);
1579 gmx_bcast(sizeof(cr
->npmenodes
),&cr
->npmenodes
,cr
);
1580 gmx_bcast(DIM
*sizeof(dd_nc
[0]),dd_nc
,cr
);
1581 gmx_bcast(sizeof(step
),&step
,cr
);
1582 gmx_bcast(sizeof(*bReadRNG
),bReadRNG
,cr
);
1583 gmx_bcast(sizeof(*bReadEkin
),bReadEkin
,cr
);
1585 ir
->bContinuation
= TRUE
;
1586 ir
->nsteps
+= ir
->init_step
- step
;
1587 ir
->init_step
= step
;
1588 ir
->simulation_part
+= 1;
1591 static void low_read_checkpoint_state(int fp
,int *simulation_part
,
1592 gmx_large_int_t
*step
,double *t
,t_state
*state
,
1596 char *version
,*btime
,*buser
,*bmach
,*fprog
,*ftime
;
1600 int flags_eks
,flags_enh
;
1602 gmx_file_position_t
*outputfiles
;
1605 do_cpt_header(gmx_fio_getxdr(fp
),TRUE
,&file_version
,
1606 &version
,&btime
,&buser
,&bmach
,&fprog
,&ftime
,
1607 &eIntegrator
,simulation_part
,step
,t
,&nppnodes
,dd_nc
,&npme
,
1608 &state
->natoms
,&state
->ngtc
,
1609 &state
->flags
,&flags_eks
,&flags_enh
,NULL
);
1611 do_cpt_state(gmx_fio_getxdr(fp
),TRUE
,state
->flags
,state
,bReadRNG
,NULL
);
1616 ret
= do_cpt_ekinstate(gmx_fio_getxdr(fp
),TRUE
,
1617 flags_eks
,&state
->ekinstate
,NULL
);
1622 ret
= do_cpt_enerhist(gmx_fio_getxdr(fp
),TRUE
,
1623 flags_enh
,&state
->enerhist
,NULL
);
1629 ret
= do_cpt_files(gmx_fio_getxdr(fp
),TRUE
,&outputfiles
,&nfiles
,NULL
,file_version
);
1637 ret
= do_cpt_footer(gmx_fio_getxdr(fp
),TRUE
,file_version
);
1651 read_checkpoint_state(const char *fn
,int *simulation_part
,
1652 gmx_large_int_t
*step
,double *t
,t_state
*state
)
1656 fp
= gmx_fio_open(fn
,"r");
1657 low_read_checkpoint_state(fp
,simulation_part
,step
,t
,state
,TRUE
);
1658 if( gmx_fio_close(fp
) != 0)
1660 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1664 void read_checkpoint_trxframe(int fp
,t_trxframe
*fr
)
1667 int simulation_part
;
1668 gmx_large_int_t step
;
1671 init_state(&state
,0,0);
1673 low_read_checkpoint_state(fp
,&simulation_part
,&step
,&t
,&state
,FALSE
);
1675 fr
->natoms
= state
.natoms
;
1678 fr
->step
= gmx_large_int_to_int(step
,
1679 "conversion of checkpoint to trajectory");
1683 fr
->lambda
= state
.lambda
;
1685 fr
->bX
= (state
.flags
& (1<<estX
));
1691 fr
->bV
= (state
.flags
& (1<<estV
));
1698 fr
->bBox
= (state
.flags
& (1<<estBOX
));
1701 copy_mat(state
.box
,fr
->box
);
1706 void list_checkpoint(const char *fn
,FILE *out
)
1710 char *version
,*btime
,*buser
,*bmach
,*fprog
,*ftime
;
1711 int eIntegrator
,simulation_part
,nppnodes
,npme
;
1712 gmx_large_int_t step
;
1716 int flags_eks
,flags_enh
;
1720 gmx_file_position_t
*outputfiles
;
1723 init_state(&state
,-1,-1);
1725 fp
= gmx_fio_open(fn
,"r");
1726 do_cpt_header(gmx_fio_getxdr(fp
),TRUE
,&file_version
,
1727 &version
,&btime
,&buser
,&bmach
,&fprog
,&ftime
,
1728 &eIntegrator
,&simulation_part
,&step
,&t
,&nppnodes
,dd_nc
,&npme
,
1729 &state
.natoms
,&state
.ngtc
,
1730 &state
.flags
,&flags_eks
,&flags_enh
,out
);
1731 ret
= do_cpt_state(gmx_fio_getxdr(fp
),TRUE
,state
.flags
,&state
,TRUE
,out
);
1736 ret
= do_cpt_ekinstate(gmx_fio_getxdr(fp
),TRUE
,
1737 flags_eks
,&state
.ekinstate
,out
);
1742 ret
= do_cpt_enerhist(gmx_fio_getxdr(fp
),TRUE
,
1743 flags_enh
,&state
.enerhist
,out
);
1747 do_cpt_files(gmx_fio_getxdr(fp
),TRUE
,&outputfiles
,&nfiles
,out
,file_version
);
1752 ret
= do_cpt_footer(gmx_fio_getxdr(fp
),TRUE
,file_version
);
1759 if( gmx_fio_close(fp
) != 0)
1761 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1768 /* This routine cannot print tons of data, since it is called before the log file is opened. */
1769 void read_checkpoint_simulation_part(const char *filename
, int *simulation_part
,
1770 gmx_large_int_t
*cpt_step
,t_commrec
*cr
)
1774 char *version
,*btime
,*buser
,*bmach
,*fprog
,*ftime
;
1775 int eIntegrator_f
,nppnodes_f
,npmenodes_f
;
1777 gmx_large_int_t step
=0;
1779 int natoms
,ngtc
,fflags
,flags_eks
,flags_enh
;
1781 if (SIMMASTER(cr
)) {
1782 if(!gmx_fexist(filename
) || ( (fp
= gmx_fio_open(filename
,"r")) < 0 ))
1784 *simulation_part
= 0;
1788 do_cpt_header(gmx_fio_getxdr(fp
),
1790 &version
,&btime
,&buser
,&bmach
,&fprog
,&ftime
,
1791 &eIntegrator_f
,simulation_part
,
1792 &step
,&t
,&nppnodes_f
,dd_nc_f
,&npmenodes_f
,
1794 &fflags
,&flags_eks
,&flags_enh
,NULL
);
1795 if( gmx_fio_close(fp
) != 0)
1797 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1809 gmx_bcast(sizeof(*simulation_part
),simulation_part
,cr
);
1811 if (NULL
!= cpt_step
)