672629960401356a4bed7da548a1d91c280a0ba2
[alexxy/gromacs.git] / src / kernel / tpbcmp.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  * 
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  * 
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <math.h>
41 #include <stdio.h>
42 #include <string.h>
43 #include "main.h"
44 #include "macros.h"
45 #include "smalloc.h"
46 #include "futil.h"
47 #include "statutil.h"
48 #include "sysstuff.h"
49 #include "txtdump.h"
50 #include "gmx_fatal.h"
51 #include "names.h"
52 #include "tpxio.h"
53 #include "enxio.h"
54 #include "mtop_util.h"
55 #include "string2.h"
56 #include "tpbcmp.h"
57
58 static void cmp_int(FILE *fp,const char *s,int index,int i1,int i2)
59 {
60   if (i1 != i2) {
61     if (index != -1)
62       fprintf(fp,"%s[%d] (%d - %d)\n",s,index,i1,i2);
63     else
64       fprintf(fp,"%s (%d - %d)\n",s,i1,i2);
65   }
66 }
67
68 static void cmp_gmx_large_int(FILE *fp,const char *s,gmx_large_int_t i1,gmx_large_int_t i2)
69 {
70   if (i1 != i2) {
71     fprintf(fp,"%s (",s);
72     fprintf(fp,gmx_large_int_pfmt,i1);
73     fprintf(fp," - ");
74     fprintf(fp,gmx_large_int_pfmt,i2);
75     fprintf(fp,")\n");
76   }
77 }
78
79 static void cmp_us(FILE *fp,const char *s,int index,unsigned short i1,unsigned short i2)
80 {
81   if (i1 != i2) {
82     if (index != -1)
83       fprintf(fp,"%s[%d] (%d - %d)\n",s,index,i1,i2);
84     else
85       fprintf(fp,"%s (%d - %d)\n",s,i1,i2);
86   }
87 }
88
89 static void cmp_uc(FILE *fp,const char *s,int index,unsigned char i1,unsigned char i2)
90 {
91   if (i1 != i2) {
92     if (index != -1)
93       fprintf(fp,"%s[%d] (%d - %d)\n",s,index,i1,i2);
94     else
95       fprintf(fp,"%s (%d - %d)\n",s,i1,i2);
96   }
97 }
98
99 static gmx_bool cmp_bool(FILE *fp, const char *s, int index, gmx_bool b1, gmx_bool b2)
100 {
101   if (b1) {
102     b1 = 1;
103   } else {
104     b1 = 0;
105   }
106   if (b2) {
107     b2 = 1;
108   } else {
109     b2 = 0;
110   }
111   if (b1 != b2) {
112     if (index != -1)
113       fprintf(fp,"%s[%d] (%s - %s)\n",s,index,
114               bool_names[b1],bool_names[b2]);
115     else
116       fprintf(fp,"%s (%s - %s)\n",s,
117               bool_names[b1],bool_names[b2]);
118   }
119   return b1 && b2;
120 }
121
122 static void cmp_str(FILE *fp, const char *s, int index,
123                     const char *s1, const char *s2)
124 {
125   if (strcmp(s1,s2) != 0) {
126     if (index != -1)
127       fprintf(fp,"%s[%d] (%s - %s)\n",s,index,s1,s2);
128     else
129       fprintf(fp,"%s (%s - %s)\n",s,s1,s2);
130   }
131 }
132
133 static gmx_bool equal_real(real i1,real i2,real ftol,real abstol)
134 {
135     return ( ( 2*fabs(i1 - i2) <= (fabs(i1) + fabs(i2))*ftol ) || fabs(i1-i2)<=abstol );
136 }
137
138 static gmx_bool equal_float(float i1,float i2,float ftol,float abstol)
139 {
140     return ( ( 2*fabs(i1 - i2) <= (fabs(i1) + fabs(i2))*ftol ) || fabs(i1-i2)<=abstol );
141 }
142
143 static gmx_bool equal_double(double i1,double i2,real ftol,real abstol)
144 {
145     return ( ( 2*fabs(i1 - i2) <= (fabs(i1) + fabs(i2))*ftol ) || fabs(i1-i2)<=abstol );
146 }
147
148 static void 
149 cmp_real(FILE *fp,const char *s,int index,real i1,real i2,real ftol,real abstol)
150 {
151   if (!equal_real(i1,i2,ftol,abstol)) {
152     if (index != -1)
153       fprintf(fp,"%s[%2d] (%e - %e)\n",s,index,i1,i2);
154     else
155       fprintf(fp,"%s (%e - %e)\n",s,i1,i2);
156   }
157 }
158
159 static void 
160 cmp_float(FILE *fp,const char *s,int index,float i1,float i2,float ftol,float abstol)
161 {
162   if (!equal_float(i1,i2,ftol,abstol)) {
163     if (index != -1)
164       fprintf(fp,"%s[%2d] (%e - %e)\n",s,index,i1,i2);
165     else
166       fprintf(fp,"%s (%e - %e)\n",s,i1,i2);
167   }
168 }
169
170
171
172 static void 
173 cmp_double(FILE *fp,const char *s,int index,double i1,double i2,double ftol,double abstol)
174 {
175   if (!equal_double(i1,i2,ftol,abstol)) {
176     if (index != -1)
177       fprintf(fp,"%s[%2d] (%16.9e - %16.9e)\n",s,index,i1,i2);
178     else
179       fprintf(fp,"%s (%16.9e - %16.9e)\n",s,i1,i2);
180   }
181 }
182
183 static void cmp_rvec(FILE *fp,const char *s,int index,rvec i1,rvec i2,real ftol,real abstol)
184 {
185     if(!equal_real(i1[XX],i2[XX],ftol,abstol) ||
186        !equal_real(i1[YY],i2[YY],ftol,abstol) ||
187        !equal_real(i1[ZZ],i2[ZZ],ftol,abstol))
188     {
189         if (index != -1)
190             fprintf(fp,"%s[%5d] (%12.5e %12.5e %12.5e) - (%12.5e %12.5e %12.5e)\n",
191                     s,index,i1[XX],i1[YY],i1[ZZ],i2[XX],i2[YY],i2[ZZ]);
192         else
193             fprintf(fp,"%s (%12.5e %12.5e %12.5e) - (%12.5e %12.5e %12.5e)\n",
194                     s,i1[XX],i1[YY],i1[ZZ],i2[XX],i2[YY],i2[ZZ]);
195     }
196 }
197
198 static void cmp_ivec(FILE *fp,const char *s,int index,ivec i1,ivec i2)
199 {
200   if ((i1[XX] != i2[XX]) || (i1[YY] != i2[YY]) || (i1[ZZ] != i2[ZZ])) {
201     if (index != -1)
202       fprintf(fp,"%s[%5d] (%8d,%8d,%8d - %8d,%8d,%8d)\n",s,index,
203               i1[XX],i1[YY],i1[ZZ],i2[XX],i2[YY],i2[ZZ]);
204     else
205       fprintf(fp,"%s (%8d,%8d,%8d - %8d,%8d,%8d)\n",s,
206               i1[XX],i1[YY],i1[ZZ],i2[XX],i2[YY],i2[ZZ]);
207   }
208 }
209
210 static void cmp_ilist(FILE *fp,int ftype,t_ilist *il1,t_ilist *il2)
211 {
212   int i;
213   char buf[256];
214  
215   fprintf(fp,"comparing ilist %s\n",interaction_function[ftype].name);
216   sprintf(buf,"%s->nr",interaction_function[ftype].name);
217   cmp_int(fp,buf,-1,il1->nr,il2->nr);
218   sprintf(buf,"%s->iatoms",interaction_function[ftype].name);
219   if (((il1->nr > 0) && (!il1->iatoms)) || 
220       ((il2->nr > 0) && (!il2->iatoms)) ||
221       ((il1->nr != il2->nr)))
222     fprintf(fp,"Comparing radically different topologies - %s is different\n",
223             buf);
224   else
225     for(i=0; (i<il1->nr); i++) 
226       cmp_int(fp,buf,i,il1->iatoms[i],il2->iatoms[i]);
227 }
228
229 void cmp_iparm(FILE *fp,const char *s,t_functype ft,
230                t_iparams ip1,t_iparams ip2,real ftol,real abstol) 
231 {
232   int i;
233   gmx_bool bDiff;
234   
235   bDiff=FALSE;
236   for(i=0; i<MAXFORCEPARAM && !bDiff; i++)
237     bDiff = !equal_real(ip1.generic.buf[i],ip2.generic.buf[i],ftol,abstol);
238   if (bDiff) {
239     fprintf(fp,"%s1: ",s);
240     pr_iparams(fp,ft,&ip1);
241     fprintf(fp,"%s2: ",s);
242     pr_iparams(fp,ft,&ip2);
243   }
244 }
245
246 void cmp_iparm_AB(FILE *fp,const char *s,t_functype ft,t_iparams ip1,real ftol,real abstol) 
247 {
248   int nrfpA,nrfpB,p0,i;
249   gmx_bool bDiff;
250   
251   /* Normally the first parameter is perturbable */
252   p0 = 0;
253   nrfpA = interaction_function[ft].nrfpA;
254   nrfpB = interaction_function[ft].nrfpB;
255   if (ft == F_PDIHS) {
256     nrfpB = 2;
257   } else if (interaction_function[ft].flags & IF_TABULATED) {
258     /* For tabulated interactions only the second parameter is perturbable */
259     p0 = 1;
260     nrfpB = 1;
261   }
262   bDiff=FALSE;
263   for(i=0; i<nrfpB && !bDiff; i++) {
264     bDiff = !equal_real(ip1.generic.buf[p0+i],ip1.generic.buf[nrfpA+i],ftol,abstol);
265   }
266   if (bDiff) {
267     fprintf(fp,"%s: ",s);
268     pr_iparams(fp,ft,&ip1);
269   }
270 }
271
272 static void cmp_idef(FILE *fp,t_idef *id1,t_idef *id2,real ftol,real abstol)
273 {
274   int i;
275   char buf1[64],buf2[64];
276   
277   fprintf(fp,"comparing idef\n");
278   if (id2) {
279     cmp_int(fp,"idef->ntypes",-1,id1->ntypes,id2->ntypes);
280     cmp_int(fp,"idef->atnr",  -1,id1->atnr,id2->atnr);
281     for(i=0; (i<id1->ntypes); i++) {
282       sprintf(buf1,"idef->functype[%d]",i);
283       sprintf(buf2,"idef->iparam[%d]",i);
284       cmp_int(fp,buf1,i,(int)id1->functype[i],(int)id2->functype[i]);
285       cmp_iparm(fp,buf2,id1->functype[i],
286                 id1->iparams[i],id2->iparams[i],ftol,abstol);
287     }
288     cmp_real(fp,"fudgeQQ",-1,id1->fudgeQQ,id2->fudgeQQ,ftol,abstol);
289     for(i=0; (i<F_NRE); i++)
290       cmp_ilist(fp,i,&(id1->il[i]),&(id2->il[i]));
291   } else {
292     for(i=0; (i<id1->ntypes); i++)
293       cmp_iparm_AB(fp,"idef->iparam",id1->functype[i],id1->iparams[i],ftol,abstol);
294   }
295 }
296
297 static void cmp_block(FILE *fp,t_block *b1,t_block *b2,const char *s)
298 {
299   int i,j,k;
300   char buf[32];
301   
302   fprintf(fp,"comparing block %s\n",s);
303   sprintf(buf,"%s.nr",s);
304   cmp_int(fp,buf,-1,b1->nr,b2->nr);
305
306
307 static void cmp_blocka(FILE *fp,t_blocka *b1,t_blocka *b2,const char *s)
308 {
309   int i,j,k;
310   char buf[32];
311   
312   fprintf(fp,"comparing blocka %s\n",s);
313   sprintf(buf,"%s.nr",s);
314   cmp_int(fp,buf,-1,b1->nr,b2->nr);
315   sprintf(buf,"%s.nra",s);
316   cmp_int(fp,buf,-1,b1->nra,b2->nra);
317
318
319 static void cmp_atom(FILE *fp,int index,t_atom *a1,t_atom *a2,real ftol,real abstol)
320 {
321   int  i;
322   char buf[256];
323
324   if (a2) {
325     cmp_us(fp,"atom.type",index,a1->type,a2->type);
326     cmp_us(fp,"atom.ptype",index,a1->ptype,a2->ptype);
327     cmp_int(fp,"atom.resind",index,a1->resind,a2->resind);
328     cmp_int(fp,"atom.atomnumber",index,a1->atomnumber,a2->atomnumber);
329     cmp_real(fp,"atom.m",index,a1->m,a2->m,ftol,abstol);
330     cmp_real(fp,"atom.q",index,a1->q,a2->q,ftol,abstol);
331     cmp_us(fp,"atom.typeB",index,a1->typeB,a2->typeB);
332     cmp_real(fp,"atom.mB",index,a1->mB,a2->mB,ftol,abstol);
333     cmp_real(fp,"atom.qB",index,a1->qB,a2->qB,ftol,abstol);
334   } else {
335     cmp_us(fp,"atom.type",index,a1->type,a1->typeB);
336     cmp_real(fp,"atom.m",index,a1->m,a1->mB,ftol,abstol);
337     cmp_real(fp,"atom.q",index,a1->q,a1->qB,ftol,abstol);
338   }
339 }
340
341 static void cmp_atoms(FILE *fp,t_atoms *a1,t_atoms *a2,real ftol, real abstol)
342 {
343   int i;
344   
345   fprintf(fp,"comparing atoms\n");
346
347   if (a2) {
348     cmp_int(fp,"atoms->nr",-1,a1->nr,a2->nr);
349     for(i=0; (i<a1->nr); i++)
350       cmp_atom(fp,i,&(a1->atom[i]),&(a2->atom[i]),ftol,abstol);
351   } else {
352     for(i=0; (i<a1->nr); i++)
353       cmp_atom(fp,i,&(a1->atom[i]),NULL,ftol,abstol);
354   }
355 }
356
357 static void cmp_top(FILE *fp,t_topology *t1,t_topology *t2,real ftol, real abstol)
358 {
359   int i;
360   
361   fprintf(fp,"comparing top\n");
362   if (t2) {
363     cmp_idef(fp,&(t1->idef),&(t2->idef),ftol,abstol);
364     cmp_atoms(fp,&(t1->atoms),&(t2->atoms),ftol,abstol);
365     cmp_block(fp,&t1->cgs,&t2->cgs,"cgs");
366     cmp_block(fp,&t1->mols,&t2->mols,"mols");
367     cmp_blocka(fp,&t1->excls,&t2->excls,"excls");
368   } else {
369     cmp_idef(fp,&(t1->idef),NULL,ftol,abstol);
370     cmp_atoms(fp,&(t1->atoms),NULL,ftol,abstol);
371   }
372 }
373
374 static void cmp_groups(FILE *fp,gmx_groups_t *g0,gmx_groups_t *g1,
375                        int natoms0,int natoms1)
376 {
377   int  i,j,ndiff;
378   char buf[32];
379
380   fprintf(fp,"comparing groups\n");
381
382   for(i=0; i<egcNR; i++) {
383     sprintf(buf,"grps[%d].nr",i);
384     cmp_int(fp,buf,-1,g0->grps[i].nr,g1->grps[i].nr);
385     if (g0->grps[i].nr == g1->grps[i].nr) {
386       for(j=0; j<g0->grps[i].nr; j++) {
387           sprintf(buf,"grps[%d].name[%d]",i,j);
388           cmp_str(fp,buf,-1,
389                   *g0->grpname[g0->grps[i].nm_ind[j]],
390                   *g1->grpname[g1->grps[i].nm_ind[j]]);
391       }
392     }
393     cmp_int(fp,"ngrpnr",i,g0->ngrpnr[i],g1->ngrpnr[i]);
394     if (g0->ngrpnr[i] == g1->ngrpnr[i] && natoms0 == natoms1 && 
395         (g0->grpnr[i] != NULL || g1->grpnr[i] != NULL)) {
396       for(j=0; j<natoms0; j++) {
397         cmp_int(fp,gtypes[i],j,ggrpnr(g0,i,j),ggrpnr(g1,i,j));
398       }
399     }
400   }
401   /* We have compared the names in the groups lists,
402    * so we can skip the grpname list comparison.
403    */
404 }
405
406 static void cmp_rvecs(FILE *fp,const char *title,int n,rvec x1[],rvec x2[],
407                       gmx_bool bRMSD,real ftol,real abstol)
408 {
409   int i,m;
410   double d,ssd;
411
412   if (bRMSD) {
413     ssd = 0;
414     for(i=0; (i<n); i++) {
415       for(m=0; m<DIM; m++) {
416         d = x1[i][m] - x2[i][m];
417         ssd += d*d;
418       }
419     }
420     fprintf(fp,"%s RMSD %g\n",title,sqrt(ssd/n));
421   } else {
422     for(i=0; (i<n); i++) {
423       cmp_rvec(fp,title,i,x1[i],x2[i],ftol,abstol);
424     }
425   }
426 }
427
428
429 /* Similar to cmp_rvecs, but this routine scales the allowed absolute tolerance
430  * by the RMS of the force components of x1.
431  */
432 static void cmp_rvecs_rmstol(FILE *fp,const char *title,int n,rvec x1[],rvec x2[],
433                              real ftol,real abstol)
434 {
435     int i,m;
436     double d;
437     double ave_x1,rms_x1;
438     
439     /* It is tricky to compare real values, in particular forces that
440      * are sums of lots of terms where the final value might be close to 0.0.
441      * To get a reference magnitude we calculate the RMS value of each
442      * component in x1, and then set the allowed absolute tolerance to the
443      * relative tolerance times this RMS magnitude.
444      */
445     ave_x1 = 0.0;
446     for(i=0;i<n;i++)
447     {
448         for(m=0;m<DIM;m++)
449         {
450             ave_x1 += x1[i][m];
451         }
452     }
453     ave_x1 /= n*DIM;
454
455     rms_x1 = 0.0;
456     for(i=0; (i<n); i++)
457     {
458         for(m=0;m<DIM;m++)
459         {
460             d       = x1[i][m] - ave_x1;
461             rms_x1 += d*d;
462         }
463     }
464     rms_x1 = sqrt(rms_x1/(DIM*n));
465     /* And now do the actual comparision with a hopefully realistic abstol. */
466     for(i=0; (i<n); i++)
467     {
468         cmp_rvec(fp,title,i,x1[i],x2[i],ftol,abstol*rms_x1);
469     }
470 }
471
472 static void cmp_grpopts(FILE *fp,t_grpopts *opt1,t_grpopts *opt2,real ftol, real abstol)
473 {
474   int i,j;
475   char buf1[256],buf2[256];
476   
477   cmp_int(fp,"inputrec->grpopts.ngtc",-1,  opt1->ngtc,opt2->ngtc);
478   cmp_int(fp,"inputrec->grpopts.ngacc",-1, opt1->ngacc,opt2->ngacc);
479   cmp_int(fp,"inputrec->grpopts.ngfrz",-1, opt1->ngfrz,opt2->ngfrz);
480   cmp_int(fp,"inputrec->grpopts.ngener",-1,opt1->ngener,opt2->ngener);
481   for(i=0; (i<min(opt1->ngtc,opt2->ngtc)); i++) {
482     cmp_real(fp,"inputrec->grpopts.nrdf",i,opt1->nrdf[i],opt2->nrdf[i],ftol,abstol);
483     cmp_real(fp,"inputrec->grpopts.ref_t",i,opt1->ref_t[i],opt2->ref_t[i],ftol,abstol);
484     cmp_real(fp,"inputrec->grpopts.tau_t",i,opt1->tau_t[i],opt2->tau_t[i],ftol,abstol);
485     cmp_int(fp,"inputrec->grpopts.annealing",i,opt1->annealing[i],opt2->annealing[i]);
486     cmp_int(fp,"inputrec->grpopts.anneal_npoints",i,
487             opt1->anneal_npoints[i],opt2->anneal_npoints[i]);
488     if(opt1->anneal_npoints[i]==opt2->anneal_npoints[i]) {
489       sprintf(buf1,"inputrec->grpopts.anneal_time[%d]",i);
490       sprintf(buf2,"inputrec->grpopts.anneal_temp[%d]",i);
491       for(j=0;j<opt1->anneal_npoints[i];j++) {
492         cmp_real(fp,buf1,j,opt1->anneal_time[i][j],opt2->anneal_time[i][j],ftol,abstol);
493         cmp_real(fp,buf2,j,opt1->anneal_temp[i][j],opt2->anneal_temp[i][j],ftol,abstol);
494       }
495     }
496   }
497   if (opt1->ngener == opt2->ngener) {
498     for(i=0; i<opt1->ngener; i++)
499       for(j=i; j<opt1->ngener; j++) {
500         sprintf(buf1,"inputrec->grpopts.egp_flags[%d]",i);
501         cmp_int(fp,buf1,j,
502                 opt1->egp_flags[opt1->ngener*i+j],
503                 opt2->egp_flags[opt1->ngener*i+j]);
504       }
505   }
506   for(i=0; (i<min(opt1->ngacc,opt2->ngacc)); i++)
507     cmp_rvec(fp,"inputrec->grpopts.acc",i,opt1->acc[i],opt2->acc[i],ftol,abstol);
508   for(i=0; (i<min(opt1->ngfrz,opt2->ngfrz)); i++)
509     cmp_ivec(fp,"inputrec->grpopts.nFreeze",i,opt1->nFreeze[i],opt2->nFreeze[i]);
510 }
511
512 static void cmp_cosines(FILE *fp,const char *s,t_cosines c1[DIM],t_cosines c2[DIM],real ftol, real abstol)
513 {
514   int i,m;
515   char buf[256];
516   
517   for(m=0; (m<DIM); m++) {
518     sprintf(buf,"inputrec->%s[%d]",s,m);
519     cmp_int(fp,buf,0,c1->n,c2->n);
520     for(i=0; (i<min(c1->n,c2->n)); i++) {
521       cmp_real(fp,buf,i,c1->a[i],c2->a[i],ftol,abstol);
522       cmp_real(fp,buf,i,c1->phi[i],c2->phi[i],ftol,abstol);
523     }
524   }
525 }
526 static void cmp_adress(FILE *fp,t_adress *ad1,t_adress *ad2,
527                        real ftol,real abstol)
528 {
529   cmp_int(fp,"ir->adress->type" ,-1,ad1->type,ad2->type);
530   cmp_real(fp,"ir->adress->const_wf" ,-1,ad1->const_wf,ad2->const_wf,ftol,abstol);
531   cmp_real(fp,"ir->adress->ex_width" ,-1,ad1->ex_width,ad2->ex_width,ftol,abstol);
532   cmp_real(fp,"ir->adress->hy_width" ,-1,ad1->hy_width,ad2->hy_width,ftol,abstol);
533   cmp_int(fp,"ir->adress->icor" ,-1,ad1->icor,ad2->icor);
534   cmp_int(fp,"ir->adress->site" ,-1,ad1->site,ad2->site);
535   cmp_rvec(fp,"ir->adress->refs" ,-1,ad1->refs,ad2->refs,ftol,abstol);
536   cmp_real(fp,"ir->adress->ex_forcecap", -1,ad1->ex_forcecap,ad2->ex_forcecap,ftol,abstol);
537 }
538
539 static void cmp_pull(FILE *fp,t_pull *pull1,t_pull *pull2,real ftol, real abstol)
540 {
541   fprintf(fp,"WARNING: Both files use COM pulling, but comparing of the pull struct is not implemented (yet). The pull parameters could be the same or different.\n");
542 }
543
544 static void cmp_simtempvals(FILE *fp,t_simtemp *simtemp1,t_simtemp *simtemp2, int n_lambda, real ftol, real abstol)
545 {
546   int i;
547   cmp_int(fp,"inputrec->simtempvals->eSimTempScale",-1,simtemp1->eSimTempScale,simtemp2->eSimTempScale);
548   cmp_real(fp,"inputrec->simtempvals->simtemp_high",-1,simtemp1->simtemp_high,simtemp2->simtemp_high,ftol,abstol);
549   cmp_real(fp,"inputrec->simtempvals->simtemp_low",-1,simtemp1->simtemp_low,simtemp2->simtemp_low,ftol,abstol);
550   for(i=0; i<n_lambda; i++)
551   {
552       cmp_real(fp,"inputrec->simtempvals->temperatures",-1,simtemp1->temperatures[i],simtemp2->temperatures[i],ftol,abstol);
553   }
554 }
555
556 static void cmp_expandedvals(FILE *fp,t_expanded *expand1,t_expanded *expand2,int n_lambda, real ftol, real abstol)
557 {
558   int i;
559
560   cmp_bool(fp,"inputrec->fepvals->bInit_weights",-1,expand1->bInit_weights,expand2->bInit_weights);
561   cmp_bool(fp,"inputrec->fepvals->bWLoneovert",-1,expand1->bWLoneovert,expand2->bWLoneovert);
562
563   for(i=0; i<n_lambda; i++)
564   {
565       cmp_real(fp,"inputrec->expandedvals->init_lambda_weights",-1,
566                expand1->init_lambda_weights[i],expand2->init_lambda_weights[i],ftol,abstol);
567   }
568
569   cmp_int(fp,"inputrec->expandedvals->lambda-stats", -1,expand1->elamstats,expand2->elamstats);
570   cmp_int(fp,"inputrec->expandedvals->lambda-mc-move", -1,expand1->elmcmove,expand2->elmcmove);
571   cmp_int(fp,"inputrec->expandedvals->lmc-repeats",-1,expand1->lmc_repeats,expand2->lmc_repeats);
572   cmp_int(fp,"inputrec->expandedvals->lmc-gibbsdelta",-1,expand1->gibbsdeltalam,expand2->gibbsdeltalam);
573   cmp_int(fp,"inputrec->expandedvals->lmc-forced-nstart",-1,expand1->lmc_forced_nstart,expand2->lmc_forced_nstart);
574   cmp_int(fp,"inputrec->expandedvals->lambda-weights-equil",-1,expand1->elmceq,expand2->elmceq);
575   cmp_int(fp,"inputrec->expandedvals->,weight-equil-number-all-lambda",-1,expand1->equil_n_at_lam,expand2->equil_n_at_lam);
576   cmp_int(fp,"inputrec->expandedvals->weight-equil-number-samples",-1,expand1->equil_samples,expand2->equil_samples);
577   cmp_int(fp,"inputrec->expandedvals->weight-equil-number-steps",-1,expand1->equil_steps,expand2->equil_steps);
578   cmp_real(fp,"inputrec->expandedvals->weight-equil-wl-delta",-1,expand1->equil_wl_delta,expand2->equil_wl_delta,ftol,abstol);
579   cmp_real(fp,"inputrec->expandedvals->weight-equil-count-ratio",-1,expand1->equil_ratio,expand2->equil_ratio,ftol,abstol);
580   cmp_bool(fp,"inputrec->expandedvals->symmetrized-transition-matrix",-1,expand1->bSymmetrizedTMatrix,expand2->bSymmetrizedTMatrix);
581   cmp_int(fp,"inputrec->expandedvals->nstTij",-1,expand1->nstTij,expand2->nstTij);
582   cmp_int(fp,"inputrec->expandedvals->mininum-var-min",-1,expand1->minvarmin,expand2->minvarmin); /*default is reasonable */
583   cmp_int(fp,"inputrec->expandedvals->weight-c-range",-1,expand1->c_range,expand2->c_range); /* default is just C=0 */
584   cmp_real(fp,"inputrec->expandedvals->wl-scale",-1,expand1->wl_scale,expand2->wl_scale,ftol,abstol);
585   cmp_real(fp,"inputrec->expandedvals->init-wl-delta",-1,expand1->init_wl_delta,expand2->init_wl_delta,ftol,abstol);
586   cmp_real(fp,"inputrec->expandedvals->wl-ratio",-1,expand1->wl_ratio,expand2->wl_ratio,ftol,abstol);
587   cmp_int(fp,"inputrec->expandedvals->nstexpanded",-1,expand1->nstexpanded,expand2->nstexpanded);
588   cmp_int(fp,"inputrec->expandedvals->lmc-seed",-1,expand1->lmc_seed,expand2->lmc_seed);
589   cmp_real(fp,"inputrec->expandedvals->mc-temperature",-1,expand1->mc_temp,expand2->mc_temp,ftol,abstol);
590 }
591
592 static void cmp_fepvals(FILE *fp,t_lambda *fep1,t_lambda *fep2,real ftol, real abstol)
593 {
594   int i,j;
595   cmp_int(fp,"inputrec->nstdhdl",-1,fep1->nstdhdl,fep2->nstdhdl);
596   cmp_double(fp,"inputrec->fepvals->init_fep_state",-1,fep1->init_fep_state,fep2->init_fep_state,ftol,abstol);
597   cmp_double(fp,"inputrec->fepvals->delta_lambda",-1,fep1->delta_lambda,fep2->delta_lambda,ftol,abstol);
598   cmp_int(fp,"inputrec->fepvals->n_lambda",-1,fep1->n_lambda,fep2->n_lambda);
599   for(i=0; i<efptNR;i++)
600   {
601       for(j=0; j<min(fep1->n_lambda,fep2->n_lambda); j++)
602       {
603           cmp_double(fp,"inputrec->fepvals->all_lambda",-1,fep1->all_lambda[i][j],fep2->all_lambda[i][j],ftol,abstol);
604       }
605   }
606   cmp_real(fp,"inputrec->fepvals->sc_alpha",-1,fep1->sc_alpha,fep2->sc_alpha,ftol,abstol);
607   cmp_int(fp,"inputrec->fepvals->sc_power",-1,fep1->sc_power,fep2->sc_power);
608   cmp_real(fp,"inputrec->fepvals->sc_r_power",-1,fep1->sc_r_power,fep2->sc_r_power,ftol,abstol);
609   cmp_real(fp,"inputrec->fepvals->sc_sigma",-1,fep1->sc_sigma,fep2->sc_sigma,ftol,abstol);
610   cmp_bool(fp,"inputrec->fepvals->bPrintEnergy",-1,fep1->bPrintEnergy,fep1->bPrintEnergy);
611   cmp_bool(fp,"inputrec->fepvals->bScCoul",-1,fep1->bScCoul,fep1->bScCoul);
612   cmp_int(fp,"inputrec->separate_dhdl_file",-1,fep1->separate_dhdl_file,fep2->separate_dhdl_file);
613   cmp_int(fp,"inputrec->dhdl_derivatives",-1,fep1->dhdl_derivatives,fep2->dhdl_derivatives);
614   cmp_int(fp,"inputrec->dh_hist_size",-1,fep1->dh_hist_size,fep2->dh_hist_size);
615   cmp_double(fp,"inputrec->dh_hist_spacing",-1,fep1->dh_hist_spacing,fep2->dh_hist_spacing,ftol,abstol);
616 }
617
618 static void cmp_inputrec(FILE *fp,t_inputrec *ir1,t_inputrec *ir2,real ftol, real abstol)
619 {
620   fprintf(fp,"comparing inputrec\n");
621
622   /* gcc 2.96 doesnt like these defines at all, but issues a huge list
623    * of warnings. Maybe it will change in future versions, but for the
624    * moment I've spelled them out instead. /EL 000820 
625    * #define CIB(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
626    * #define CII(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
627    * #define CIR(s) cmp_real(fp,"inputrec->"#s,0,ir1->##s,ir2->##s,ftol)
628    */
629   cmp_int(fp,"inputrec->eI",-1,ir1->eI,ir2->eI);
630   cmp_gmx_large_int(fp,"inputrec->nsteps",ir1->nsteps,ir2->nsteps);
631   cmp_gmx_large_int(fp,"inputrec->init_step",ir1->init_step,ir2->init_step);
632   cmp_int(fp,"inputrec->simulation_part",-1,ir1->simulation_part,ir2->simulation_part);
633   cmp_int(fp,"inputrec->ePBC",-1,ir1->ePBC,ir2->ePBC);
634   cmp_int(fp,"inputrec->bPeriodicMols",-1,ir1->bPeriodicMols,ir2->bPeriodicMols);
635   cmp_int(fp,"inputrec->cutoff_scheme",-1,ir1->cutoff_scheme,ir2->cutoff_scheme);
636   cmp_int(fp,"inputrec->ns_type",-1,ir1->ns_type,ir2->ns_type);
637   cmp_int(fp,"inputrec->nstlist",-1,ir1->nstlist,ir2->nstlist);
638   cmp_int(fp,"inputrec->ndelta",-1,ir1->ndelta,ir2->ndelta);
639   cmp_int(fp,"inputrec->nstcomm",-1,ir1->nstcomm,ir2->nstcomm);
640   cmp_int(fp,"inputrec->comm_mode",-1,ir1->comm_mode,ir2->comm_mode);
641   cmp_int(fp,"inputrec->nstcheckpoint",-1,ir1->nstcheckpoint,ir2->nstcheckpoint);
642   cmp_int(fp,"inputrec->nstlog",-1,ir1->nstlog,ir2->nstlog);
643   cmp_int(fp,"inputrec->nstxout",-1,ir1->nstxout,ir2->nstxout);
644   cmp_int(fp,"inputrec->nstvout",-1,ir1->nstvout,ir2->nstvout);
645   cmp_int(fp,"inputrec->nstfout",-1,ir1->nstfout,ir2->nstfout);
646   cmp_int(fp,"inputrec->nstcalcenergy",-1,ir1->nstcalcenergy,ir2->nstcalcenergy);
647   cmp_int(fp,"inputrec->nstenergy",-1,ir1->nstenergy,ir2->nstenergy);
648   cmp_int(fp,"inputrec->nstxtcout",-1,ir1->nstxtcout,ir2->nstxtcout);
649   cmp_double(fp,"inputrec->init_t",-1,ir1->init_t,ir2->init_t,ftol,abstol);
650   cmp_double(fp,"inputrec->delta_t",-1,ir1->delta_t,ir2->delta_t,ftol,abstol);
651   cmp_real(fp,"inputrec->xtcprec",-1,ir1->xtcprec,ir2->xtcprec,ftol,abstol);
652   cmp_real(fp,"inputrec->fourierspacing",-1,ir1->fourier_spacing,ir2->fourier_spacing,ftol,abstol);
653   cmp_int(fp,"inputrec->nkx",-1,ir1->nkx,ir2->nkx);
654   cmp_int(fp,"inputrec->nky",-1,ir1->nky,ir2->nky);
655   cmp_int(fp,"inputrec->nkz",-1,ir1->nkz,ir2->nkz);
656   cmp_int(fp,"inputrec->pme_order",-1,ir1->pme_order,ir2->pme_order);
657   cmp_real(fp,"inputrec->ewald_rtol",-1,ir1->ewald_rtol,ir2->ewald_rtol,ftol,abstol);
658   cmp_int(fp,"inputrec->ewald_geometry",-1,ir1->ewald_geometry,ir2->ewald_geometry);
659   cmp_real(fp,"inputrec->epsilon_surface",-1,ir1->epsilon_surface,ir2->epsilon_surface,ftol,abstol);
660   cmp_int(fp,"inputrec->bOptFFT",-1,ir1->bOptFFT,ir2->bOptFFT);
661   cmp_int(fp,"inputrec->bContinuation",-1,ir1->bContinuation,ir2->bContinuation);
662   cmp_int(fp,"inputrec->bShakeSOR",-1,ir1->bShakeSOR,ir2->bShakeSOR);
663   cmp_int(fp,"inputrec->etc",-1,ir1->etc,ir2->etc);
664   cmp_int(fp,"inputrec->bPrintNHChains",-1,ir1->bPrintNHChains,ir2->bPrintNHChains);
665   cmp_int(fp,"inputrec->epc",-1,ir1->epc,ir2->epc);
666   cmp_int(fp,"inputrec->epct",-1,ir1->epct,ir2->epct);
667   cmp_real(fp,"inputrec->tau_p",-1,ir1->tau_p,ir2->tau_p,ftol,abstol);
668   cmp_rvec(fp,"inputrec->ref_p(x)",-1,ir1->ref_p[XX],ir2->ref_p[XX],ftol,abstol);
669   cmp_rvec(fp,"inputrec->ref_p(y)",-1,ir1->ref_p[YY],ir2->ref_p[YY],ftol,abstol);
670   cmp_rvec(fp,"inputrec->ref_p(z)",-1,ir1->ref_p[ZZ],ir2->ref_p[ZZ],ftol,abstol);
671   cmp_rvec(fp,"inputrec->compress(x)",-1,ir1->compress[XX],ir2->compress[XX],ftol,abstol);
672   cmp_rvec(fp,"inputrec->compress(y)",-1,ir1->compress[YY],ir2->compress[YY],ftol,abstol);
673   cmp_rvec(fp,"inputrec->compress(z)",-1,ir1->compress[ZZ],ir2->compress[ZZ],ftol,abstol);
674   cmp_int(fp,"refcoord_scaling",-1,ir1->refcoord_scaling,ir2->refcoord_scaling);
675   cmp_rvec(fp,"inputrec->posres_com",-1,ir1->posres_com,ir2->posres_com,ftol,abstol);
676   cmp_rvec(fp,"inputrec->posres_comB",-1,ir1->posres_comB,ir2->posres_comB,ftol,abstol);
677   cmp_real(fp,"inputrec->verletbuf_drift",-1,ir1->verletbuf_drift,ir2->verletbuf_drift,ftol,abstol);
678   cmp_real(fp,"inputrec->rlist",-1,ir1->rlist,ir2->rlist,ftol,abstol);
679   cmp_real(fp,"inputrec->rlistlong",-1,ir1->rlistlong,ir2->rlistlong,ftol,abstol);
680   cmp_int(fp,"inputrec->nstcalclr",-1,ir1->nstcalclr,ir2->nstcalclr);
681   cmp_real(fp,"inputrec->rtpi",-1,ir1->rtpi,ir2->rtpi,ftol,abstol);
682   cmp_int(fp,"inputrec->coulombtype",-1,ir1->coulombtype,ir2->coulombtype);
683   cmp_int(fp,"inputrec->coulomb_modifier",-1,ir1->coulomb_modifier,ir2->coulomb_modifier);
684   cmp_real(fp,"inputrec->rcoulomb_switch",-1,ir1->rcoulomb_switch,ir2->rcoulomb_switch,ftol,abstol);
685   cmp_real(fp,"inputrec->rcoulomb",-1,ir1->rcoulomb,ir2->rcoulomb,ftol,abstol);
686   cmp_int(fp,"inputrec->vdwtype",-1,ir1->vdwtype,ir2->vdwtype);
687   cmp_int(fp,"inputrec->vdw_modifier",-1,ir1->vdw_modifier,ir2->vdw_modifier);  cmp_real(fp,"inputrec->rvdw_switch",-1,ir1->rvdw_switch,ir2->rvdw_switch,ftol,abstol);
688   cmp_real(fp,"inputrec->rvdw",-1,ir1->rvdw,ir2->rvdw,ftol,abstol);
689   cmp_real(fp,"inputrec->epsilon_r",-1,ir1->epsilon_r,ir2->epsilon_r,ftol,abstol);
690   cmp_real(fp,"inputrec->epsilon_rf",-1,ir1->epsilon_rf,ir2->epsilon_rf,ftol,abstol);
691   cmp_real(fp,"inputrec->tabext",-1,ir1->tabext,ir2->tabext,ftol,abstol);
692   cmp_int(fp,"inputrec->implicit_solvent",-1,ir1->implicit_solvent,ir2->implicit_solvent);
693   cmp_int(fp,"inputrec->gb_algorithm",-1,ir1->gb_algorithm,ir2->gb_algorithm);
694   cmp_int(fp,"inputrec->nstgbradii",-1,ir1->nstgbradii,ir2->nstgbradii);
695   cmp_real(fp,"inputrec->rgbradii",-1,ir1->rgbradii,ir2->rgbradii,ftol,abstol);
696   cmp_real(fp,"inputrec->gb_saltconc",-1,ir1->gb_saltconc,ir2->gb_saltconc,ftol,abstol);
697   cmp_real(fp,"inputrec->gb_epsilon_solvent",-1,ir1->gb_epsilon_solvent,ir2->gb_epsilon_solvent,ftol,abstol);
698   cmp_real(fp,"inputrec->gb_obc_alpha",-1,ir1->gb_obc_alpha,ir2->gb_obc_alpha,ftol,abstol);
699   cmp_real(fp,"inputrec->gb_obc_beta",-1,ir1->gb_obc_beta,ir2->gb_obc_beta,ftol,abstol);
700   cmp_real(fp,"inputrec->gb_obc_gamma",-1,ir1->gb_obc_gamma,ir2->gb_obc_gamma,ftol,abstol);
701   cmp_real(fp,"inputrec->gb_dielectric_offset",-1,ir1->gb_dielectric_offset,ir2->gb_dielectric_offset,ftol,abstol);
702   cmp_int(fp,"inputrec->sa_algorithm",-1,ir1->sa_algorithm,ir2->sa_algorithm);
703   cmp_real(fp,"inputrec->sa_surface_tension",-1,ir1->sa_surface_tension,ir2->sa_surface_tension,ftol,abstol);   
704
705   cmp_int(fp,"inputrec->eDispCorr",-1,ir1->eDispCorr,ir2->eDispCorr);
706   cmp_real(fp,"inputrec->shake_tol",-1,ir1->shake_tol,ir2->shake_tol,ftol,abstol);
707   cmp_int(fp,"inputrec->efep",-1,ir1->efep,ir2->efep);
708   cmp_fepvals(fp,ir1->fepvals,ir2->fepvals,ftol,abstol);
709   cmp_int(fp,"inputrec->bSimTemp",-1,ir1->bSimTemp,ir2->bSimTemp);
710   if ((ir1->bSimTemp == ir2->bSimTemp) && (ir1->bSimTemp))
711   {
712       cmp_simtempvals(fp,ir1->simtempvals,ir2->simtempvals,min(ir1->fepvals->n_lambda,ir2->fepvals->n_lambda),ftol,abstol);
713   }
714   cmp_int(fp,"inputrec->bExpanded",-1,ir1->bExpanded,ir2->bExpanded);
715   if ((ir1->bExpanded == ir2->bExpanded) && (ir1->bExpanded))
716   {
717       cmp_expandedvals(fp,ir1->expandedvals,ir2->expandedvals,min(ir1->fepvals->n_lambda,ir2->fepvals->n_lambda),ftol,abstol);
718   }
719   cmp_int(fp,"inputrec->nwall",-1,ir1->nwall,ir2->nwall);
720   cmp_int(fp,"inputrec->wall_type",-1,ir1->wall_type,ir2->wall_type);
721   cmp_int(fp,"inputrec->wall_atomtype[0]",-1,ir1->wall_atomtype[0],ir2->wall_atomtype[0]);
722   cmp_int(fp,"inputrec->wall_atomtype[1]",-1,ir1->wall_atomtype[1],ir2->wall_atomtype[1]);
723   cmp_real(fp,"inputrec->wall_density[0]",-1,ir1->wall_density[0],ir2->wall_density[0],ftol,abstol);
724   cmp_real(fp,"inputrec->wall_density[1]",-1,ir1->wall_density[1],ir2->wall_density[1],ftol,abstol);
725   cmp_real(fp,"inputrec->wall_ewald_zfac",-1,ir1->wall_ewald_zfac,ir2->wall_ewald_zfac,ftol,abstol);
726
727   cmp_int(fp,"inputrec->ePull",-1,ir1->ePull,ir2->ePull);
728   if (ir1->ePull == ir2->ePull && ir1->ePull != epullNO)
729     cmp_pull(fp,ir1->pull,ir2->pull,ftol,abstol);
730   
731   cmp_int(fp,"inputrec->eDisre",-1,ir1->eDisre,ir2->eDisre);
732   cmp_real(fp,"inputrec->dr_fc",-1,ir1->dr_fc,ir2->dr_fc,ftol,abstol);
733   cmp_int(fp,"inputrec->eDisreWeighting",-1,ir1->eDisreWeighting,ir2->eDisreWeighting);
734   cmp_int(fp,"inputrec->bDisreMixed",-1,ir1->bDisreMixed,ir2->bDisreMixed);
735   cmp_int(fp,"inputrec->nstdisreout",-1,ir1->nstdisreout,ir2->nstdisreout);
736   cmp_real(fp,"inputrec->dr_tau",-1,ir1->dr_tau,ir2->dr_tau,ftol,abstol);
737   cmp_real(fp,"inputrec->orires_fc",-1,ir1->orires_fc,ir2->orires_fc,ftol,abstol);
738   cmp_real(fp,"inputrec->orires_tau",-1,ir1->orires_tau,ir2->orires_tau,ftol,abstol);
739   cmp_int(fp,"inputrec->nstorireout",-1,ir1->nstorireout,ir2->nstorireout);
740   cmp_real(fp,"inputrec->dihre_fc",-1,ir1->dihre_fc,ir2->dihre_fc,ftol,abstol);
741   cmp_real(fp,"inputrec->em_stepsize",-1,ir1->em_stepsize,ir2->em_stepsize,ftol,abstol);
742   cmp_real(fp,"inputrec->em_tol",-1,ir1->em_tol,ir2->em_tol,ftol,abstol);
743   cmp_int(fp,"inputrec->niter",-1,ir1->niter,ir2->niter);
744   cmp_real(fp,"inputrec->fc_stepsize",-1,ir1->fc_stepsize,ir2->fc_stepsize,ftol,abstol);
745   cmp_int(fp,"inputrec->nstcgsteep",-1,ir1->nstcgsteep,ir2->nstcgsteep);
746   cmp_int(fp,"inputrec->nbfgscorr",0,ir1->nbfgscorr,ir2->nbfgscorr);
747   cmp_int(fp,"inputrec->eConstrAlg",-1,ir1->eConstrAlg,ir2->eConstrAlg);
748   cmp_int(fp,"inputrec->nProjOrder",-1,ir1->nProjOrder,ir2->nProjOrder);
749   cmp_real(fp,"inputrec->LincsWarnAngle",-1,ir1->LincsWarnAngle,ir2->LincsWarnAngle,ftol,abstol);
750   cmp_int(fp,"inputrec->nLincsIter",-1,ir1->nLincsIter,ir2->nLincsIter);
751   cmp_real(fp,"inputrec->bd_fric",-1,ir1->bd_fric,ir2->bd_fric,ftol,abstol);
752   cmp_int(fp,"inputrec->ld_seed",-1,ir1->ld_seed,ir2->ld_seed);
753   cmp_real(fp,"inputrec->cos_accel",-1,ir1->cos_accel,ir2->cos_accel,ftol,abstol);
754   cmp_rvec(fp,"inputrec->deform(a)",-1,ir1->deform[XX],ir2->deform[XX],ftol,abstol);
755   cmp_rvec(fp,"inputrec->deform(b)",-1,ir1->deform[YY],ir2->deform[YY],ftol,abstol);
756   cmp_rvec(fp,"inputrec->deform(c)",-1,ir1->deform[ZZ],ir2->deform[ZZ],ftol,abstol);
757
758   
759   cmp_bool(fp,"ir->bAdress->type" ,-1,ir1->bAdress,ir2->bAdress);
760   if (ir1->bAdress && ir2->bAdress) {
761       cmp_adress(fp,ir1->adress,ir2->adress,ftol,abstol);
762   }
763
764   cmp_int(fp,"inputrec->userint1",-1,ir1->userint1,ir2->userint1);
765   cmp_int(fp,"inputrec->userint2",-1,ir1->userint2,ir2->userint2);
766   cmp_int(fp,"inputrec->userint3",-1,ir1->userint3,ir2->userint3);
767   cmp_int(fp,"inputrec->userint4",-1,ir1->userint4,ir2->userint4);
768   cmp_real(fp,"inputrec->userreal1",-1,ir1->userreal1,ir2->userreal1,ftol,abstol);
769   cmp_real(fp,"inputrec->userreal2",-1,ir1->userreal2,ir2->userreal2,ftol,abstol);
770   cmp_real(fp,"inputrec->userreal3",-1,ir1->userreal3,ir2->userreal3,ftol,abstol);
771   cmp_real(fp,"inputrec->userreal4",-1,ir1->userreal4,ir2->userreal4,ftol,abstol);
772   cmp_grpopts(fp,&(ir1->opts),&(ir2->opts),ftol,abstol);
773   cmp_cosines(fp,"ex",ir1->ex,ir2->ex,ftol,abstol);
774   cmp_cosines(fp,"et",ir1->et,ir2->et,ftol,abstol);
775 }
776
777 static void comp_pull_AB(FILE *fp,t_pull *pull,real ftol,real abstol)
778 {
779   int i;
780
781   for(i=0; i<pull->ngrp+1; i++) {
782     fprintf(fp,"comparing pull group %d\n",i);
783     cmp_real(fp,"pullgrp->k",-1,pull->grp[i].k,pull->grp[i].kB,ftol,abstol);
784   }
785 }
786
787 static void comp_state(t_state *st1, t_state *st2,
788                        gmx_bool bRMSD,real ftol,real abstol)
789 {
790   int i,j,nc;
791
792   fprintf(stdout,"comparing flags\n");
793   cmp_int(stdout,"flags",-1,st1->flags,st2->flags);
794   fprintf(stdout,"comparing box\n");
795   cmp_rvecs(stdout,"box",DIM,st1->box,st2->box,FALSE,ftol,abstol);
796   fprintf(stdout,"comparing box_rel\n");
797   cmp_rvecs(stdout,"box_rel",DIM,st1->box_rel,st2->box_rel,FALSE,ftol,abstol);
798   fprintf(stdout,"comparing boxv\n");
799   cmp_rvecs(stdout,"boxv",DIM,st1->boxv,st2->boxv,FALSE,ftol,abstol);
800   if (st1->flags & (1<<estSVIR_PREV)) {
801     fprintf(stdout,"comparing shake vir_prev\n");
802     cmp_rvecs_rmstol(stdout,"svir_prev",DIM,st1->svir_prev,st2->svir_prev,ftol,abstol);
803   }
804   if (st1->flags & (1<<estFVIR_PREV)) {
805     fprintf(stdout,"comparing force vir_prev\n");
806     cmp_rvecs_rmstol(stdout,"fvir_prev",DIM,st1->fvir_prev,st2->fvir_prev,ftol,abstol);
807   }
808   if (st1->flags & (1<<estPRES_PREV)) {
809     fprintf(stdout,"comparing prev_pres\n");
810     cmp_rvecs_rmstol(stdout,"pres_prev",DIM,st1->pres_prev,st2->pres_prev,ftol,abstol);
811   }
812   cmp_int(stdout,"ngtc",-1,st1->ngtc,st2->ngtc);
813   cmp_int(stdout,"nhchainlength",-1,st1->nhchainlength,st2->nhchainlength);
814   if (st1->ngtc == st2->ngtc && st1->nhchainlength == st2->nhchainlength){
815     for(i=0; i<st1->ngtc; i++) {
816       nc = i*st1->nhchainlength;
817       for(j=0; j<nc; j++) {
818         cmp_real(stdout,"nosehoover_xi",
819                  i,st1->nosehoover_xi[nc+j],st2->nosehoover_xi[nc+j],ftol,abstol);
820       }
821     }
822   }
823   cmp_int(stdout,"nnhpres",-1,st1->nnhpres,st2->nnhpres);
824   if (st1->nnhpres == st2->nnhpres && st1->nhchainlength == st2->nhchainlength) {
825     for(i=0; i<st1->nnhpres; i++) {
826       nc = i*st1->nhchainlength;
827       for(j=0; j<nc; j++) {
828         cmp_real(stdout,"nosehoover_xi",
829                  i,st1->nhpres_xi[nc+j],st2->nhpres_xi[nc+j],ftol,abstol);
830       }
831     }
832   }
833
834   cmp_int(stdout,"natoms",-1,st1->natoms,st2->natoms);
835   if (st1->natoms == st2->natoms) {
836     if ((st1->flags & (1<<estX)) && (st2->flags & (1<<estX))) {
837       fprintf(stdout,"comparing x\n");
838       cmp_rvecs(stdout,"x",st1->natoms,st1->x,st2->x,bRMSD,ftol,abstol);
839     }
840     if ((st1->flags & (1<<estV)) && (st2->flags & (1<<estV))) {
841       fprintf(stdout,"comparing v\n");
842       cmp_rvecs(stdout,"v",st1->natoms,st1->v,st2->v,bRMSD,ftol,abstol);
843     }
844   }
845 }
846
847 void comp_tpx(const char *fn1,const char *fn2,
848               gmx_bool bRMSD,real ftol,real abstol)
849 {
850   const char  *ff[2];
851   t_tpxheader sh[2];
852   t_inputrec  ir[2];
853   t_state     state[2];
854   gmx_mtop_t  mtop[2];
855   t_topology  top[2];
856   int         i;
857
858   ff[0]=fn1;
859   ff[1]=fn2;
860   for(i=0; i<(fn2 ? 2 : 1); i++) {
861     read_tpx_state(ff[i],&(ir[i]),&state[i],NULL,&(mtop[i]));
862   }
863   if (fn2) {
864     cmp_inputrec(stdout,&ir[0],&ir[1],ftol,abstol);
865     /* Convert gmx_mtop_t to t_topology.
866      * We should implement direct mtop comparison,
867      * but it might be useful to keep t_topology comparison as an option.
868      */
869     top[0] = gmx_mtop_t_to_t_topology(&mtop[0]);
870     top[1] = gmx_mtop_t_to_t_topology(&mtop[1]);
871     cmp_top(stdout,&top[0],&top[1],ftol,abstol);
872     cmp_groups(stdout,&mtop[0].groups,&mtop[1].groups,
873                mtop[0].natoms,mtop[1].natoms);
874     comp_state(&state[0],&state[1],bRMSD,ftol,abstol);
875   } else {
876     if (ir[0].efep == efepNO) {
877       fprintf(stdout,"inputrec->efep = %s\n",efep_names[ir[0].efep]);
878     } else {
879       if (ir[0].ePull != epullNO) {
880         comp_pull_AB(stdout,ir->pull,ftol,abstol);
881       }
882       /* Convert gmx_mtop_t to t_topology.
883        * We should implement direct mtop comparison,
884        * but it might be useful to keep t_topology comparison as an option.
885        */
886       top[0] = gmx_mtop_t_to_t_topology(&mtop[0]);
887       cmp_top(stdout,&top[0],NULL,ftol,abstol);
888     }
889   }
890 }
891
892 void comp_frame(FILE *fp, t_trxframe *fr1, t_trxframe *fr2,
893                 gmx_bool bRMSD, real ftol,real abstol)
894 {
895   fprintf(fp,"\n");
896   cmp_int(fp,"flags",-1,fr1->flags,fr2->flags);
897   cmp_int(fp,"not_ok",-1,fr1->not_ok,fr2->not_ok);
898   cmp_int(fp,"natoms",-1,fr1->natoms,fr2->natoms);
899   cmp_real(fp,"t0",-1,fr1->t0,fr2->t0,ftol,abstol);
900   if (cmp_bool(fp,"bTitle",-1,fr1->bTitle,fr2->bTitle))
901     cmp_str(fp,"title", -1, fr1->title, fr2->title);
902   if (cmp_bool(fp,"bStep",-1,fr1->bStep,fr2->bStep))
903     cmp_int(fp,"step",-1,fr1->step,fr2->step);
904   cmp_int(fp,"step",-1,fr1->step,fr2->step);
905   if (cmp_bool(fp,"bTime",-1,fr1->bTime,fr2->bTime))   
906     cmp_real(fp,"time",-1,fr1->time,fr2->time,ftol,abstol);
907   if (cmp_bool(fp,"bLambda",-1,fr1->bLambda,fr2->bLambda)) 
908     cmp_real(fp,"lambda",-1,fr1->lambda,fr2->lambda,ftol,abstol);
909   if (cmp_bool(fp,"bAtoms",-1,fr1->bAtoms,fr2->bAtoms))
910     cmp_atoms(fp,fr1->atoms,fr2->atoms,ftol,abstol);
911   if (cmp_bool(fp,"bPrec",-1,fr1->bPrec,fr2->bPrec))
912     cmp_real(fp,"prec",-1,fr1->prec,fr2->prec,ftol,abstol);
913   if (cmp_bool(fp,"bX",-1,fr1->bX,fr2->bX))
914     cmp_rvecs(fp,"x",min(fr1->natoms,fr2->natoms),fr1->x,fr2->x,bRMSD,ftol,abstol);
915   if (cmp_bool(fp,"bV",-1,fr1->bV,fr2->bV))
916     cmp_rvecs(fp,"v",min(fr1->natoms,fr2->natoms),fr1->v,fr2->v,bRMSD,ftol,abstol);
917   if (cmp_bool(fp,"bF",-1,fr1->bF,fr2->bF))
918     cmp_rvecs_rmstol(fp,"f",min(fr1->natoms,fr2->natoms),fr1->f,fr2->f,ftol,abstol);
919   if (cmp_bool(fp,"bBox",-1,fr1->bBox,fr2->bBox))
920     cmp_rvecs(fp,"box",3,fr1->box,fr2->box,FALSE,ftol,abstol);
921 }
922
923 void comp_trx(const output_env_t oenv,const char *fn1, const char *fn2, 
924               gmx_bool bRMSD,real ftol,real abstol)
925 {
926   int i;
927   const char *fn[2];
928   t_trxframe fr[2];
929   t_trxstatus *status[2];
930   gmx_bool b[2];
931   
932   fn[0]=fn1;
933   fn[1]=fn2;
934   fprintf(stderr,"Comparing trajectory files %s and %s\n",fn1,fn2);
935   for (i=0; i<2; i++)
936     b[i] = read_first_frame(oenv,&status[i],fn[i],&fr[i],TRX_READ_X|TRX_READ_V|TRX_READ_F);
937   
938   if (b[0] && b[1]) { 
939     do {
940       comp_frame(stdout, &(fr[0]), &(fr[1]), bRMSD, ftol, abstol);
941       
942       for (i=0; i<2; i++)
943         b[i] = read_next_frame(oenv,status[i],&fr[i]);
944     } while (b[0] && b[1]);
945     
946     for (i=0; i<2; i++) {
947       if (b[i] && !b[1-i])
948         fprintf(stdout,"\nEnd of file on %s but not on %s\n",fn[1-i],fn[i]);
949       close_trj(status[i]);
950     }
951   }
952   if (!b[0] && !b[1])
953     fprintf(stdout,"\nBoth files read correctly\n");
954 }
955
956 static real ener_tensor_diag(int n,int *ind1,int *ind2,
957                              gmx_enxnm_t *enm1,
958                              int *tensi,int i,
959                              t_energy e1[],t_energy e2[])
960 {
961   int  d1,d2;
962   int  len;
963   int  j;
964   real prod1,prod2;
965   int  nfound;
966
967   d1 = tensi[i]/DIM;
968   d2 = tensi[i] - d1*DIM;
969   
970   /* Find the diagonal elements d1 and d2 */
971   len = strlen(enm1[ind1[i]].name);
972   prod1 = 1;
973   prod2 = 1;
974   nfound = 0;
975   for(j=0; j<n; j++) {
976     if (tensi[j] >= 0 &&
977         strlen(enm1[ind1[j]].name) == len &&
978         strncmp(enm1[ind1[i]].name,enm1[ind1[j]].name,len-2) == 0 &&
979         (tensi[j] == d1*DIM+d1 || tensi[j] == d2*DIM+d2)) {
980       prod1 *= fabs(e1[ind1[j]].e);
981       prod2 *= fabs(e2[ind2[j]].e);
982       nfound++;
983     }
984   }
985
986   if (nfound == 2) {
987     return 0.5*(sqrt(prod1) + sqrt(prod2));
988   } else {
989     return 0;
990   }
991 }
992
993 static gmx_bool enernm_equal(const char *nm1,const char *nm2)
994 {
995   int len1,len2;
996
997   len1 = strlen(nm1);
998   len2 = strlen(nm2);
999
1000   /* Remove " (bar)" at the end of a name */
1001   if (len1 > 6 && strcmp(nm1+len1-6," (bar)") == 0) {
1002     len1 -= 6;
1003   }
1004   if (len2 > 6 && strcmp(nm2+len2-6," (bar)") == 0) {
1005     len2 -= 6;
1006   }
1007
1008   return (len1 == len2 && gmx_strncasecmp(nm1,nm2,len1) == 0);
1009 }
1010
1011 static void cmp_energies(FILE *fp,int step1,int step2,
1012                          t_energy e1[],t_energy e2[],
1013                          gmx_enxnm_t *enm1,gmx_enxnm_t *enm2,
1014                          real ftol,real abstol,
1015                          int nre,int *ind1,int *ind2,int maxener)
1016 {
1017   int  i,ii;
1018   int  *tensi,len,d1,d2;
1019   real ftol_i,abstol_i;
1020
1021   snew(tensi,maxener);
1022   /* Check for tensor elements ending on "-XX", "-XY", ... , "-ZZ" */
1023   for(i=0; (i<maxener); i++) {
1024     ii = ind1[i];
1025     tensi[i] = -1;
1026     len = strlen(enm1[ii].name);
1027     if (len > 3 && enm1[ii].name[len-3] == '-') {
1028       d1 = enm1[ii].name[len-2] - 'X';
1029       d2 = enm1[ii].name[len-1] - 'X';
1030       if (d1 >= 0 && d1 < DIM &&
1031           d2 >= 0 && d2 < DIM) {
1032         tensi[i] = d1*DIM + d2;
1033       }
1034     }
1035   }
1036   
1037   for(i=0; (i<maxener); i++) {
1038     /* Check if this is an off-diagonal tensor element */
1039     if (tensi[i] >= 0 && tensi[i] != 0 && tensi[i] != 4 && tensi[i] != 8) {
1040       /* Turn on the relative tolerance check (4 is maximum relative diff.) */
1041       ftol_i = 5;
1042       /* Do the relative tolerance through an absolute tolerance times
1043        * the size of diagonal components of the tensor.
1044        */
1045       abstol_i = ftol*ener_tensor_diag(nre,ind1,ind2,enm1,tensi,i,e1,e2);
1046       if (debug) {
1047         fprintf(debug,"tensor '%s' val %f diag %f\n",
1048                 enm1[i].name,e1[i].e,abstol_i/ftol);
1049       }
1050       if (abstol_i > 0) {
1051         /* We found a diagonal, we need to check with the minimum tolerance */
1052         abstol_i = min(abstol_i,abstol);
1053       } else {
1054         /* We did not find a diagonal, ignore the relative tolerance check */
1055         abstol_i = abstol;
1056       }
1057     } else {
1058       ftol_i   = ftol;
1059       abstol_i = abstol;
1060     }
1061     if (!equal_real(e1[ind1[i]].e,e2[ind2[i]].e,ftol_i,abstol_i)) {
1062       fprintf(fp,"%-15s  step %3d:  %12g,  step %3d: %12g\n",
1063               enm1[ind1[i]].name,
1064               step1,e1[ind1[i]].e,
1065               step2,e2[ind2[i]].e);
1066     }
1067   }
1068
1069   sfree(tensi);
1070 }
1071
1072 #if 0
1073 static void cmp_disres(t_enxframe *fr1,t_enxframe *fr2,real ftol, real abstol)
1074 {
1075   int i;
1076   char bav[64],bt[64],bs[22];
1077     
1078   cmp_int(stdout,"ndisre",-1,fr1->ndisre,fr2->ndisre);
1079   if ((fr1->ndisre == fr2->ndisre) && (fr1->ndisre > 0)) {
1080     sprintf(bav,"step %s: disre rav",gmx_step_str(fr1->step,bs));
1081     sprintf(bt, "step %s: disre  rt",gmx_step_str(fr1->step,bs));
1082     for(i=0; (i<fr1->ndisre); i++) {
1083       cmp_real(stdout,bav,i,fr1->disre_rm3tav[i],fr2->disre_rm3tav[i],ftol,abstol);
1084       cmp_real(stdout,bt ,i,fr1->disre_rt[i]    ,fr2->disre_rt[i]    ,ftol,abstol);
1085     }
1086   }
1087 }
1088 #endif
1089
1090 static void cmp_eblocks(t_enxframe *fr1,t_enxframe *fr2,real ftol, real abstol)
1091 {
1092     int i,j,k;
1093     char buf[64],bs[22];
1094
1095     cmp_int(stdout,"nblock",-1,fr1->nblock,fr2->nblock);  
1096     if ((fr1->nblock == fr2->nblock) && (fr1->nblock > 0)) 
1097     {
1098         for(j=0; (j<fr1->nblock); j++) 
1099         {
1100             t_enxblock *b1, *b2; /* convenience vars */
1101
1102             b1=&(fr1->block[j]);
1103             b2=&(fr2->block[j]);
1104
1105             sprintf(buf,"step %s: block[%d]",gmx_step_str(fr1->step,bs),j);
1106             cmp_int(stdout,buf,-1,b1->nsub,b2->nsub);
1107             cmp_int(stdout,buf,-1,b1->id,b2->id);
1108
1109             if ( (b1->nsub==b2->nsub) && (b1->id==b2->id) )
1110             {
1111                 for(i=0;i<b1->nsub;i++)
1112                 {
1113                     t_enxsubblock *s1, *s2;
1114
1115                     s1=&(b1->sub[i]);
1116                     s2=&(b2->sub[i]);
1117
1118                     cmp_int(stdout, buf, -1, (int)s1->type, (int)s2->type);
1119                     cmp_gmx_large_int(stdout, buf, s1->nr, s2->nr);
1120
1121                     if ((s1->type == s2->type) && (s1->nr == s2->nr))
1122                     {
1123                         switch(s1->type)
1124                         {
1125                             case xdr_datatype_float:
1126                                 for(k=0;k<s1->nr;k++)
1127                                 {
1128                                     cmp_float(stdout, buf, i,
1129                                              s1->fval[k], s2->fval[k], 
1130                                              ftol, abstol);
1131                                 }
1132                                 break;
1133                             case xdr_datatype_double:
1134                                 for(k=0;k<s1->nr;k++)
1135                                 {
1136                                     cmp_double(stdout, buf, i,
1137                                              s1->dval[k], s2->dval[k], 
1138                                              ftol, abstol);
1139                                 }
1140                                 break;
1141                             case xdr_datatype_int:
1142                                 for(k=0;k<s1->nr;k++)
1143                                 {
1144                                     cmp_int(stdout, buf, i,
1145                                             s1->ival[k], s2->ival[k]);
1146                                 }
1147                                 break;
1148                             case xdr_datatype_large_int:
1149                                 for(k=0;k<s1->nr;k++)
1150                                 {
1151                                     cmp_gmx_large_int(stdout, buf, 
1152                                                       s1->lval[k], s2->lval[k]);
1153                                 }
1154                                 break;
1155                             case xdr_datatype_char:
1156                                 for(k=0;k<s1->nr;k++)
1157                                 {
1158                                     cmp_uc(stdout, buf, i,
1159                                            s1->cval[k], s2->cval[k]);
1160                                 }
1161                                 break;
1162                             case xdr_datatype_string:
1163                                 for(k=0;k<s1->nr;k++)
1164                                 {
1165                                     cmp_str(stdout, buf, i,
1166                                             s1->sval[k], s2->sval[k]);
1167                                 }
1168                                 break;
1169                             default:
1170                                 gmx_incons("Unknown data type!!");
1171                         }
1172                     }
1173                 }
1174             }
1175         }
1176     }
1177 }
1178
1179 void comp_enx(const char *fn1,const char *fn2,real ftol,real abstol,const char *lastener)
1180 {
1181   int        nre,nre1,nre2,block;
1182   ener_file_t in1, in2;
1183   int        i,j,maxener,*ind1,*ind2,*have;
1184   char       buf[256];
1185   gmx_enxnm_t *enm1=NULL,*enm2=NULL;
1186   t_enxframe *fr1,*fr2;
1187   gmx_bool       b1,b2;
1188   
1189   fprintf(stdout,"comparing energy file %s and %s\n\n",fn1,fn2);
1190
1191   in1 = open_enx(fn1,"r");
1192   in2 = open_enx(fn2,"r");
1193   do_enxnms(in1,&nre1,&enm1);
1194   do_enxnms(in2,&nre2,&enm2);
1195   if (nre1 != nre2) {
1196     fprintf(stdout,"There are %d and %d terms in the energy files\n\n",
1197             nre1,nre2);
1198   } else {
1199     fprintf(stdout,"There are %d terms in the energy files\n\n",nre1);
1200   }
1201
1202   snew(ind1,nre1);
1203   snew(ind2,nre2);
1204   snew(have,nre2);
1205   nre = 0;
1206   for(i=0; i<nre1; i++) {
1207     for(j=0; j<nre2; j++) {
1208       if (enernm_equal(enm1[i].name,enm2[j].name)) {
1209         ind1[nre] = i;
1210         ind2[nre] = j;
1211         have[j] = 1;
1212         nre++;
1213         break;
1214       }
1215     }
1216     if (nre == 0 || ind1[nre-1] != i) {
1217       cmp_str(stdout,"enm",i,enm1[i].name,"-");
1218     }
1219   }
1220   for(i=0; i<nre2; i++) {
1221     if (have[i] == 0) {
1222       cmp_str(stdout,"enm",i,"-",enm2[i].name);
1223     }
1224   }
1225
1226   maxener = nre;
1227   for(i=0; i<nre; i++) {
1228     if ((lastener != NULL) && (strstr(enm1[i].name,lastener) != NULL)) {
1229       maxener=i+1;
1230       break;
1231     }
1232   }
1233
1234   fprintf(stdout,"There are %d terms to compare in the energy files\n\n",
1235           maxener);
1236
1237   for(i=0; i<maxener; i++) {
1238     cmp_str(stdout,"unit",i,enm1[ind1[i]].unit,enm2[ind2[i]].unit);
1239   }
1240   
1241   snew(fr1,1);
1242   snew(fr2,1);
1243   do { 
1244     b1 = do_enx(in1,fr1);
1245     b2 = do_enx(in2,fr2);
1246     if (b1 && !b2)
1247       fprintf(stdout,"\nEnd of file on %s but not on %s\n",fn2,fn1);
1248     else if (!b1 && b2) 
1249       fprintf(stdout,"\nEnd of file on %s but not on %s\n",fn1,fn2);
1250     else if (!b1 && !b2)
1251       fprintf(stdout,"\nFiles read successfully\n");
1252     else {
1253       cmp_real(stdout,"t",-1,fr1->t,fr2->t,ftol,abstol);
1254       cmp_int(stdout,"step",-1,fr1->step,fr2->step);
1255       /* We don't want to print the nre mismatch for every frame */
1256       /* cmp_int(stdout,"nre",-1,fr1->nre,fr2->nre); */
1257       if ((fr1->nre >= nre) && (fr2->nre >= nre))
1258         cmp_energies(stdout,fr1->step,fr1->step,fr1->ener,fr2->ener,
1259                      enm1,enm2,ftol,abstol,nre,ind1,ind2,maxener);
1260       /*cmp_disres(fr1,fr2,ftol,abstol);*/
1261       cmp_eblocks(fr1,fr2,ftol,abstol);
1262     }
1263   } while (b1 && b2);
1264     
1265   close_enx(in1);
1266   close_enx(in2);
1267   
1268   free_enxframe(fr2);
1269   sfree(fr2);
1270   free_enxframe(fr1);
1271   sfree(fr1);
1272 }