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