Removed the unnecessary hacks from some boolean statements.
[alexxy/gromacs.git] / src / gmxlib / checkpoint.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  * This file is part of Gromacs        Copyright (c) 1991-2008
5  * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
6  *
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.
11  *
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
14  * 
15  * And Hey:
16  * Gnomes, ROck Monsters And Chili Sauce
17  */
18
19 /* The source code in this file should be thread-safe. 
20  Please keep it that way. */
21
22
23 #ifdef HAVE_CONFIG_H
24 #include <config.h>
25 #endif
26
27 #include <string.h>
28 #include <time.h>
29
30 #ifdef HAVE_SYS_TIME_H
31 #include <sys/time.h>
32 #endif
33
34
35 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
36 /* _chsize_s */
37 #include <io.h>
38 #include <sys/locking.h>
39 #endif
40
41
42 #include "filenm.h"
43 #include "names.h"
44 #include "typedefs.h"
45 #include "smalloc.h"
46 #include "gmxfio.h"
47 #include "xdrf.h"
48 #include "statutil.h"
49 #include "txtdump.h"
50 #include "vec.h"
51 #include "network.h"
52 #include "gmx_random.h"
53 #include "checkpoint.h"
54 #include "futil.h"
55 #include "string2.h"
56 #include <fcntl.h>
57
58
59 #ifdef GMX_FAHCORE
60 #include "corewrap.h"
61 #endif
62
63
64 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
65 char *
66 gmx_ctime_r(const time_t *clock,char *buf, int n);
67
68
69 #define CPT_MAGIC1 171817
70 #define CPT_MAGIC2 171819
71
72 /* cpt_version should normally only be changed
73  * when the header of footer format changes.
74  * The state data format itself is backward and forward compatible.
75  * But old code can not read a new entry that is present in the file
76  * (but can read a new format when new entries are not present).
77  */
78 static const int cpt_version = 12;
79
80
81 const char *est_names[estNR]=
82 {
83     "FE-lambda",
84     "box", "box-rel", "box-v", "pres_prev",
85     "nosehoover-xi", "thermostat-integral",
86     "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
87     "disre_initf", "disre_rm3tav",
88     "orire_initf", "orire_Dtav",
89     "svir_prev", "nosehoover-vxi", "v_eta", "vol0", "nhpres_xi", "nhpres_vxi", "fvir_prev",
90 };
91
92 enum { eeksEKIN_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksEKINF, eeksEKINO, eeksEKINSCALEF, eeksEKINSCALEH, eeksVSCALE, eeksEKINTOTAL, eeksNR };
93
94 const char *eeks_names[eeksNR]=
95 {
96     "Ekin_n", "Ekinh", "dEkindlambda", "mv_cos",
97     "Ekinf", "Ekinh_old", "EkinScaleF_NHC", "EkinScaleH_NHC","Vscale_NHC","Ekin_Total"
98 };
99
100 enum { eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
101        eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
102        eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM, 
103        eenhENERGY_DELTA_H_NN,
104        eenhENERGY_DELTA_H_LIST, eenhENERGY_DELTA_H_STARTTIME, 
105        eenhNR };
106
107 const char *eenh_names[eenhNR]=
108 {
109     "energy_n", "energy_aver", "energy_sum", "energy_nsum",
110     "energy_sum_sim", "energy_nsum_sim",
111     "energy_nsteps", "energy_nsteps_sim", 
112     "energy_delta_h_nn",
113     "energy_delta_h_list", "energy_delta_h_starttime"
114 };
115
116
117
118 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
119 static int
120 gmx_wintruncate(const char *filename, __int64 size)
121 {
122 #ifdef GMX_FAHCORE
123     /*we do this elsewhere*/
124     return 0;
125 #else
126     FILE *fp;
127     int   rc;
128     
129     fp=fopen(filename,"rb+");
130     
131     if(fp==NULL)
132     {
133         return -1;
134     }
135     
136     return _chsize_s( fileno(fp), size);
137 #endif
138 }
139 #endif
140
141
142 enum { ecprREAL, ecprRVEC, ecprMATRIX };
143
144 static const char *st_names(int cptp,int ecpt)
145 {
146     switch (cptp)
147     {
148     case 0: return est_names [ecpt]; break;
149     case 1: return eeks_names[ecpt]; break;
150     case 2: return eenh_names[ecpt]; break;
151     }
152
153     return NULL;
154 }
155
156 static void cp_warning(FILE *fp)
157 {
158     fprintf(fp,"\nWARNING: Checkpoint file is corrupted or truncated\n\n");
159 }
160
161 static void cp_error()
162 {
163     gmx_fatal(FARGS,"Checkpoint file corrupted/truncated, or maybe you are out of quota?");
164 }
165
166 static void do_cpt_string_err(XDR *xd,gmx_bool bRead,const char *desc,char **s,FILE *list)
167 {
168 #define CPTSTRLEN 1024
169     bool_t res=0;
170     
171     if (bRead)
172     {
173         snew(*s,CPTSTRLEN);
174     }
175     res = xdr_string(xd,s,CPTSTRLEN);
176     if (res == 0)
177     {
178         cp_error();
179     }
180     if (list)
181     {
182         fprintf(list,"%s = %s\n",desc,*s);
183         sfree(*s);
184     }
185 }
186
187 static int do_cpt_int(XDR *xd,const char *desc,int *i,FILE *list)
188 {
189     bool_t res=0;
190     
191     res = xdr_int(xd,i);
192     if (res == 0)
193     {
194         return -1;
195     }
196     if (list)
197     {
198         fprintf(list,"%s = %d\n",desc,*i);
199     }
200     return 0;
201 }
202
203 static int do_cpt_u_chars(XDR *xd,const char *desc,int n,unsigned char *i,FILE *list)
204 {
205     bool_t res=1;
206     int j;
207     if (list)
208     {
209         fprintf(list,"%s = ",desc);
210     }
211     for (j=0; j<n && res; j++)
212     {
213         res &= xdr_u_char(xd,&i[j]);
214         if (list)
215         {
216             fprintf(list,"%02x",i[j]);
217         }
218     }
219     if (list)
220     {
221         fprintf(list,"\n");
222     }
223     if (res == 0)
224     {
225         return -1;
226     }
227
228     return 0;
229 }
230
231 static void do_cpt_int_err(XDR *xd,const char *desc,int *i,FILE *list)
232 {
233     if (do_cpt_int(xd,desc,i,list) < 0)
234     {
235         cp_error();
236     }
237 }
238
239 static void do_cpt_step_err(XDR *xd,const char *desc,gmx_large_int_t *i,FILE *list)
240 {
241     bool_t res=0;
242     char   buf[STEPSTRSIZE];
243
244     res = xdr_gmx_large_int(xd,i,"reading checkpoint file");
245     if (res == 0)
246     {
247         cp_error();
248     }
249     if (list)
250     {
251         fprintf(list,"%s = %s\n",desc,gmx_step_str(*i,buf));
252     }
253 }
254
255 static void do_cpt_double_err(XDR *xd,const char *desc,double *f,FILE *list)
256 {
257     bool_t res=0;
258     
259     res = xdr_double(xd,f);
260     if (res == 0)
261     {
262         cp_error();
263     }
264     if (list)
265     {
266         fprintf(list,"%s = %f\n",desc,*f);
267     }
268 }
269
270 /* If nval >= 0, nval is used; on read this should match the passed value.
271  * If nval n<0, *nptr is used; on read the value is stored in nptr
272  */
273 static int do_cpte_reals_low(XDR *xd,int cptp,int ecpt,int sflags,
274                              int nval,int *nptr,real **v,
275                              FILE *list,int erealtype)
276 {
277     bool_t res=0;
278 #ifndef GMX_DOUBLE
279     int  dtc=xdr_datatype_float; 
280 #else
281     int  dtc=xdr_datatype_double;
282 #endif
283     real *vp,*va=NULL;
284     float  *vf;
285     double *vd;
286     int  nf,dt,i;
287     
288     if (list != NULL)
289     {
290         if (nval >= 0)
291         {
292             nf = nval;
293         }
294         else
295         {
296         if (nptr == NULL)
297         {
298             gmx_incons("*ntpr=NULL in do_cpte_reals_low");
299         }
300         nf = *nptr;
301         }
302     }
303     res = xdr_int(xd,&nf);
304     if (res == 0)
305     {
306         return -1;
307     }
308     if (list != NULL)
309     {
310         if (nval >= 0)
311         {
312             if (nf != nval)
313             {
314                 gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),nval,nf);
315             }
316         }
317         else
318         {
319             *nptr = nf;
320         }
321     }
322     dt = dtc;
323     res = xdr_int(xd,&dt);
324     if (res == 0)
325     {
326         return -1;
327     }
328     if (dt != dtc)
329     {
330         fprintf(stderr,"Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
331                 st_names(cptp,ecpt),xdr_datatype_names[dtc],
332                 xdr_datatype_names[dt]);
333     }
334     if (list || !(sflags & (1<<ecpt)))
335     {
336         snew(va,nf);
337         vp = va;
338     }
339     else
340     {
341         if (*v == NULL)
342         {
343             snew(*v,nf);
344         }
345         vp = *v;
346     }
347     if (dt == xdr_datatype_float)
348     {
349         if (dtc == xdr_datatype_float)
350         {
351             vf = (float *)vp;
352         }
353         else
354         {
355             snew(vf,nf);
356         }
357         res = xdr_vector(xd,(char *)vf,nf,
358                          (unsigned int)sizeof(float),(xdrproc_t)xdr_float);
359         if (res == 0)
360         {
361             return -1;
362         }
363         if (dtc != xdr_datatype_float)
364         {
365             for(i=0; i<nf; i++)
366             {
367                 vp[i] = vf[i];
368             }
369             sfree(vf);
370         }
371     }
372     else
373     {
374         if (dtc == xdr_datatype_double)
375         {
376             vd = (double *)vp;
377         }
378         else
379         {
380             snew(vd,nf);
381         }
382         res = xdr_vector(xd,(char *)vd,nf,
383                          (unsigned int)sizeof(double),(xdrproc_t)xdr_double);
384         if (res == 0)
385         {
386             return -1;
387         }
388         if (dtc != xdr_datatype_double)
389         {
390             for(i=0; i<nf; i++)
391             {
392                 vp[i] = vd[i];
393             }
394             sfree(vd);
395         }
396     }
397     
398     if (list)
399     {
400         switch (erealtype)
401         {
402         case ecprREAL:
403             pr_reals(list,0,st_names(cptp,ecpt),vp,nf);
404             break;
405         case ecprRVEC:
406             pr_rvecs(list,0,st_names(cptp,ecpt),(rvec *)vp,nf/3);
407             break;
408         default:
409             gmx_incons("Unknown checkpoint real type");
410         }
411     }
412     if (va)
413     {
414         sfree(va);
415     }
416
417     return 0;
418 }
419
420
421 /* This function stores n along with the reals for reading,
422  * but on reading it assumes that n matches the value in the checkpoint file,
423  * a fatal error is generated when this is not the case.
424  */
425 static int do_cpte_reals(XDR *xd,int cptp,int ecpt,int sflags,
426                          int n,real **v,FILE *list)
427 {
428     return do_cpte_reals_low(xd,cptp,ecpt,sflags,n,NULL,v,list,ecprREAL);
429 }
430
431 /* This function does the same as do_cpte_reals,
432  * except that on reading it ignores the passed value of *n
433  * and stored the value read from the checkpoint file in *n.
434  */
435 static int do_cpte_n_reals(XDR *xd,int cptp,int ecpt,int sflags,
436                            int *n,real **v,FILE *list)
437 {
438     return do_cpte_reals_low(xd,cptp,ecpt,sflags,-1,n,v,list,ecprREAL);
439 }
440
441 static int do_cpte_real(XDR *xd,int cptp,int ecpt,int sflags,
442                         real *r,FILE *list)
443 {
444     int n;
445
446     return do_cpte_reals_low(xd,cptp,ecpt,sflags,1,NULL,&r,list,ecprREAL);
447 }
448
449 static int do_cpte_ints(XDR *xd,int cptp,int ecpt,int sflags,
450                         int n,int **v,FILE *list)
451 {
452     bool_t res=0;
453     int  dtc=xdr_datatype_int;
454     int *vp,*va=NULL;
455     int  nf,dt,i;
456     
457     nf = n;
458     res = xdr_int(xd,&nf);
459     if (res == 0)
460     {
461         return -1;
462     }
463     if (list == NULL && v != NULL && nf != n)
464     {
465         gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),n,nf);
466     }
467     dt = dtc;
468     res = xdr_int(xd,&dt);
469     if (res == 0)
470     {
471         return -1;
472     }
473     if (dt != dtc)
474     {
475         gmx_fatal(FARGS,"Type mismatch for state entry %s, code type is %s, file type is %s\n",
476                   st_names(cptp,ecpt),xdr_datatype_names[dtc],
477                   xdr_datatype_names[dt]);
478     }
479     if (list || !(sflags & (1<<ecpt)) || v == NULL)
480     {
481         snew(va,nf);
482         vp = va;
483     }
484     else
485     {
486         if (*v == NULL)
487         {
488             snew(*v,nf);
489         }
490         vp = *v;
491     }
492     res = xdr_vector(xd,(char *)vp,nf,
493                      (unsigned int)sizeof(int),(xdrproc_t)xdr_int);
494     if (res == 0)
495     {
496         return -1;
497     }
498     if (list)
499     {
500         pr_ivec(list,0,st_names(cptp,ecpt),vp,nf,TRUE);
501     }
502     if (va)
503     {
504         sfree(va);
505     }
506
507     return 0;
508 }
509
510 static int do_cpte_int(XDR *xd,int cptp,int ecpt,int sflags,
511                        int *i,FILE *list)
512 {
513     return do_cpte_ints(xd,cptp,ecpt,sflags,1,&i,list);
514 }
515
516 static int do_cpte_doubles(XDR *xd,int cptp,int ecpt,int sflags,
517                            int n,double **v,FILE *list)
518 {
519     bool_t res=0;
520     int  dtc=xdr_datatype_double;
521     double *vp,*va=NULL;
522     int  nf,dt,i;
523     
524     nf = n;
525     res = xdr_int(xd,&nf);
526     if (res == 0)
527     {
528         return -1;
529     }
530     if (list == NULL && nf != n)
531     {
532         gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),n,nf);
533     }
534     dt = dtc;
535     res = xdr_int(xd,&dt);
536     if (res == 0)
537     {
538         return -1;
539     }
540     if (dt != dtc)
541     {
542         gmx_fatal(FARGS,"Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
543                   st_names(cptp,ecpt),xdr_datatype_names[dtc],
544                   xdr_datatype_names[dt]);
545     }
546     if (list || !(sflags & (1<<ecpt)))
547     {
548         snew(va,nf);
549         vp = va;
550     }
551     else
552     {
553         if (*v == NULL)
554         {
555             snew(*v,nf);
556         }
557         vp = *v;
558     }
559     res = xdr_vector(xd,(char *)vp,nf,
560                      (unsigned int)sizeof(double),(xdrproc_t)xdr_double);
561     if (res == 0)
562     {
563         return -1;
564     }
565     if (list)
566     {
567         pr_doubles(list,0,st_names(cptp,ecpt),vp,nf);
568     }
569     if (va)
570     {
571         sfree(va);
572     }
573
574     return 0;
575 }
576
577
578 static int do_cpte_rvecs(XDR *xd,int cptp,int ecpt,int sflags,
579                          int n,rvec **v,FILE *list)
580 {
581     int n3;
582
583     return do_cpte_reals_low(xd,cptp,ecpt,sflags,
584                              n*DIM,NULL,(real **)v,list,ecprRVEC);
585 }
586
587 static int do_cpte_matrix(XDR *xd,int cptp,int ecpt,int sflags,
588                           matrix v,FILE *list)
589 {
590     real *vr;
591     real ret;
592
593     vr = (real *)&(v[0][0]);
594     ret = do_cpte_reals_low(xd,cptp,ecpt,sflags,
595                             DIM*DIM,NULL,&vr,NULL,ecprMATRIX);
596     
597     if (list && ret == 0)
598     {
599         pr_rvecs(list,0,st_names(cptp,ecpt),v,DIM);
600     }
601     
602     return ret;
603 }
604
605 static int do_cpte_matrices(XDR *xd,int cptp,int ecpt,int sflags,
606                             int n,matrix **v,FILE *list)
607 {
608     bool_t res=0;
609     matrix *vp,*va=NULL;
610     real *vr;
611     int  nf,i,j,k;
612     int  ret;
613
614     nf = n;
615     res = xdr_int(xd,&nf);
616     if (res == 0)
617     {
618         return -1;
619     }
620     if (list == NULL && nf != n)
621     {
622         gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),n,nf);
623     }
624     if (list || !(sflags & (1<<ecpt)))
625     {
626         snew(va,nf);
627         vp = va;
628     }
629     else
630     {
631         if (*v == NULL)
632         {
633             snew(*v,nf);
634         }
635         vp = *v;
636     }
637     snew(vr,nf*DIM*DIM);
638     for(i=0; i<nf; i++)
639     {
640         for(j=0; j<DIM; j++)
641         {
642             for(k=0; k<DIM; k++)
643             {
644                 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
645             }
646         }
647     }
648     ret = do_cpte_reals_low(xd,cptp,ecpt,sflags,
649                             nf*DIM*DIM,NULL,&vr,NULL,ecprMATRIX);
650     for(i=0; i<nf; i++)
651     {
652         for(j=0; j<DIM; j++)
653         {
654             for(k=0; k<DIM; k++)
655             {
656                 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
657             }
658         }
659     }
660     sfree(vr);
661     
662     if (list && ret == 0)
663     {
664         for(i=0; i<nf; i++)
665         {
666             pr_rvecs(list,0,st_names(cptp,ecpt),vp[i],DIM);
667         }
668     }
669     if (va)
670     {
671         sfree(va);
672     }
673     
674     return ret;
675 }
676
677 static void do_cpt_header(XDR *xd,gmx_bool bRead,int *file_version,
678                           char **version,char **btime,char **buser,char **bmach,
679                           char **fprog,char **ftime,
680                           int *eIntegrator,int *simulation_part,
681                           gmx_large_int_t *step,double *t,
682                           int *nnodes,int *dd_nc,int *npme,
683                           int *natoms,int *ngtc, int *nnhpres, int *nhchainlength,
684                           int *flags_state,int *flags_eks,int *flags_enh,
685                           FILE *list)
686 {
687     bool_t res=0;
688     int  magic;
689     int  idum=0;
690     int  i;
691     char *fhost;
692
693     if (bRead)
694     {
695         magic = -1;
696     }
697     else
698     {
699         magic = CPT_MAGIC1;
700     }
701     res = xdr_int(xd,&magic);
702     if (res == 0)
703     {
704         gmx_fatal(FARGS,"The checkpoint file is empty/corrupted, or maybe you are out of quota?");
705     }
706     if (magic != CPT_MAGIC1)
707     {
708         gmx_fatal(FARGS,"Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
709                   "The checkpoint file is corrupted or not a checkpoint file",
710                   magic,CPT_MAGIC1);
711     }
712     if (!bRead)
713     {
714         snew(fhost,255);
715 #ifdef HAVE_UNISTD_H
716         if (gethostname(fhost,255) != 0)
717         {
718             sprintf(fhost,"unknown");
719         }
720 #else
721         sprintf(fhost,"unknown");
722 #endif  
723     }
724     do_cpt_string_err(xd,bRead,"GROMACS version"           ,version,list);
725     do_cpt_string_err(xd,bRead,"GROMACS build time"        ,btime,list);
726     do_cpt_string_err(xd,bRead,"GROMACS build user"        ,buser,list);
727     do_cpt_string_err(xd,bRead,"GROMACS build machine"     ,bmach,list);
728     do_cpt_string_err(xd,bRead,"generating program"        ,fprog,list);
729     do_cpt_string_err(xd,bRead,"generation time"           ,ftime,list);
730     *file_version = cpt_version;
731     do_cpt_int_err(xd,"checkpoint file version",file_version,list);
732     if (*file_version > cpt_version)
733     {
734         gmx_fatal(FARGS,"Attempting to read a checkpoint file of version %d with code of version %d\n",*file_version,cpt_version);
735     }
736     if (*file_version >= 12)
737     {
738         do_cpt_string_err(xd,bRead,"generating host"           ,&fhost,list);
739         if (list == NULL)
740         {
741             sfree(fhost);
742         }
743     }
744     do_cpt_int_err(xd,"#atoms"            ,natoms     ,list);
745     do_cpt_int_err(xd,"#T-coupling groups",ngtc       ,list);
746     if (*file_version >= 10) 
747     {
748         do_cpt_int_err(xd,"#Nose-Hoover T-chains",nhchainlength,list);
749     }
750     else
751     {
752         *nhchainlength = 1;
753     }
754     if (*file_version >= 11)
755     {
756         do_cpt_int_err(xd,"#Nose-Hoover T-chains for barostat ",nnhpres,list);
757     }
758     else
759     {
760         *nnhpres = 0;
761     }
762     do_cpt_int_err(xd,"integrator"        ,eIntegrator,list);
763         if (*file_version >= 3)
764         {
765                 do_cpt_int_err(xd,"simulation part #", simulation_part,list);
766         }
767         else
768         {
769                 *simulation_part = 1;
770         }
771     if (*file_version >= 5)
772     {
773         do_cpt_step_err(xd,"step"         ,step       ,list);
774     }
775     else
776     {
777         do_cpt_int_err(xd,"step"          ,&idum      ,list);
778         *step = idum;
779     }
780     do_cpt_double_err(xd,"t"              ,t          ,list);
781     do_cpt_int_err(xd,"#PP-nodes"         ,nnodes     ,list);
782     idum = 1;
783     do_cpt_int_err(xd,"dd_nc[x]",dd_nc ? &(dd_nc[0]) : &idum,list);
784     do_cpt_int_err(xd,"dd_nc[y]",dd_nc ? &(dd_nc[1]) : &idum,list);
785     do_cpt_int_err(xd,"dd_nc[z]",dd_nc ? &(dd_nc[2]) : &idum,list);
786     do_cpt_int_err(xd,"#PME-only nodes",npme,list);
787     do_cpt_int_err(xd,"state flags",flags_state,list);
788         if (*file_version >= 4)
789     {
790         do_cpt_int_err(xd,"ekin data flags",flags_eks,list);
791         do_cpt_int_err(xd,"energy history flags",flags_enh,list);
792     }
793     else
794     {
795         *flags_eks  = 0;
796         *flags_enh   = (*flags_state >> (estORIRE_DTAV+1));
797         *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
798                                          (1<<(estORIRE_DTAV+2)) |
799                                          (1<<(estORIRE_DTAV+3))));
800     }
801 }
802
803 static int do_cpt_footer(XDR *xd,gmx_bool bRead,int file_version)
804 {
805     bool_t res=0;
806     int  magic;
807     
808     if (file_version >= 2)
809     {
810         magic = CPT_MAGIC2;
811         res = xdr_int(xd,&magic);
812         if (res == 0)
813         {
814             cp_error();
815         }
816         if (magic != CPT_MAGIC2)
817         {
818             return -1;
819         }
820     }
821
822     return 0;
823 }
824
825 static int do_cpt_state(XDR *xd,gmx_bool bRead,
826                         int fflags,t_state *state,
827                         gmx_bool bReadRNG,FILE *list)
828 {
829     int  sflags;
830     int  **rng_p,**rngi_p;
831     int  i;
832     int  ret;
833     int  nnht,nnhtp;
834
835     ret = 0;
836     
837     nnht = state->nhchainlength*state->ngtc;
838     nnhtp = state->nhchainlength*state->nnhpres;
839
840     if (bReadRNG)
841     {
842         rng_p  = (int **)&state->ld_rng;
843         rngi_p = &state->ld_rngi;
844     }
845     else
846     {
847         /* Do not read the RNG data */
848         rng_p  = NULL;
849         rngi_p = NULL;
850     }
851
852     sflags = state->flags;
853     for(i=0; (i<estNR && ret == 0); i++)
854     {
855         if (fflags & (1<<i))
856         {
857             switch (i)
858             {
859             case estLAMBDA:  ret = do_cpte_real(xd,0,i,sflags,&state->lambda,list); break;
860             case estBOX:     ret = do_cpte_matrix(xd,0,i,sflags,state->box,list); break;
861             case estBOX_REL: ret = do_cpte_matrix(xd,0,i,sflags,state->box_rel,list); break;
862             case estBOXV:    ret = do_cpte_matrix(xd,0,i,sflags,state->boxv,list); break;
863             case estPRES_PREV: ret = do_cpte_matrix(xd,0,i,sflags,state->pres_prev,list); break;
864             case estSVIR_PREV:  ret = do_cpte_matrix(xd,0,i,sflags,state->svir_prev,list); break;
865             case estFVIR_PREV:  ret = do_cpte_matrix(xd,0,i,sflags,state->fvir_prev,list); break;
866             case estNH_XI:   ret = do_cpte_doubles(xd,0,i,sflags,nnht,&state->nosehoover_xi,list); break;
867             case estNH_VXI:  ret = do_cpte_doubles(xd,0,i,sflags,nnht,&state->nosehoover_vxi,list); break;
868             case estNHPRES_XI:   ret = do_cpte_doubles(xd,0,i,sflags,nnhtp,&state->nhpres_xi,list); break;
869             case estNHPRES_VXI:  ret = do_cpte_doubles(xd,0,i,sflags,nnhtp,&state->nhpres_vxi,list); break;
870             case estTC_INT:  ret = do_cpte_doubles(xd,0,i,sflags,state->ngtc,&state->therm_integral,list); break;
871             case estVETA:    ret = do_cpte_real(xd,0,i,sflags,&state->veta,list); break;
872             case estVOL0:    ret = do_cpte_real(xd,0,i,sflags,&state->vol0,list); break;
873             case estX:       ret = do_cpte_rvecs(xd,0,i,sflags,state->natoms,&state->x,list); break;
874             case estV:       ret = do_cpte_rvecs(xd,0,i,sflags,state->natoms,&state->v,list); break;
875             case estSDX:     ret = do_cpte_rvecs(xd,0,i,sflags,state->natoms,&state->sd_X,list); break;
876             case estLD_RNG:  ret = do_cpte_ints(xd,0,i,sflags,state->nrng,rng_p,list); break;
877             case estLD_RNGI: ret = do_cpte_ints(xd,0,i,sflags,state->nrngi,rngi_p,list); break;
878             case estDISRE_INITF:  ret = do_cpte_real (xd,0,i,sflags,&state->hist.disre_initf,list); break;
879             case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd,0,i,sflags,&state->hist.ndisrepairs,&state->hist.disre_rm3tav,list); break;
880             case estORIRE_INITF:  ret = do_cpte_real (xd,0,i,sflags,&state->hist.orire_initf,list); break;
881             case estORIRE_DTAV:   ret = do_cpte_n_reals(xd,0,i,sflags,&state->hist.norire_Dtav,&state->hist.orire_Dtav,list); break;
882             default:
883                 gmx_fatal(FARGS,"Unknown state entry %d\n"
884                           "You are probably reading a new checkpoint file with old code",i);
885             }
886         }
887     }
888     
889     return ret;
890 }
891
892 static int do_cpt_ekinstate(XDR *xd,gmx_bool bRead,
893                             int fflags,ekinstate_t *ekins,
894                             FILE *list)
895 {
896     int  i;
897     int  ret;
898
899     ret = 0;
900
901     for(i=0; (i<eeksNR && ret == 0); i++)
902     {
903         if (fflags & (1<<i))
904         {
905             switch (i)
906             {
907                 
908                         case eeksEKIN_N:     ret = do_cpte_int(xd,1,i,fflags,&ekins->ekin_n,list); break;
909                         case eeksEKINH :     ret = do_cpte_matrices(xd,1,i,fflags,ekins->ekin_n,&ekins->ekinh,list); break;
910                         case eeksEKINF:      ret = do_cpte_matrices(xd,1,i,fflags,ekins->ekin_n,&ekins->ekinf,list); break;
911                         case eeksEKINO:      ret = do_cpte_matrices(xd,1,i,fflags,ekins->ekin_n,&ekins->ekinh_old,list); break;
912             case eeksEKINTOTAL:  ret = do_cpte_matrix(xd,1,i,fflags,ekins->ekin_total,list); break;
913             case eeksEKINSCALEF: ret = do_cpte_doubles(xd,1,i,fflags,ekins->ekin_n,&ekins->ekinscalef_nhc,list); break;
914             case eeksVSCALE:     ret = do_cpte_doubles(xd,1,i,fflags,ekins->ekin_n,&ekins->vscale_nhc,list); break;
915             case eeksEKINSCALEH: ret = do_cpte_doubles(xd,1,i,fflags,ekins->ekin_n,&ekins->ekinscaleh_nhc,list); break;
916                         case eeksDEKINDL :   ret = do_cpte_real(xd,1,i,fflags,&ekins->dekindl,list); break;
917             case eeksMVCOS:      ret = do_cpte_real(xd,1,i,fflags,&ekins->mvcos,list); break;                   
918             default:
919                 gmx_fatal(FARGS,"Unknown ekin data state entry %d\n"
920                           "You are probably reading a new checkpoint file with old code",i);
921             }
922         }
923     }
924     
925     return ret;
926 }
927
928
929 static int do_cpt_enerhist(XDR *xd,gmx_bool bRead,
930                            int fflags,energyhistory_t *enerhist,
931                            FILE *list)
932 {
933     int  i;
934     int  j;
935     int  ret;
936
937     ret = 0;
938
939     if (bRead)
940     {
941         enerhist->nsteps     = 0;
942         enerhist->nsum       = 0;
943         enerhist->nsteps_sim = 0;
944         enerhist->nsum_sim   = 0;
945         enerhist->dht        = NULL;
946
947         if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
948         {
949             snew(enerhist->dht,1);
950             enerhist->dht->ndh = NULL;
951             enerhist->dht->dh = NULL;
952         }
953     }
954
955     for(i=0; (i<eenhNR && ret == 0); i++)
956     {
957         if (fflags & (1<<i))
958         {
959             switch (i)
960             {
961                 case eenhENERGY_N:     ret = do_cpte_int(xd,2,i,fflags,&enerhist->nener,list); break;
962                 case eenhENERGY_AVER:  ret = do_cpte_doubles(xd,2,i,fflags,enerhist->nener,&enerhist->ener_ave,list); break;
963                 case eenhENERGY_SUM:   ret = do_cpte_doubles(xd,2,i,fflags,enerhist->nener,&enerhist->ener_sum,list); break;
964                 case eenhENERGY_NSUM:  do_cpt_step_err(xd,eenh_names[i],&enerhist->nsum,list); break;
965                 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd,2,i,fflags,enerhist->nener,&enerhist->ener_sum_sim,list); break;
966                 case eenhENERGY_NSUM_SIM:   do_cpt_step_err(xd,eenh_names[i],&enerhist->nsum_sim,list); break;
967                 case eenhENERGY_NSTEPS:     do_cpt_step_err(xd,eenh_names[i],&enerhist->nsteps,list); break;
968                 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd,eenh_names[i],&enerhist->nsteps_sim,list); break;
969                 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd,eenh_names[i], &(enerhist->dht->nndh), list); 
970                     if (bRead) /* now allocate memory for it */
971                     {
972                         snew(enerhist->dht->dh, enerhist->dht->nndh);
973                         snew(enerhist->dht->ndh, enerhist->dht->nndh);
974                         for(j=0;j<enerhist->dht->nndh;j++)
975                         {
976                             enerhist->dht->ndh[j] = 0;
977                             enerhist->dht->dh[j] = NULL;
978                         }
979                     }
980                 break;
981                 case eenhENERGY_DELTA_H_LIST: 
982                     for(j=0;j<enerhist->dht->nndh;j++)
983                     {
984                         ret=do_cpte_n_reals(xd, 2, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list); 
985                     }
986                     break;
987                 case eenhENERGY_DELTA_H_STARTTIME: 
988                     ret=do_cpte_real(xd, 2, i, fflags, &(enerhist->dht->starttime), list); break;
989                 default:
990                     gmx_fatal(FARGS,"Unknown energy history entry %d\n"
991                               "You are probably reading a new checkpoint file with old code",i);
992             }
993         }
994     }
995
996     if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
997     {
998         /* Assume we have an old file format and copy sum to sum_sim */
999         srenew(enerhist->ener_sum_sim,enerhist->nener);
1000         for(i=0; i<enerhist->nener; i++)
1001         {
1002             enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1003         }
1004         fflags |= (1<<eenhENERGY_SUM_SIM);
1005     }
1006     
1007     if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1008         !(fflags & (1<<eenhENERGY_NSTEPS)))
1009     {
1010         /* Assume we have an old file format and copy nsum to nsteps */
1011         enerhist->nsteps = enerhist->nsum;
1012         fflags |= (1<<eenhENERGY_NSTEPS);
1013     }
1014     if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1015         !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1016     {
1017         /* Assume we have an old file format and copy nsum to nsteps */
1018         enerhist->nsteps_sim = enerhist->nsum_sim;
1019         fflags |= (1<<eenhENERGY_NSTEPS_SIM);
1020     }
1021
1022     return ret;
1023 }
1024
1025 static int do_cpt_files(XDR *xd, gmx_bool bRead, 
1026                         gmx_file_position_t **p_outputfiles, int *nfiles, 
1027                         FILE *list, int file_version)
1028 {
1029     int    i,j;
1030     gmx_off_t  offset;
1031     gmx_off_t  mask = 0xFFFFFFFFL;
1032     int    offset_high,offset_low;
1033     char   *buf;
1034     gmx_file_position_t *outputfiles;
1035
1036     if (do_cpt_int(xd,"number of output files",nfiles,list) != 0)
1037     {
1038         return -1;
1039     }
1040
1041     if(bRead)
1042     {
1043         snew(*p_outputfiles,*nfiles);
1044     }
1045
1046     outputfiles = *p_outputfiles;
1047
1048     for(i=0;i<*nfiles;i++)
1049     {
1050         /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1051         if(bRead)
1052         {
1053             do_cpt_string_err(xd,bRead,"output filename",&buf,list);
1054             strncpy(outputfiles[i].filename,buf,CPTSTRLEN-1);
1055             if(list==NULL)
1056             {
1057                 sfree(buf);                     
1058             }
1059
1060             if (do_cpt_int(xd,"file_offset_high",&offset_high,list) != 0)
1061             {
1062                 return -1;
1063             }
1064             if (do_cpt_int(xd,"file_offset_low",&offset_low,list) != 0)
1065             {
1066                 return -1;
1067             }
1068 #if (SIZEOF_GMX_OFF_T > 4)
1069             outputfiles[i].offset = ( ((gmx_off_t) offset_high) << 32 ) | ( (gmx_off_t) offset_low & mask );
1070 #else
1071             outputfiles[i].offset = offset_low;
1072 #endif
1073         }
1074         else
1075         {
1076             buf = outputfiles[i].filename;
1077             do_cpt_string_err(xd,bRead,"output filename",&buf,list);
1078             /* writing */
1079             offset      = outputfiles[i].offset;
1080             if (offset == -1)
1081             {
1082                 offset_low  = -1;
1083                 offset_high = -1;
1084             }
1085             else
1086             {
1087 #if (SIZEOF_GMX_OFF_T > 4)
1088                 offset_low  = (int) (offset & mask);
1089                 offset_high = (int) ((offset >> 32) & mask);
1090 #else
1091                 offset_low  = offset;
1092                 offset_high = 0;
1093 #endif
1094             }
1095             if (do_cpt_int(xd,"file_offset_high",&offset_high,list) != 0)
1096             {
1097                 return -1;
1098             }
1099             if (do_cpt_int(xd,"file_offset_low",&offset_low,list) != 0)
1100             {
1101                 return -1;
1102             }
1103         }
1104         if (file_version >= 8)
1105         {
1106             if (do_cpt_int(xd,"file_checksum_size",&(outputfiles[i].chksum_size),
1107                            list) != 0)
1108             {
1109                 return -1;
1110             }
1111             if (do_cpt_u_chars(xd,"file_checksum",16,outputfiles[i].chksum,list) != 0)
1112             {
1113                 return -1;
1114             }
1115         } 
1116         else 
1117         {
1118             outputfiles[i].chksum_size = -1;
1119         }
1120     }
1121     return 0;
1122 }
1123
1124
1125 void write_checkpoint(const char *fn,gmx_bool bNumberAndKeep,
1126                       FILE *fplog,t_commrec *cr,
1127                       int eIntegrator,int simulation_part,
1128                       gmx_large_int_t step,double t,t_state *state)
1129 {
1130     t_fileio *fp;
1131     int  file_version;
1132     char *version;
1133     char *btime;
1134     char *buser;
1135     char *bmach;
1136     char *fprog;
1137     char *fntemp; /* the temporary checkpoint file name */
1138     time_t now;
1139     char timebuf[STRLEN];
1140     int  nppnodes,npmenodes,flag_64bit;
1141     char buf[1024],suffix[5+STEPSTRSIZE],sbuf[STEPSTRSIZE];
1142     gmx_file_position_t *outputfiles;
1143     int  noutputfiles;
1144     char *ftime;
1145     int  flags_eks,flags_enh,i;
1146     t_fileio *ret;
1147                 
1148     if (PAR(cr))
1149     {
1150         if (DOMAINDECOMP(cr))
1151         {
1152             nppnodes  = cr->dd->nnodes;
1153             npmenodes = cr->npmenodes;
1154         }
1155         else
1156         {
1157             nppnodes  = cr->nnodes;
1158             npmenodes = 0;
1159         }
1160     }
1161     else
1162     {
1163         nppnodes  = 1;
1164         npmenodes = 0;
1165     }
1166
1167     /* make the new temporary filename */
1168     snew(fntemp, strlen(fn)+5+STEPSTRSIZE);
1169     strcpy(fntemp,fn);
1170     fntemp[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1171     sprintf(suffix,"_%s%s","step",gmx_step_str(step,sbuf));
1172     strcat(fntemp,suffix);
1173     strcat(fntemp,fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1174    
1175     time(&now);
1176     gmx_ctime_r(&now,timebuf,STRLEN);
1177
1178     if (fplog)
1179     { 
1180         fprintf(fplog,"Writing checkpoint, step %s at %s\n\n",
1181                 gmx_step_str(step,buf),timebuf);
1182     }
1183     
1184     /* Get offsets for open files */
1185     gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1186
1187     fp = gmx_fio_open(fntemp,"w");
1188         
1189     if (state->ekinstate.bUpToDate)
1190     {
1191         flags_eks =
1192             ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) | 
1193              (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) | 
1194              (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1195     }
1196     else
1197     {
1198         flags_eks = 0;
1199     }
1200
1201     flags_enh = 0;
1202     if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1203     {
1204         flags_enh |= (1<<eenhENERGY_N);
1205         if (state->enerhist.nsum > 0)
1206         {
1207             flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1208                           (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSUM));
1209         }
1210         if (state->enerhist.nsum_sim > 0)
1211         {
1212             flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSTEPS_SIM) |
1213                           (1<<eenhENERGY_NSUM_SIM));
1214         }
1215         if (state->enerhist.dht)
1216         {
1217             flags_enh |= ( (1<< eenhENERGY_DELTA_H_NN) |
1218                            (1<< eenhENERGY_DELTA_H_LIST) | 
1219                            (1<< eenhENERGY_DELTA_H_STARTTIME) );
1220         }
1221     }
1222
1223     
1224     version = strdup(VERSION);
1225     btime   = strdup(BUILD_TIME);
1226     buser   = strdup(BUILD_USER);
1227     bmach   = strdup(BUILD_MACHINE);
1228     fprog   = strdup(Program());
1229
1230     ftime   = &(timebuf[0]);
1231     
1232     do_cpt_header(gmx_fio_getxdr(fp),FALSE,&file_version,
1233                   &version,&btime,&buser,&bmach,&fprog,&ftime,
1234                   &eIntegrator,&simulation_part,&step,&t,&nppnodes,
1235                   DOMAINDECOMP(cr) ? cr->dd->nc : NULL,&npmenodes,
1236                   &state->natoms,&state->ngtc,&state->nnhpres,
1237                   &state->nhchainlength, &state->flags,&flags_eks,&flags_enh,
1238                   NULL);
1239     
1240     sfree(version);
1241     sfree(btime);
1242     sfree(buser);
1243     sfree(bmach);
1244     sfree(fprog);
1245
1246     if((do_cpt_state(gmx_fio_getxdr(fp),FALSE,state->flags,state,TRUE,NULL) < 0)        ||
1247        (do_cpt_ekinstate(gmx_fio_getxdr(fp),FALSE,flags_eks,&state->ekinstate,NULL) < 0)||
1248        (do_cpt_enerhist(gmx_fio_getxdr(fp),FALSE,flags_enh,&state->enerhist,NULL) < 0)  ||
1249        (do_cpt_files(gmx_fio_getxdr(fp),FALSE,&outputfiles,&noutputfiles,NULL,
1250                      file_version) < 0))
1251     {
1252         gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1253     }
1254
1255     do_cpt_footer(gmx_fio_getxdr(fp),FALSE,file_version);
1256
1257     /* we really, REALLY, want to make sure to physically write the checkpoint, 
1258        and all the files it depends on, out to disk. Because we've
1259        opened the checkpoint with gmx_fio_open(), it's in our list
1260        of open files.  */
1261     ret=gmx_fio_all_output_fsync();
1262
1263     if (ret)
1264     {
1265         char buf[STRLEN];
1266         sprintf(buf,
1267                 "Cannot fsync '%s'; maybe you are out of disk space or quota?",
1268                 gmx_fio_getname(ret));
1269
1270         if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV)==NULL)
1271         {
1272             gmx_file(buf);
1273         }
1274         else
1275         {
1276             gmx_warning(buf);
1277         }
1278     }
1279
1280     if( gmx_fio_close(fp) != 0)
1281     {
1282         gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1283     }
1284
1285     /* we don't move the checkpoint if the user specified they didn't want it,
1286        or if the fsyncs failed */
1287     if (!bNumberAndKeep && !ret)
1288     {
1289         if (gmx_fexist(fn))
1290         {
1291             /* Rename the previous checkpoint file */
1292             strcpy(buf,fn);
1293             buf[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1294             strcat(buf,"_prev");
1295             strcat(buf,fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1296 #ifndef GMX_FAHCORE
1297             /* we copy here so that if something goes wrong between now and
1298              * the rename below, there's always a state.cpt.
1299              * If renames are atomic (such as in POSIX systems),
1300              * this copying should be unneccesary.
1301              */
1302             gmx_file_copy(fn, buf, FALSE);
1303             /* We don't really care if this fails: 
1304              * there's already a new checkpoint.
1305              */
1306 #else
1307             gmx_file_rename(fn, buf);
1308 #endif
1309         }
1310         if (gmx_file_rename(fntemp, fn) != 0)
1311         {
1312             gmx_file("Cannot rename checkpoint file; maybe you are out of quota?");
1313         }
1314     }
1315
1316     sfree(outputfiles);
1317     sfree(fntemp);
1318
1319 #ifdef GMX_FAHCORE
1320     /*code for alternate checkpointing scheme.  moved from top of loop over 
1321       steps */
1322     fcRequestCheckPoint();
1323     if ( fcCheckPointParallel( cr->nodeid, NULL,0) == 0 ) {
1324         gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", step );
1325     }
1326 #endif /* end GMX_FAHCORE block */
1327 }
1328
1329 static void print_flag_mismatch(FILE *fplog,int sflags,int fflags)
1330 {
1331     int i;
1332     
1333     fprintf(fplog,"\nState entry mismatch between the simulation and the checkpoint file\n");
1334     fprintf(fplog,"Entries which are not present in the checkpoint file will not be updated\n");
1335     fprintf(fplog,"  %24s    %11s    %11s\n","","simulation","checkpoint");
1336     for(i=0; i<estNR; i++)
1337     {
1338         if ((sflags & (1<<i)) || (fflags & (1<<i)))
1339         {
1340             fprintf(fplog,"  %24s    %11s    %11s\n",
1341                     est_names[i],
1342                     (sflags & (1<<i)) ? "  present  " : "not present",
1343                     (fflags & (1<<i)) ? "  present  " : "not present");
1344         }
1345     }
1346 }
1347
1348 static void check_int(FILE *fplog,const char *type,int p,int f,gmx_bool *mm)
1349 {
1350         FILE *fp = fplog ? fplog : stderr;
1351
1352     if (p != f)
1353     {
1354                 fprintf(fp,"  %s mismatch,\n",type);
1355                 fprintf(fp,"    current program: %d\n",p);
1356                 fprintf(fp,"    checkpoint file: %d\n",f);
1357                 fprintf(fp,"\n");
1358         *mm = TRUE;
1359     }
1360 }
1361
1362 static void check_string(FILE *fplog,const char *type,const char *p,
1363                          const char *f,gmx_bool *mm)
1364 {
1365         FILE *fp = fplog ? fplog : stderr;
1366         
1367     if (strcmp(p,f) != 0)
1368     {
1369                 fprintf(fp,"  %s mismatch,\n",type);
1370                 fprintf(fp,"    current program: %s\n",p);
1371                 fprintf(fp,"    checkpoint file: %s\n",f);
1372                 fprintf(fp,"\n");
1373         *mm = TRUE;
1374     }
1375 }
1376
1377 static void check_match(FILE *fplog,
1378                         char *version,
1379                         char *btime,char *buser,char *bmach,char *fprog,
1380                         t_commrec *cr,gmx_bool bPartDecomp,int npp_f,int npme_f,
1381                         ivec dd_nc,ivec dd_nc_f)
1382 {
1383     int  npp;
1384     gmx_bool mm;
1385     
1386     mm = FALSE;
1387     
1388     check_string(fplog,"Version"      ,VERSION      ,version,&mm);
1389     check_string(fplog,"Build time"   ,BUILD_TIME   ,btime  ,&mm);
1390     check_string(fplog,"Build user"   ,BUILD_USER   ,buser  ,&mm);
1391     check_string(fplog,"Build machine",BUILD_MACHINE,bmach  ,&mm);
1392     check_string(fplog,"Program name" ,Program()    ,fprog  ,&mm);
1393     
1394     npp = cr->nnodes - cr->npmenodes;
1395     check_int   (fplog,"#nodes"       ,cr->nnodes   ,npp_f+npme_f ,&mm);
1396     if (bPartDecomp)
1397     {
1398         dd_nc[XX] = 1;
1399         dd_nc[YY] = 1;
1400         dd_nc[ZZ] = 1;
1401     }
1402     if (npp > 1)
1403     {
1404         check_int (fplog,"#PME-nodes"  ,cr->npmenodes,npme_f     ,&mm);
1405         if (npp == npp_f)
1406         {
1407             check_int (fplog,"#DD-cells[x]",dd_nc[XX]    ,dd_nc_f[XX],&mm);
1408             check_int (fplog,"#DD-cells[y]",dd_nc[YY]    ,dd_nc_f[YY],&mm);
1409             check_int (fplog,"#DD-cells[z]",dd_nc[ZZ]    ,dd_nc_f[ZZ],&mm);
1410         }
1411     }
1412     
1413     if (mm)
1414     {
1415                 fprintf(stderr,
1416                                 "Gromacs binary or parallel settings not identical to previous run.\n"
1417                                 "Continuation is exact, but is not guaranteed to be binary identical%s.\n\n",
1418                                 fplog ? ",\n see the log file for details" : "");
1419                 
1420         if (fplog)
1421         {
1422                         fprintf(fplog,
1423                                         "Gromacs binary or parallel settings not identical to previous run.\n"
1424                                         "Continuation is exact, but is not guaranteed to be binary identical.\n\n");
1425                 }
1426     }
1427 }
1428
1429 static void read_checkpoint(const char *fn,FILE **pfplog,
1430                             t_commrec *cr,gmx_bool bPartDecomp,ivec dd_nc,
1431                             int eIntegrator,gmx_large_int_t *step,double *t,
1432                             t_state *state,gmx_bool *bReadRNG,gmx_bool *bReadEkin,
1433                             int *simulation_part,gmx_bool bAppendOutputFiles)
1434 {
1435     t_fileio *fp;
1436     int  i,j,rc;
1437     int  file_version;
1438     char *version,*btime,*buser,*bmach,*fprog,*ftime;
1439         char filename[STRLEN],buf[STEPSTRSIZE];
1440     int  nppnodes,eIntegrator_f,nppnodes_f,npmenodes_f;
1441     ivec dd_nc_f;
1442     int  natoms,ngtc,nnhpres,nhchainlength,fflags,flags_eks,flags_enh;
1443     int  d;
1444     int  ret;
1445         gmx_file_position_t *outputfiles;
1446         int  nfiles;
1447         t_fileio *chksum_file;
1448         FILE* fplog = *pfplog;
1449         unsigned char digest[16];
1450 #if !((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
1451         struct flock fl = { 0, SEEK_SET, 0,       F_WRLCK,     0 }; 
1452 #endif
1453         
1454     const char *int_warn=
1455         "WARNING: The checkpoint file was generator with integrator %s,\n"
1456         "         while the simulation uses integrator %s\n\n";
1457     const char *sd_note=
1458         "NOTE: The checkpoint file was for %d nodes doing SD or BD,\n"
1459         "      while the simulation uses %d SD or BD nodes,\n"
1460         "      continuation will be exact, except for the random state\n\n";
1461     
1462     if (PARTDECOMP(cr))
1463     {
1464         gmx_fatal(FARGS,
1465                   "read_checkpoint not (yet) supported with particle decomposition");
1466     }
1467     
1468     fp = gmx_fio_open(fn,"r");
1469     do_cpt_header(gmx_fio_getxdr(fp),TRUE,&file_version,
1470                   &version,&btime,&buser,&bmach,&fprog,&ftime,
1471                   &eIntegrator_f,simulation_part,step,t,
1472                   &nppnodes_f,dd_nc_f,&npmenodes_f,
1473                   &natoms,&ngtc,&nnhpres,&nhchainlength,
1474                   &fflags,&flags_eks,&flags_enh,NULL);
1475     
1476     if (cr == NULL || MASTER(cr))
1477     {
1478         fprintf(stderr,"\nReading checkpoint file %s generated: %s\n\n",
1479                 fn,ftime);
1480     }
1481         
1482         /* This will not be written if we do appending, since fplog is still NULL then */
1483     if (fplog)
1484     {
1485         fprintf(fplog,"\n");
1486         fprintf(fplog,"Reading checkpoint file %s\n",fn);
1487         fprintf(fplog,"  file generated by:     %s\n",fprog);  
1488         fprintf(fplog,"  file generated at:     %s\n",ftime);  
1489         fprintf(fplog,"  GROMACS build time:    %s\n",btime);  
1490         fprintf(fplog,"  GROMACS build user:    %s\n",buser);  
1491         fprintf(fplog,"  GROMACS build machine: %s\n",bmach);  
1492         fprintf(fplog,"  simulation part #:     %d\n",*simulation_part);
1493         fprintf(fplog,"  step:                  %s\n",gmx_step_str(*step,buf));
1494         fprintf(fplog,"  time:                  %f\n",*t);  
1495         fprintf(fplog,"\n");
1496     }
1497     
1498     if (natoms != state->natoms)
1499     {
1500         gmx_fatal(FARGS,"Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms",natoms,state->natoms);
1501     }
1502     if (ngtc != state->ngtc)
1503     {
1504         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);
1505     }
1506     if (nnhpres != state->nnhpres)
1507     {
1508         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);
1509     }
1510
1511     init_gtc_state(state,state->ngtc,state->nnhpres,nhchainlength); /* need to keep this here to keep the tpr format working */
1512     /* write over whatever was read; we use the number of Nose-Hoover chains from the checkpoint */
1513     
1514     if (eIntegrator_f != eIntegrator)
1515     {
1516         if (MASTER(cr))
1517         {
1518             fprintf(stderr,int_warn,EI(eIntegrator_f),EI(eIntegrator));
1519         }
1520                 if(bAppendOutputFiles)
1521                 {
1522                         gmx_fatal(FARGS,
1523                                           "Output file appending requested, but input/checkpoint integrators do not match.\n"
1524                                           "Stopping the run to prevent you from ruining all your data...\n"
1525                                           "If you _really_ know what you are doing, try with the -noappend option.\n");
1526                 }
1527         if (fplog)
1528         {
1529             fprintf(fplog,int_warn,EI(eIntegrator_f),EI(eIntegrator));
1530         }
1531     }
1532
1533     if (!PAR(cr))
1534     {
1535         nppnodes = 1;
1536     }
1537     else if (bPartDecomp)
1538     {
1539         nppnodes = cr->nnodes;
1540     }
1541     else if (cr->nnodes == nppnodes_f + npmenodes_f)
1542     {
1543         if (cr->npmenodes < 0)
1544         {
1545             cr->npmenodes = npmenodes_f;
1546         }
1547         nppnodes = cr->nnodes - cr->npmenodes;
1548         if (nppnodes == nppnodes_f)
1549         {
1550             for(d=0; d<DIM; d++)
1551             {
1552                 if (dd_nc[d] == 0)
1553                 {
1554                     dd_nc[d] = dd_nc_f[d];
1555                 }
1556             }
1557         }
1558     }
1559     else
1560     {
1561         /* The number of PP nodes has not been set yet */
1562         nppnodes = -1;
1563     }
1564
1565     if ((EI_SD(eIntegrator) || eIntegrator == eiBD) && nppnodes > 0)
1566     {
1567         /* Correct the RNG state size for the number of PP nodes.
1568          * Such assignments should all be moved to one central function.
1569          */
1570         state->nrng  = nppnodes*gmx_rng_n();
1571         state->nrngi = nppnodes;
1572     }
1573     
1574     *bReadRNG = TRUE;
1575     if (fflags != state->flags)
1576     {
1577                 
1578         if (MASTER(cr))
1579         {
1580                         if(bAppendOutputFiles)
1581                         {
1582                                 gmx_fatal(FARGS,
1583                                                   "Output file appending requested, but input and checkpoint states are not identical.\n"
1584                                                   "Stopping the run to prevent you from ruining all your data...\n"
1585                                                   "You can try with the -noappend option, and get more info in the log file.\n");
1586                         }
1587                         
1588             fprintf(stderr,
1589                     "WARNING: The checkpoint state entries do not match the simulation,\n"
1590                     "         see the log file for details\n\n");
1591         }
1592                 
1593                 if(fplog)
1594                 {
1595                         print_flag_mismatch(fplog,state->flags,fflags);
1596                 }
1597     }
1598     else
1599     {
1600         if ((EI_SD(eIntegrator) || eIntegrator == eiBD) &&
1601             nppnodes != nppnodes_f)
1602         {
1603             *bReadRNG = FALSE;
1604             if (MASTER(cr))
1605             {
1606                 fprintf(stderr,sd_note,nppnodes_f,nppnodes);
1607             }
1608             if (fplog)
1609             {
1610                 fprintf(fplog ,sd_note,nppnodes_f,nppnodes);
1611             }
1612         }
1613         if (MASTER(cr))
1614         {
1615             check_match(fplog,version,btime,buser,bmach,fprog,
1616                         cr,bPartDecomp,nppnodes_f,npmenodes_f,dd_nc,dd_nc_f);
1617         }
1618     }
1619     ret = do_cpt_state(gmx_fio_getxdr(fp),TRUE,fflags,state,*bReadRNG,NULL);
1620     if (ret)
1621     {
1622         cp_error();
1623     }
1624     ret = do_cpt_ekinstate(gmx_fio_getxdr(fp),TRUE,
1625                            flags_eks,&state->ekinstate,NULL);
1626     if (ret)
1627     {
1628         cp_error();
1629     }
1630     *bReadEkin = ((flags_eks & (1<<eeksEKINH)) || (flags_eks & (1<<eeksEKINF)) || (flags_eks & (1<<eeksEKINO)) ||
1631                   (flags_eks & (1<<eeksEKINSCALEF)) | (flags_eks & (1<<eeksEKINSCALEH)) | (flags_eks & (1<<eeksVSCALE)));
1632     
1633     ret = do_cpt_enerhist(gmx_fio_getxdr(fp),TRUE,
1634                           flags_enh,&state->enerhist,NULL);
1635     if (ret)
1636     {
1637         cp_error();
1638     }
1639
1640     if (file_version < 6)
1641     {
1642         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.";
1643
1644         fprintf(stderr,"\nWARNING: %s\n\n",warn);
1645         if (fplog)
1646         {
1647             fprintf(fplog,"\nWARNING: %s\n\n",warn);
1648         }
1649         state->enerhist.nsum     = *step;
1650         state->enerhist.nsum_sim = *step;
1651     }
1652
1653         ret = do_cpt_files(gmx_fio_getxdr(fp),TRUE,&outputfiles,&nfiles,NULL,file_version);
1654         if (ret)
1655         {
1656                 cp_error();
1657         }
1658                                            
1659     ret = do_cpt_footer(gmx_fio_getxdr(fp),TRUE,file_version);
1660     if (ret)
1661     {
1662         cp_error();
1663     }
1664     if( gmx_fio_close(fp) != 0)
1665         {
1666                 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1667         }
1668     
1669     sfree(fprog);
1670     sfree(ftime);
1671     sfree(btime);
1672     sfree(buser);
1673     sfree(bmach);
1674         
1675         /* If the user wants to append to output files,
1676      * we use the file pointer positions of the output files stored
1677      * in the checkpoint file and truncate the files such that any frames
1678      * written after the checkpoint time are removed.
1679      * All files are md5sum checked such that we can be sure that
1680      * we do not truncate other (maybe imprortant) files.
1681          */
1682     if (bAppendOutputFiles)
1683     {
1684         if (fn2ftp(outputfiles[0].filename)!=efLOG)
1685         {
1686             /* make sure first file is log file so that it is OK to use it for 
1687              * locking
1688              */
1689             gmx_fatal(FARGS,"The first output file should always be the log "
1690                       "file but instead is: %s", outputfiles[0].filename);
1691         }
1692         for(i=0;i<nfiles;i++)
1693         {
1694             if (outputfiles[i].offset < 0)
1695             {
1696                 gmx_fatal(FARGS,"The original run wrote a file called '%s' which "
1697                     "is larger than 2 GB, but mdrun did not support large file"
1698                     " offsets. Can not append. Run mdrun with -noappend",
1699                     outputfiles[i].filename);
1700             }
1701 #ifdef GMX_FAHCORE
1702             chksum_file=gmx_fio_open(outputfiles[i].filename,"a");
1703
1704 #else
1705             chksum_file=gmx_fio_open(outputfiles[i].filename,"r+");
1706
1707             /* lock log file */                
1708             if (i==0)
1709             {
1710 #if !((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__) 
1711                 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl)
1712                     ==-1)
1713 #else
1714                 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX)==-1)
1715 #endif
1716                 {
1717                     gmx_fatal(FARGS,"Failed to lock: %s. Already running "
1718                         "simulation?", outputfiles[i].filename);
1719                 }
1720             }
1721             
1722             /* compute md5 chksum */ 
1723             if (outputfiles[i].chksum_size != -1)
1724             {
1725                 if (gmx_fio_get_file_md5(chksum_file,outputfiles[i].offset,
1726                                      digest) != outputfiles[i].chksum_size)
1727                 {
1728                     gmx_fatal(FARGS,"Can't read %d bytes of '%s' to compute checksum. The file has been replaced or its contents has been modified.",
1729                               outputfiles[i].chksum_size, 
1730                               outputfiles[i].filename);
1731                 }
1732             } 
1733             else if (i==0)  /*log file need to be seeked even when not reading md5*/
1734             {
1735                 gmx_fio_seek(chksum_file,outputfiles[i].offset);
1736             }
1737 #endif
1738             
1739             if (i==0) /*open log file here - so that lock is never lifted 
1740                         after chksum is calculated */
1741             {
1742                 *pfplog = gmx_fio_getfp(chksum_file);
1743             }
1744             else
1745             {
1746                 gmx_fio_close(chksum_file);
1747             }
1748 #ifndef GMX_FAHCORE            
1749             /* compare md5 chksum */
1750             if (outputfiles[i].chksum_size != -1 &&
1751                 memcmp(digest,outputfiles[i].chksum,16)!=0) 
1752             {
1753                 if (debug)
1754                 {
1755                     fprintf(debug,"chksum for %s: ",outputfiles[i].filename);
1756                     for (j=0; j<16; j++)
1757                     {
1758                         fprintf(debug,"%02x",digest[j]);
1759                     }
1760                     fprintf(debug,"\n");
1761                 }
1762                 gmx_fatal(FARGS,"Checksum wrong for '%s'. The file has been replaced or its contents has been modified.",
1763                           outputfiles[i].filename);
1764             }
1765 #endif        
1766
1767               
1768             if (i!=0) /*log file is already seeked to correct position */
1769             {
1770 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
1771                 rc = gmx_wintruncate(outputfiles[i].filename,outputfiles[i].offset);
1772 #else            
1773                 rc = truncate(outputfiles[i].filename,outputfiles[i].offset);
1774 #endif
1775                 if(rc!=0)
1776                 {
1777                     gmx_fatal(FARGS,"Truncation of file %s failed.",outputfiles[i].filename);
1778                 }
1779             }
1780         }
1781     }
1782
1783     sfree(outputfiles);
1784 }
1785
1786
1787 void load_checkpoint(const char *fn,FILE **fplog,
1788                      t_commrec *cr,gmx_bool bPartDecomp,ivec dd_nc,
1789                      t_inputrec *ir,t_state *state,
1790                      gmx_bool *bReadRNG,gmx_bool *bReadEkin,gmx_bool bAppend)
1791 {
1792     gmx_large_int_t step;
1793     double t;
1794
1795     if (SIMMASTER(cr)) {
1796       /* Read the state from the checkpoint file */
1797       read_checkpoint(fn,fplog,
1798                       cr,bPartDecomp,dd_nc,
1799                       ir->eI,&step,&t,state,bReadRNG,bReadEkin,
1800                       &ir->simulation_part,bAppend);
1801     }
1802     if (PAR(cr)) {
1803       gmx_bcast(sizeof(cr->npmenodes),&cr->npmenodes,cr);
1804       gmx_bcast(DIM*sizeof(dd_nc[0]),dd_nc,cr);
1805       gmx_bcast(sizeof(step),&step,cr);
1806       gmx_bcast(sizeof(*bReadRNG),bReadRNG,cr);
1807       gmx_bcast(sizeof(*bReadEkin),bReadEkin,cr);
1808     }
1809     ir->bContinuation    = TRUE;
1810     if (ir->nsteps >= 0)
1811     {
1812         ir->nsteps          += ir->init_step - step;
1813     }
1814     ir->init_step        = step;
1815         ir->simulation_part += 1;
1816 }
1817
1818 static void read_checkpoint_data(t_fileio *fp,int *simulation_part,
1819                                  gmx_large_int_t *step,double *t,t_state *state,
1820                                  gmx_bool bReadRNG,
1821                                  int *nfiles,gmx_file_position_t **outputfiles)
1822 {
1823     int  file_version;
1824     char *version,*btime,*buser,*bmach,*fprog,*ftime;
1825     int  eIntegrator;
1826     int  nppnodes,npme;
1827     ivec dd_nc;
1828     int  flags_eks,flags_enh;
1829     int  nfiles_loc;
1830     gmx_file_position_t *files_loc=NULL;
1831     int  ret;
1832         
1833     do_cpt_header(gmx_fio_getxdr(fp),TRUE,&file_version,
1834                   &version,&btime,&buser,&bmach,&fprog,&ftime,
1835                   &eIntegrator,simulation_part,step,t,&nppnodes,dd_nc,&npme,
1836                   &state->natoms,&state->ngtc,&state->nnhpres,&state->nhchainlength,
1837                   &state->flags,&flags_eks,&flags_enh,NULL);
1838     ret =
1839         do_cpt_state(gmx_fio_getxdr(fp),TRUE,state->flags,state,bReadRNG,NULL);
1840     if (ret)
1841     {
1842         cp_error();
1843     }
1844     ret = do_cpt_ekinstate(gmx_fio_getxdr(fp),TRUE,
1845                            flags_eks,&state->ekinstate,NULL);
1846     if (ret)
1847     {
1848         cp_error();
1849     }
1850     ret = do_cpt_enerhist(gmx_fio_getxdr(fp),TRUE,
1851                           flags_enh,&state->enerhist,NULL);
1852     if (ret)
1853     {
1854         cp_error();
1855     }
1856
1857     ret = do_cpt_files(gmx_fio_getxdr(fp),TRUE,
1858                        outputfiles != NULL ? outputfiles : &files_loc,
1859                        outputfiles != NULL ? nfiles : &nfiles_loc,
1860                        NULL,file_version);
1861     if (files_loc != NULL)
1862     {
1863         sfree(files_loc);
1864     }
1865         
1866     if (ret)
1867     {
1868         cp_error();
1869     }
1870         
1871     ret = do_cpt_footer(gmx_fio_getxdr(fp),TRUE,file_version);
1872     if (ret)
1873     {
1874         cp_error();
1875     }
1876
1877     sfree(fprog);
1878     sfree(ftime);
1879     sfree(btime);
1880     sfree(buser);
1881     sfree(bmach);
1882 }
1883
1884 void 
1885 read_checkpoint_state(const char *fn,int *simulation_part,
1886                       gmx_large_int_t *step,double *t,t_state *state)
1887 {
1888     t_fileio *fp;
1889     
1890     fp = gmx_fio_open(fn,"r");
1891     read_checkpoint_data(fp,simulation_part,step,t,state,TRUE,NULL,NULL);
1892     if( gmx_fio_close(fp) != 0)
1893         {
1894                 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1895         }
1896 }
1897
1898 void read_checkpoint_trxframe(t_fileio *fp,t_trxframe *fr)
1899 {
1900     t_state state;
1901     int simulation_part;
1902     gmx_large_int_t step;
1903     double t;
1904     
1905     init_state(&state,0,0,0,0);
1906     
1907     read_checkpoint_data(fp,&simulation_part,&step,&t,&state,FALSE,NULL,NULL);
1908     
1909     fr->natoms  = state.natoms;
1910     fr->bTitle  = FALSE;
1911     fr->bStep   = TRUE;
1912     fr->step    = gmx_large_int_to_int(step,
1913                                     "conversion of checkpoint to trajectory");
1914     fr->bTime   = TRUE;
1915     fr->time    = t;
1916     fr->bLambda = TRUE;
1917     fr->lambda  = state.lambda;
1918     fr->bAtoms  = FALSE;
1919     fr->bX      = (state.flags & (1<<estX));
1920     if (fr->bX)
1921     {
1922         fr->x     = state.x;
1923         state.x   = NULL;
1924     }
1925     fr->bV      = (state.flags & (1<<estV));
1926     if (fr->bV)
1927     {
1928         fr->v     = state.v;
1929         state.v   = NULL;
1930     }
1931     fr->bF      = FALSE;
1932     fr->bBox    = (state.flags & (1<<estBOX));
1933     if (fr->bBox)
1934     {
1935         copy_mat(state.box,fr->box);
1936     }
1937     done_state(&state);
1938 }
1939
1940 void list_checkpoint(const char *fn,FILE *out)
1941 {
1942     t_fileio *fp;
1943     int  file_version;
1944     char *version,*btime,*buser,*bmach,*fprog,*ftime;
1945     int  eIntegrator,simulation_part,nppnodes,npme;
1946     gmx_large_int_t step;
1947     double t;
1948     ivec dd_nc;
1949     t_state state;
1950     int  flags_eks,flags_enh;
1951     int  indent;
1952     int  i,j;
1953     int  ret;
1954     gmx_file_position_t *outputfiles;
1955         int  nfiles;
1956         
1957     init_state(&state,-1,-1,-1,-1);
1958
1959     fp = gmx_fio_open(fn,"r");
1960     do_cpt_header(gmx_fio_getxdr(fp),TRUE,&file_version,
1961                   &version,&btime,&buser,&bmach,&fprog,&ftime,
1962                   &eIntegrator,&simulation_part,&step,&t,&nppnodes,dd_nc,&npme,
1963                   &state.natoms,&state.ngtc,&state.nnhpres,&state.nhchainlength,
1964                   &state.flags,&flags_eks,&flags_enh,out);
1965     ret = do_cpt_state(gmx_fio_getxdr(fp),TRUE,state.flags,&state,TRUE,out);
1966     if (ret)
1967     {
1968         cp_error();
1969     }
1970     ret = do_cpt_ekinstate(gmx_fio_getxdr(fp),TRUE,
1971                            flags_eks,&state.ekinstate,out);
1972     if (ret)
1973     {
1974         cp_error();
1975     }
1976     ret = do_cpt_enerhist(gmx_fio_getxdr(fp),TRUE,
1977                           flags_enh,&state.enerhist,out);
1978
1979     if (ret == 0)
1980     {
1981                 do_cpt_files(gmx_fio_getxdr(fp),TRUE,&outputfiles,&nfiles,out,file_version);
1982         }
1983         
1984     if (ret == 0)
1985     {
1986         ret = do_cpt_footer(gmx_fio_getxdr(fp),TRUE,file_version);
1987     }
1988         
1989     if (ret)
1990     {
1991         cp_warning(out);
1992     }
1993     if( gmx_fio_close(fp) != 0)
1994         {
1995                 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1996         }
1997     
1998     done_state(&state);
1999 }
2000
2001
2002 static gmx_bool exist_output_file(const char *fnm_cp,int nfile,const t_filenm fnm[])
2003 {
2004     int i;
2005
2006     /* Check if the output file name stored in the checkpoint file
2007      * is one of the output file names of mdrun.
2008      */
2009     i = 0;
2010     while (i < nfile &&
2011            !(is_output(&fnm[i]) && strcmp(fnm_cp,fnm[i].fns[0]) == 0))
2012     {
2013         i++;
2014     }
2015     
2016     return (i < nfile && gmx_fexist(fnm_cp));
2017 }
2018
2019 /* This routine cannot print tons of data, since it is called before the log file is opened. */
2020 gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
2021                                      gmx_large_int_t *cpt_step,t_commrec *cr,
2022                                      gmx_bool bAppendReq,
2023                                      int nfile,const t_filenm fnm[],
2024                                      const char *part_suffix,gmx_bool *bAddPart)
2025 {
2026     t_fileio *fp;
2027     gmx_large_int_t step=0;
2028         double t;
2029     t_state state;
2030     int  nfiles;
2031     gmx_file_position_t *outputfiles;
2032     int  nexist,f;
2033     gmx_bool bAppend;
2034     char *fn,suf_up[STRLEN];
2035
2036     bAppend = FALSE;
2037
2038     if (SIMMASTER(cr)) {
2039         if(!gmx_fexist(filename) || (!(fp = gmx_fio_open(filename,"r")) ))
2040         {
2041             *simulation_part = 0;
2042         }
2043         else 
2044         {
2045             init_state(&state,0,0,0,0);
2046
2047             read_checkpoint_data(fp,simulation_part,&step,&t,&state,FALSE,
2048                                  &nfiles,&outputfiles);
2049             if( gmx_fio_close(fp) != 0)
2050             {
2051                 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
2052             }
2053             done_state(&state);
2054
2055             if (bAppendReq)
2056             {
2057                 nexist = 0;
2058                 for(f=0; f<nfiles; f++)
2059                 {
2060                     if (exist_output_file(outputfiles[f].filename,nfile,fnm))
2061                     {
2062                         nexist++;
2063                     }
2064                 }
2065                 if (nexist == nfiles)
2066                 {
2067                     bAppend = bAppendReq;
2068                 }
2069                 else if (nexist > 0)
2070                 {
2071                     fprintf(stderr,
2072                             "Output file appending has been requested,\n"
2073                             "but some output files listed in the checkpoint file %s\n"
2074                             "are not present or are named differently by the current program:\n",
2075                             filename);
2076                     fprintf(stderr,"output files present:");
2077                     for(f=0; f<nfiles; f++)
2078                     {
2079                         if (exist_output_file(outputfiles[f].filename,
2080                                               nfile,fnm))
2081                         {
2082                             fprintf(stderr," %s",outputfiles[f].filename);
2083                         }
2084                     }
2085                     fprintf(stderr,"\n");
2086                     fprintf(stderr,"output files not present or named differently:");
2087                     for(f=0; f<nfiles; f++)
2088                     {
2089                         if (!exist_output_file(outputfiles[f].filename,
2090                                                nfile,fnm))
2091                         {
2092                             fprintf(stderr," %s",outputfiles[f].filename);
2093                         }
2094                     }
2095                     fprintf(stderr,"\n");
2096                     
2097                     gmx_fatal(FARGS,"File appending requested, but only %d of the %d output files are present",nexist,nfiles);
2098                 }
2099             }
2100             
2101             if (bAppend)
2102             {
2103                 if (nfiles == 0)
2104                 {
2105                     gmx_fatal(FARGS,"File appending requested, but no output file information is stored in the checkpoint file");
2106                 }
2107                 fn = outputfiles[0].filename;
2108                 if (strlen(fn) < 4 ||
2109                     gmx_strcasecmp(fn+strlen(fn)-4,ftp2ext(efLOG)) == 0)
2110                 {
2111                     gmx_fatal(FARGS,"File appending requested, but the log file is not the first file listed in the checkpoint file");
2112                 }
2113                 /* Set bAddPart to whether the suffix string '.part' is present
2114                  * in the log file name.
2115                  */
2116                 strcpy(suf_up,part_suffix);
2117                 upstring(suf_up);
2118                 *bAddPart = (strstr(fn,part_suffix) != NULL ||
2119                              strstr(fn,suf_up) != NULL);
2120             }
2121
2122             sfree(outputfiles);
2123         }
2124     }
2125     if (PAR(cr))
2126     {
2127         gmx_bcast(sizeof(*simulation_part),simulation_part,cr);
2128
2129         if (*simulation_part > 0 && bAppendReq)
2130         {
2131             gmx_bcast(sizeof(bAppend),&bAppend,cr);
2132             gmx_bcast(sizeof(*bAddPart),bAddPart,cr);
2133         }
2134     }
2135     if (NULL != cpt_step)
2136     {
2137         *cpt_step = step;
2138     }
2139
2140     return bAppend;
2141 }