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