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