Remove utility/ dependencies on network.h and main.h
[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/fileio/futil.h"
48 #include "sysstuff.h"
49 #include "txtdump.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "names.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "gromacs/fileio/trxio.h"
54 #include "gromacs/fileio/enxio.h"
55 #include "mtop_util.h"
56 #include "gromacs/utility/cstringutil.h"
57
58 static void cmp_int(FILE *fp, const char *s, int index, int i1, int i2)
59 {
60     if (i1 != i2)
61     {
62         if (index != -1)
63         {
64             fprintf(fp, "%s[%d] (%d - %d)\n", s, index, i1, i2);
65         }
66         else
67         {
68             fprintf(fp, "%s (%d - %d)\n", s, i1, i2);
69         }
70     }
71 }
72
73 static void cmp_int64(FILE *fp, const char *s, gmx_int64_t i1, gmx_int64_t i2)
74 {
75     if (i1 != i2)
76     {
77         fprintf(fp, "%s (", s);
78         fprintf(fp, "%"GMX_PRId64, i1);
79         fprintf(fp, " - ");
80         fprintf(fp, "%"GMX_PRId64, i2);
81         fprintf(fp, ")\n");
82     }
83 }
84
85 static void cmp_us(FILE *fp, const char *s, int index, unsigned short i1, unsigned short i2)
86 {
87     if (i1 != i2)
88     {
89         if (index != -1)
90         {
91             fprintf(fp, "%s[%d] (%u - %u)\n", s, index, i1, i2);
92         }
93         else
94         {
95             fprintf(fp, "%s (%u - %u)\n", s, i1, i2);
96         }
97     }
98 }
99
100 static void cmp_uc(FILE *fp, const char *s, int index, unsigned char i1, unsigned char i2)
101 {
102     if (i1 != i2)
103     {
104         if (index != -1)
105         {
106             fprintf(fp, "%s[%d] (%d - %d)\n", s, index, i1, i2);
107         }
108         else
109         {
110             fprintf(fp, "%s (%d - %d)\n", s, i1, i2);
111         }
112     }
113 }
114
115 static gmx_bool cmp_bool(FILE *fp, const char *s, int index, gmx_bool b1, gmx_bool b2)
116 {
117     if (b1)
118     {
119         b1 = 1;
120     }
121     else
122     {
123         b1 = 0;
124     }
125     if (b2)
126     {
127         b2 = 1;
128     }
129     else
130     {
131         b2 = 0;
132     }
133     if (b1 != b2)
134     {
135         if (index != -1)
136         {
137             fprintf(fp, "%s[%d] (%s - %s)\n", s, index,
138                     bool_names[b1], bool_names[b2]);
139         }
140         else
141         {
142             fprintf(fp, "%s (%s - %s)\n", s,
143                     bool_names[b1], bool_names[b2]);
144         }
145     }
146     return b1 && b2;
147 }
148
149 static void cmp_str(FILE *fp, const char *s, int index,
150                     const char *s1, const char *s2)
151 {
152     if (strcmp(s1, s2) != 0)
153     {
154         if (index != -1)
155         {
156             fprintf(fp, "%s[%d] (%s - %s)\n", s, index, s1, s2);
157         }
158         else
159         {
160             fprintf(fp, "%s (%s - %s)\n", s, s1, s2);
161         }
162     }
163 }
164
165 static gmx_bool equal_real(real i1, real i2, real ftol, real abstol)
166 {
167     return ( ( 2*fabs(i1 - i2) <= (fabs(i1) + fabs(i2))*ftol ) || fabs(i1-i2) <= abstol );
168 }
169
170 static gmx_bool equal_float(float i1, float i2, float ftol, float abstol)
171 {
172     return ( ( 2*fabs(i1 - i2) <= (fabs(i1) + fabs(i2))*ftol ) || fabs(i1-i2) <= abstol );
173 }
174
175 static gmx_bool equal_double(double i1, double i2, real ftol, real abstol)
176 {
177     return ( ( 2*fabs(i1 - i2) <= (fabs(i1) + fabs(i2))*ftol ) || fabs(i1-i2) <= abstol );
178 }
179
180 static void
181 cmp_real(FILE *fp, const char *s, int index, real i1, real i2, real ftol, real abstol)
182 {
183     if (!equal_real(i1, i2, ftol, abstol))
184     {
185         if (index != -1)
186         {
187             fprintf(fp, "%s[%2d] (%e - %e)\n", s, index, i1, i2);
188         }
189         else
190         {
191             fprintf(fp, "%s (%e - %e)\n", s, i1, i2);
192         }
193     }
194 }
195
196 static void
197 cmp_float(FILE *fp, const char *s, int index, float i1, float i2, float ftol, float abstol)
198 {
199     if (!equal_float(i1, i2, ftol, abstol))
200     {
201         if (index != -1)
202         {
203             fprintf(fp, "%s[%2d] (%e - %e)\n", s, index, i1, i2);
204         }
205         else
206         {
207             fprintf(fp, "%s (%e - %e)\n", s, i1, i2);
208         }
209     }
210 }
211
212
213
214 static void
215 cmp_double(FILE *fp, const char *s, int index, double i1, double i2, double ftol, double abstol)
216 {
217     if (!equal_double(i1, i2, ftol, abstol))
218     {
219         if (index != -1)
220         {
221             fprintf(fp, "%s[%2d] (%16.9e - %16.9e)\n", s, index, i1, i2);
222         }
223         else
224         {
225             fprintf(fp, "%s (%16.9e - %16.9e)\n", s, i1, i2);
226         }
227     }
228 }
229
230 static void cmp_rvec(FILE *fp, const char *s, int index, rvec i1, rvec i2, real ftol, real abstol)
231 {
232     if (!equal_real(i1[XX], i2[XX], ftol, abstol) ||
233         !equal_real(i1[YY], i2[YY], ftol, abstol) ||
234         !equal_real(i1[ZZ], i2[ZZ], ftol, abstol))
235     {
236         if (index != -1)
237         {
238             fprintf(fp, "%s[%5d] (%12.5e %12.5e %12.5e) - (%12.5e %12.5e %12.5e)\n",
239                     s, index, i1[XX], i1[YY], i1[ZZ], i2[XX], i2[YY], i2[ZZ]);
240         }
241         else
242         {
243             fprintf(fp, "%s (%12.5e %12.5e %12.5e) - (%12.5e %12.5e %12.5e)\n",
244                     s, i1[XX], i1[YY], i1[ZZ], i2[XX], i2[YY], i2[ZZ]);
245         }
246     }
247 }
248
249 static void cmp_ivec(FILE *fp, const char *s, int index, ivec i1, ivec i2)
250 {
251     if ((i1[XX] != i2[XX]) || (i1[YY] != i2[YY]) || (i1[ZZ] != i2[ZZ]))
252     {
253         if (index != -1)
254         {
255             fprintf(fp, "%s[%5d] (%8d,%8d,%8d - %8d,%8d,%8d)\n", s, index,
256                     i1[XX], i1[YY], i1[ZZ], i2[XX], i2[YY], i2[ZZ]);
257         }
258         else
259         {
260             fprintf(fp, "%s (%8d,%8d,%8d - %8d,%8d,%8d)\n", s,
261                     i1[XX], i1[YY], i1[ZZ], i2[XX], i2[YY], i2[ZZ]);
262         }
263     }
264 }
265
266 static void cmp_ilist(FILE *fp, int ftype, t_ilist *il1, t_ilist *il2)
267 {
268     int  i;
269     char buf[256];
270
271     fprintf(fp, "comparing ilist %s\n", interaction_function[ftype].name);
272     sprintf(buf, "%s->nr", interaction_function[ftype].name);
273     cmp_int(fp, buf, -1, il1->nr, il2->nr);
274     sprintf(buf, "%s->iatoms", interaction_function[ftype].name);
275     if (((il1->nr > 0) && (!il1->iatoms)) ||
276         ((il2->nr > 0) && (!il2->iatoms)) ||
277         ((il1->nr != il2->nr)))
278     {
279         fprintf(fp, "Comparing radically different topologies - %s is different\n",
280                 buf);
281     }
282     else
283     {
284         for (i = 0; (i < il1->nr); i++)
285         {
286             cmp_int(fp, buf, i, il1->iatoms[i], il2->iatoms[i]);
287         }
288     }
289 }
290
291 void cmp_iparm(FILE *fp, const char *s, t_functype ft,
292                t_iparams ip1, t_iparams ip2, real ftol, real abstol)
293 {
294     int      i;
295     gmx_bool bDiff;
296
297     bDiff = FALSE;
298     for (i = 0; i < MAXFORCEPARAM && !bDiff; i++)
299     {
300         bDiff = !equal_real(ip1.generic.buf[i], ip2.generic.buf[i], ftol, abstol);
301     }
302     if (bDiff)
303     {
304         fprintf(fp, "%s1: ", s);
305         pr_iparams(fp, ft, &ip1);
306         fprintf(fp, "%s2: ", s);
307         pr_iparams(fp, ft, &ip2);
308     }
309 }
310
311 void cmp_iparm_AB(FILE *fp, const char *s, t_functype ft, t_iparams ip1, real ftol, real abstol)
312 {
313     int      nrfpA, nrfpB, p0, i;
314     gmx_bool bDiff;
315
316     /* Normally the first parameter is perturbable */
317     p0    = 0;
318     nrfpA = interaction_function[ft].nrfpA;
319     nrfpB = interaction_function[ft].nrfpB;
320     if (ft == F_PDIHS)
321     {
322         nrfpB = 2;
323     }
324     else if (interaction_function[ft].flags & IF_TABULATED)
325     {
326         /* For tabulated interactions only the second parameter is perturbable */
327         p0    = 1;
328         nrfpB = 1;
329     }
330     bDiff = FALSE;
331     for (i = 0; i < nrfpB && !bDiff; i++)
332     {
333         bDiff = !equal_real(ip1.generic.buf[p0+i], ip1.generic.buf[nrfpA+i], ftol, abstol);
334     }
335     if (bDiff)
336     {
337         fprintf(fp, "%s: ", s);
338         pr_iparams(fp, ft, &ip1);
339     }
340 }
341
342 static void cmp_idef(FILE *fp, t_idef *id1, t_idef *id2, real ftol, real abstol)
343 {
344     int  i;
345     char buf1[64], buf2[64];
346
347     fprintf(fp, "comparing idef\n");
348     if (id2)
349     {
350         cmp_int(fp, "idef->ntypes", -1, id1->ntypes, id2->ntypes);
351         cmp_int(fp, "idef->atnr",  -1, id1->atnr, id2->atnr);
352         for (i = 0; (i < min(id1->ntypes, id2->ntypes)); i++)
353         {
354             sprintf(buf1, "idef->functype[%d]", i);
355             sprintf(buf2, "idef->iparam[%d]", i);
356             cmp_int(fp, buf1, i, (int)id1->functype[i], (int)id2->functype[i]);
357             cmp_iparm(fp, buf2, id1->functype[i],
358                       id1->iparams[i], id2->iparams[i], ftol, abstol);
359         }
360         cmp_real(fp, "fudgeQQ", -1, id1->fudgeQQ, id2->fudgeQQ, ftol, abstol);
361         for (i = 0; (i < F_NRE); i++)
362         {
363             cmp_ilist(fp, i, &(id1->il[i]), &(id2->il[i]));
364         }
365     }
366     else
367     {
368         for (i = 0; (i < id1->ntypes); i++)
369         {
370             cmp_iparm_AB(fp, "idef->iparam", id1->functype[i], id1->iparams[i], ftol, abstol);
371         }
372     }
373 }
374
375 static void cmp_block(FILE *fp, t_block *b1, t_block *b2, const char *s)
376 {
377     int  i, j, k;
378     char buf[32];
379
380     fprintf(fp, "comparing block %s\n", s);
381     sprintf(buf, "%s.nr", s);
382     cmp_int(fp, buf, -1, b1->nr, b2->nr);
383 }
384
385 static void cmp_blocka(FILE *fp, t_blocka *b1, t_blocka *b2, const char *s)
386 {
387     int  i, j, k;
388     char buf[32];
389
390     fprintf(fp, "comparing blocka %s\n", s);
391     sprintf(buf, "%s.nr", s);
392     cmp_int(fp, buf, -1, b1->nr, b2->nr);
393     sprintf(buf, "%s.nra", s);
394     cmp_int(fp, buf, -1, b1->nra, b2->nra);
395 }
396
397 static void cmp_atom(FILE *fp, int index, t_atom *a1, t_atom *a2, real ftol, real abstol)
398 {
399     int  i;
400     char buf[256];
401
402     if (a2)
403     {
404         cmp_us(fp, "atom.type", index, a1->type, a2->type);
405         cmp_us(fp, "atom.ptype", index, a1->ptype, a2->ptype);
406         cmp_int(fp, "atom.resind", index, a1->resind, a2->resind);
407         cmp_int(fp, "atom.atomnumber", index, a1->atomnumber, a2->atomnumber);
408         cmp_real(fp, "atom.m", index, a1->m, a2->m, ftol, abstol);
409         cmp_real(fp, "atom.q", index, a1->q, a2->q, ftol, abstol);
410         cmp_us(fp, "atom.typeB", index, a1->typeB, a2->typeB);
411         cmp_real(fp, "atom.mB", index, a1->mB, a2->mB, ftol, abstol);
412         cmp_real(fp, "atom.qB", index, a1->qB, a2->qB, ftol, abstol);
413     }
414     else
415     {
416         cmp_us(fp, "atom.type", index, a1->type, a1->typeB);
417         cmp_real(fp, "atom.m", index, a1->m, a1->mB, ftol, abstol);
418         cmp_real(fp, "atom.q", index, a1->q, a1->qB, ftol, abstol);
419     }
420 }
421
422 static void cmp_atoms(FILE *fp, t_atoms *a1, t_atoms *a2, real ftol, real abstol)
423 {
424     int i;
425
426     fprintf(fp, "comparing atoms\n");
427
428     if (a2)
429     {
430         cmp_int(fp, "atoms->nr", -1, a1->nr, a2->nr);
431         for (i = 0; (i < a1->nr); i++)
432         {
433             cmp_atom(fp, i, &(a1->atom[i]), &(a2->atom[i]), ftol, abstol);
434         }
435     }
436     else
437     {
438         for (i = 0; (i < a1->nr); i++)
439         {
440             cmp_atom(fp, i, &(a1->atom[i]), NULL, ftol, abstol);
441         }
442     }
443 }
444
445 static void cmp_top(FILE *fp, t_topology *t1, t_topology *t2, real ftol, real abstol)
446 {
447     int i;
448
449     fprintf(fp, "comparing top\n");
450     if (t2)
451     {
452         cmp_idef(fp, &(t1->idef), &(t2->idef), ftol, abstol);
453         cmp_atoms(fp, &(t1->atoms), &(t2->atoms), ftol, abstol);
454         cmp_block(fp, &t1->cgs, &t2->cgs, "cgs");
455         cmp_block(fp, &t1->mols, &t2->mols, "mols");
456         cmp_blocka(fp, &t1->excls, &t2->excls, "excls");
457     }
458     else
459     {
460         cmp_idef(fp, &(t1->idef), NULL, ftol, abstol);
461         cmp_atoms(fp, &(t1->atoms), NULL, ftol, abstol);
462     }
463 }
464
465 static void cmp_groups(FILE *fp, gmx_groups_t *g0, gmx_groups_t *g1,
466                        int natoms0, int natoms1)
467 {
468     int  i, j, ndiff;
469     char buf[32];
470
471     fprintf(fp, "comparing groups\n");
472
473     for (i = 0; i < egcNR; i++)
474     {
475         sprintf(buf, "grps[%d].nr", i);
476         cmp_int(fp, buf, -1, g0->grps[i].nr, g1->grps[i].nr);
477         if (g0->grps[i].nr == g1->grps[i].nr)
478         {
479             for (j = 0; j < g0->grps[i].nr; j++)
480             {
481                 sprintf(buf, "grps[%d].name[%d]", i, j);
482                 cmp_str(fp, buf, -1,
483                         *g0->grpname[g0->grps[i].nm_ind[j]],
484                         *g1->grpname[g1->grps[i].nm_ind[j]]);
485             }
486         }
487         cmp_int(fp, "ngrpnr", i, g0->ngrpnr[i], g1->ngrpnr[i]);
488         if (g0->ngrpnr[i] == g1->ngrpnr[i] && natoms0 == natoms1 &&
489             (g0->grpnr[i] != NULL || g1->grpnr[i] != NULL))
490         {
491             for (j = 0; j < natoms0; j++)
492             {
493                 cmp_int(fp, gtypes[i], j, ggrpnr(g0, i, j), ggrpnr(g1, i, j));
494             }
495         }
496     }
497     /* We have compared the names in the groups lists,
498      * so we can skip the grpname list comparison.
499      */
500 }
501
502 static void cmp_rvecs(FILE *fp, const char *title, int n, rvec x1[], rvec x2[],
503                       gmx_bool bRMSD, real ftol, real abstol)
504 {
505     int    i, m;
506     double d, ssd;
507
508     if (bRMSD)
509     {
510         ssd = 0;
511         for (i = 0; (i < n); i++)
512         {
513             for (m = 0; m < DIM; m++)
514             {
515                 d    = x1[i][m] - x2[i][m];
516                 ssd += d*d;
517             }
518         }
519         fprintf(fp, "%s RMSD %g\n", title, sqrt(ssd/n));
520     }
521     else
522     {
523         for (i = 0; (i < n); i++)
524         {
525             cmp_rvec(fp, title, i, x1[i], x2[i], ftol, abstol);
526         }
527     }
528 }
529
530
531 /* Similar to cmp_rvecs, but this routine scales the allowed absolute tolerance
532  * by the RMS of the force components of x1.
533  */
534 static void cmp_rvecs_rmstol(FILE *fp, const char *title, int n, rvec x1[], rvec x2[],
535                              real ftol, real abstol)
536 {
537     int    i, m;
538     double d;
539     double ave_x1, rms_x1;
540
541     /* It is tricky to compare real values, in particular forces that
542      * are sums of lots of terms where the final value might be close to 0.0.
543      * To get a reference magnitude we calculate the RMS value of each
544      * component in x1, and then set the allowed absolute tolerance to the
545      * relative tolerance times this RMS magnitude.
546      */
547     ave_x1 = 0.0;
548     for (i = 0; i < n; i++)
549     {
550         for (m = 0; m < DIM; m++)
551         {
552             ave_x1 += x1[i][m];
553         }
554     }
555     ave_x1 /= n*DIM;
556
557     rms_x1 = 0.0;
558     for (i = 0; (i < n); i++)
559     {
560         for (m = 0; m < DIM; m++)
561         {
562             d       = x1[i][m] - ave_x1;
563             rms_x1 += d*d;
564         }
565     }
566     rms_x1 = sqrt(rms_x1/(DIM*n));
567     /* And now do the actual comparision with a hopefully realistic abstol. */
568     for (i = 0; (i < n); i++)
569     {
570         cmp_rvec(fp, title, i, x1[i], x2[i], ftol, abstol*rms_x1);
571     }
572 }
573
574 static void cmp_grpopts(FILE *fp, t_grpopts *opt1, t_grpopts *opt2, real ftol, real abstol)
575 {
576     int  i, j;
577     char buf1[256], buf2[256];
578
579     cmp_int(fp, "inputrec->grpopts.ngtc", -1,  opt1->ngtc, opt2->ngtc);
580     cmp_int(fp, "inputrec->grpopts.ngacc", -1, opt1->ngacc, opt2->ngacc);
581     cmp_int(fp, "inputrec->grpopts.ngfrz", -1, opt1->ngfrz, opt2->ngfrz);
582     cmp_int(fp, "inputrec->grpopts.ngener", -1, opt1->ngener, opt2->ngener);
583     for (i = 0; (i < min(opt1->ngtc, opt2->ngtc)); i++)
584     {
585         cmp_real(fp, "inputrec->grpopts.nrdf", i, opt1->nrdf[i], opt2->nrdf[i], ftol, abstol);
586         cmp_real(fp, "inputrec->grpopts.ref_t", i, opt1->ref_t[i], opt2->ref_t[i], ftol, abstol);
587         cmp_real(fp, "inputrec->grpopts.tau_t", i, opt1->tau_t[i], opt2->tau_t[i], ftol, abstol);
588         cmp_int(fp, "inputrec->grpopts.annealing", i, opt1->annealing[i], opt2->annealing[i]);
589         cmp_int(fp, "inputrec->grpopts.anneal_npoints", i,
590                 opt1->anneal_npoints[i], opt2->anneal_npoints[i]);
591         if (opt1->anneal_npoints[i] == opt2->anneal_npoints[i])
592         {
593             sprintf(buf1, "inputrec->grpopts.anneal_time[%d]", i);
594             sprintf(buf2, "inputrec->grpopts.anneal_temp[%d]", i);
595             for (j = 0; j < opt1->anneal_npoints[i]; j++)
596             {
597                 cmp_real(fp, buf1, j, opt1->anneal_time[i][j], opt2->anneal_time[i][j], ftol, abstol);
598                 cmp_real(fp, buf2, j, opt1->anneal_temp[i][j], opt2->anneal_temp[i][j], ftol, abstol);
599             }
600         }
601     }
602     if (opt1->ngener == opt2->ngener)
603     {
604         for (i = 0; i < opt1->ngener; i++)
605         {
606             for (j = i; j < opt1->ngener; j++)
607             {
608                 sprintf(buf1, "inputrec->grpopts.egp_flags[%d]", i);
609                 cmp_int(fp, buf1, j,
610                         opt1->egp_flags[opt1->ngener*i+j],
611                         opt2->egp_flags[opt1->ngener*i+j]);
612             }
613         }
614     }
615     for (i = 0; (i < min(opt1->ngacc, opt2->ngacc)); i++)
616     {
617         cmp_rvec(fp, "inputrec->grpopts.acc", i, opt1->acc[i], opt2->acc[i], ftol, abstol);
618     }
619     for (i = 0; (i < min(opt1->ngfrz, opt2->ngfrz)); i++)
620     {
621         cmp_ivec(fp, "inputrec->grpopts.nFreeze", i, opt1->nFreeze[i], opt2->nFreeze[i]);
622     }
623 }
624
625 static void cmp_cosines(FILE *fp, const char *s, t_cosines c1[DIM], t_cosines c2[DIM], real ftol, real abstol)
626 {
627     int  i, m;
628     char buf[256];
629
630     for (m = 0; (m < DIM); m++)
631     {
632         sprintf(buf, "inputrec->%s[%d]", s, m);
633         cmp_int(fp, buf, 0, c1->n, c2->n);
634         for (i = 0; (i < min(c1->n, c2->n)); i++)
635         {
636             cmp_real(fp, buf, i, c1->a[i], c2->a[i], ftol, abstol);
637             cmp_real(fp, buf, i, c1->phi[i], c2->phi[i], ftol, abstol);
638         }
639     }
640 }
641 static void cmp_adress(FILE *fp, t_adress *ad1, t_adress *ad2,
642                        real ftol, real abstol)
643 {
644     cmp_int(fp, "ir->adress->type", -1, ad1->type, ad2->type);
645     cmp_real(fp, "ir->adress->const_wf", -1, ad1->const_wf, ad2->const_wf, ftol, abstol);
646     cmp_real(fp, "ir->adress->ex_width", -1, ad1->ex_width, ad2->ex_width, ftol, abstol);
647     cmp_real(fp, "ir->adress->hy_width", -1, ad1->hy_width, ad2->hy_width, ftol, abstol);
648     cmp_int(fp, "ir->adress->icor", -1, ad1->icor, ad2->icor);
649     cmp_int(fp, "ir->adress->site", -1, ad1->site, ad2->site);
650     cmp_rvec(fp, "ir->adress->refs", -1, ad1->refs, ad2->refs, ftol, abstol);
651     cmp_real(fp, "ir->adress->ex_forcecap", -1, ad1->ex_forcecap, ad2->ex_forcecap, ftol, abstol);
652 }
653
654 static void cmp_pull(FILE *fp)
655 {
656     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");
657 }
658
659 static void cmp_simtempvals(FILE *fp, t_simtemp *simtemp1, t_simtemp *simtemp2, int n_lambda, real ftol, real abstol)
660 {
661     int i;
662     cmp_int(fp, "inputrec->simtempvals->eSimTempScale", -1, simtemp1->eSimTempScale, simtemp2->eSimTempScale);
663     cmp_real(fp, "inputrec->simtempvals->simtemp_high", -1, simtemp1->simtemp_high, simtemp2->simtemp_high, ftol, abstol);
664     cmp_real(fp, "inputrec->simtempvals->simtemp_low", -1, simtemp1->simtemp_low, simtemp2->simtemp_low, ftol, abstol);
665     for (i = 0; i < n_lambda; i++)
666     {
667         cmp_real(fp, "inputrec->simtempvals->temperatures", -1, simtemp1->temperatures[i], simtemp2->temperatures[i], ftol, abstol);
668     }
669 }
670
671 static void cmp_expandedvals(FILE *fp, t_expanded *expand1, t_expanded *expand2, int n_lambda, real ftol, real abstol)
672 {
673     int i;
674
675     cmp_bool(fp, "inputrec->fepvals->bInit_weights", -1, expand1->bInit_weights, expand2->bInit_weights);
676     cmp_bool(fp, "inputrec->fepvals->bWLoneovert", -1, expand1->bWLoneovert, expand2->bWLoneovert);
677
678     for (i = 0; i < n_lambda; i++)
679     {
680         cmp_real(fp, "inputrec->expandedvals->init_lambda_weights", -1,
681                  expand1->init_lambda_weights[i], expand2->init_lambda_weights[i], ftol, abstol);
682     }
683
684     cmp_int(fp, "inputrec->expandedvals->lambda-stats", -1, expand1->elamstats, expand2->elamstats);
685     cmp_int(fp, "inputrec->expandedvals->lambda-mc-move", -1, expand1->elmcmove, expand2->elmcmove);
686     cmp_int(fp, "inputrec->expandedvals->lmc-repeats", -1, expand1->lmc_repeats, expand2->lmc_repeats);
687     cmp_int(fp, "inputrec->expandedvals->lmc-gibbsdelta", -1, expand1->gibbsdeltalam, expand2->gibbsdeltalam);
688     cmp_int(fp, "inputrec->expandedvals->lmc-forced-nstart", -1, expand1->lmc_forced_nstart, expand2->lmc_forced_nstart);
689     cmp_int(fp, "inputrec->expandedvals->lambda-weights-equil", -1, expand1->elmceq, expand2->elmceq);
690     cmp_int(fp, "inputrec->expandedvals->,weight-equil-number-all-lambda", -1, expand1->equil_n_at_lam, expand2->equil_n_at_lam);
691     cmp_int(fp, "inputrec->expandedvals->weight-equil-number-samples", -1, expand1->equil_samples, expand2->equil_samples);
692     cmp_int(fp, "inputrec->expandedvals->weight-equil-number-steps", -1, expand1->equil_steps, expand2->equil_steps);
693     cmp_real(fp, "inputrec->expandedvals->weight-equil-wl-delta", -1, expand1->equil_wl_delta, expand2->equil_wl_delta, ftol, abstol);
694     cmp_real(fp, "inputrec->expandedvals->weight-equil-count-ratio", -1, expand1->equil_ratio, expand2->equil_ratio, ftol, abstol);
695     cmp_bool(fp, "inputrec->expandedvals->symmetrized-transition-matrix", -1, expand1->bSymmetrizedTMatrix, expand2->bSymmetrizedTMatrix);
696     cmp_int(fp, "inputrec->expandedvals->nstTij", -1, expand1->nstTij, expand2->nstTij);
697     cmp_int(fp, "inputrec->expandedvals->mininum-var-min", -1, expand1->minvarmin, expand2->minvarmin); /*default is reasonable */
698     cmp_int(fp, "inputrec->expandedvals->weight-c-range", -1, expand1->c_range, expand2->c_range);      /* default is just C=0 */
699     cmp_real(fp, "inputrec->expandedvals->wl-scale", -1, expand1->wl_scale, expand2->wl_scale, ftol, abstol);
700     cmp_real(fp, "inputrec->expandedvals->init-wl-delta", -1, expand1->init_wl_delta, expand2->init_wl_delta, ftol, abstol);
701     cmp_real(fp, "inputrec->expandedvals->wl-ratio", -1, expand1->wl_ratio, expand2->wl_ratio, ftol, abstol);
702     cmp_int(fp, "inputrec->expandedvals->nstexpanded", -1, expand1->nstexpanded, expand2->nstexpanded);
703     cmp_int(fp, "inputrec->expandedvals->lmc-seed", -1, expand1->lmc_seed, expand2->lmc_seed);
704     cmp_real(fp, "inputrec->expandedvals->mc-temperature", -1, expand1->mc_temp, expand2->mc_temp, ftol, abstol);
705 }
706
707 static void cmp_fepvals(FILE *fp, t_lambda *fep1, t_lambda *fep2, real ftol, real abstol)
708 {
709     int i, j;
710     cmp_int(fp, "inputrec->nstdhdl", -1, fep1->nstdhdl, fep2->nstdhdl);
711     cmp_double(fp, "inputrec->fepvals->init_fep_state", -1, fep1->init_fep_state, fep2->init_fep_state, ftol, abstol);
712     cmp_double(fp, "inputrec->fepvals->delta_lambda", -1, fep1->delta_lambda, fep2->delta_lambda, ftol, abstol);
713     cmp_int(fp, "inputrec->fepvals->n_lambda", -1, fep1->n_lambda, fep2->n_lambda);
714     for (i = 0; i < efptNR; i++)
715     {
716         for (j = 0; j < min(fep1->n_lambda, fep2->n_lambda); j++)
717         {
718             cmp_double(fp, "inputrec->fepvals->all_lambda", -1, fep1->all_lambda[i][j], fep2->all_lambda[i][j], ftol, abstol);
719         }
720     }
721     cmp_int(fp, "inputrec->fepvals->lambda_neighbors", 1, fep1->lambda_neighbors,
722             fep2->lambda_neighbors);
723     cmp_real(fp, "inputrec->fepvals->sc_alpha", -1, fep1->sc_alpha, fep2->sc_alpha, ftol, abstol);
724     cmp_int(fp, "inputrec->fepvals->sc_power", -1, fep1->sc_power, fep2->sc_power);
725     cmp_real(fp, "inputrec->fepvals->sc_r_power", -1, fep1->sc_r_power, fep2->sc_r_power, ftol, abstol);
726     cmp_real(fp, "inputrec->fepvals->sc_sigma", -1, fep1->sc_sigma, fep2->sc_sigma, ftol, abstol);
727     cmp_bool(fp, "inputrec->fepvals->bPrintEnergy", -1, fep1->bPrintEnergy, fep1->bPrintEnergy);
728     cmp_bool(fp, "inputrec->fepvals->bScCoul", -1, fep1->bScCoul, fep1->bScCoul);
729     cmp_int(fp, "inputrec->separate_dhdl_file", -1, fep1->separate_dhdl_file, fep2->separate_dhdl_file);
730     cmp_int(fp, "inputrec->dhdl_derivatives", -1, fep1->dhdl_derivatives, fep2->dhdl_derivatives);
731     cmp_int(fp, "inputrec->dh_hist_size", -1, fep1->dh_hist_size, fep2->dh_hist_size);
732     cmp_double(fp, "inputrec->dh_hist_spacing", -1, fep1->dh_hist_spacing, fep2->dh_hist_spacing, ftol, abstol);
733 }
734
735 static void cmp_inputrec(FILE *fp, t_inputrec *ir1, t_inputrec *ir2, real ftol, real abstol)
736 {
737     fprintf(fp, "comparing inputrec\n");
738
739     /* gcc 2.96 doesnt like these defines at all, but issues a huge list
740      * of warnings. Maybe it will change in future versions, but for the
741      * moment I've spelled them out instead. /EL 000820
742      * #define CIB(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
743      * #define CII(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
744      * #define CIR(s) cmp_real(fp,"inputrec->"#s,0,ir1->##s,ir2->##s,ftol)
745      */
746     cmp_int(fp, "inputrec->eI", -1, ir1->eI, ir2->eI);
747     cmp_int64(fp, "inputrec->nsteps", ir1->nsteps, ir2->nsteps);
748     cmp_int64(fp, "inputrec->init_step", ir1->init_step, ir2->init_step);
749     cmp_int(fp, "inputrec->simulation_part", -1, ir1->simulation_part, ir2->simulation_part);
750     cmp_int(fp, "inputrec->ePBC", -1, ir1->ePBC, ir2->ePBC);
751     cmp_int(fp, "inputrec->bPeriodicMols", -1, ir1->bPeriodicMols, ir2->bPeriodicMols);
752     cmp_int(fp, "inputrec->cutoff_scheme", -1, ir1->cutoff_scheme, ir2->cutoff_scheme);
753     cmp_int(fp, "inputrec->ns_type", -1, ir1->ns_type, ir2->ns_type);
754     cmp_int(fp, "inputrec->nstlist", -1, ir1->nstlist, ir2->nstlist);
755     cmp_int(fp, "inputrec->ndelta", -1, ir1->ndelta, ir2->ndelta);
756     cmp_int(fp, "inputrec->nstcomm", -1, ir1->nstcomm, ir2->nstcomm);
757     cmp_int(fp, "inputrec->comm_mode", -1, ir1->comm_mode, ir2->comm_mode);
758     cmp_int(fp, "inputrec->nstcheckpoint", -1, ir1->nstcheckpoint, ir2->nstcheckpoint);
759     cmp_int(fp, "inputrec->nstlog", -1, ir1->nstlog, ir2->nstlog);
760     cmp_int(fp, "inputrec->nstxout", -1, ir1->nstxout, ir2->nstxout);
761     cmp_int(fp, "inputrec->nstvout", -1, ir1->nstvout, ir2->nstvout);
762     cmp_int(fp, "inputrec->nstfout", -1, ir1->nstfout, ir2->nstfout);
763     cmp_int(fp, "inputrec->nstcalcenergy", -1, ir1->nstcalcenergy, ir2->nstcalcenergy);
764     cmp_int(fp, "inputrec->nstenergy", -1, ir1->nstenergy, ir2->nstenergy);
765     cmp_int(fp, "inputrec->nstxout_compressed", -1, ir1->nstxout_compressed, ir2->nstxout_compressed);
766     cmp_double(fp, "inputrec->init_t", -1, ir1->init_t, ir2->init_t, ftol, abstol);
767     cmp_double(fp, "inputrec->delta_t", -1, ir1->delta_t, ir2->delta_t, ftol, abstol);
768     cmp_real(fp, "inputrec->x_compression_precision", -1, ir1->x_compression_precision, ir2->x_compression_precision, ftol, abstol);
769     cmp_real(fp, "inputrec->fourierspacing", -1, ir1->fourier_spacing, ir2->fourier_spacing, ftol, abstol);
770     cmp_int(fp, "inputrec->nkx", -1, ir1->nkx, ir2->nkx);
771     cmp_int(fp, "inputrec->nky", -1, ir1->nky, ir2->nky);
772     cmp_int(fp, "inputrec->nkz", -1, ir1->nkz, ir2->nkz);
773     cmp_int(fp, "inputrec->pme_order", -1, ir1->pme_order, ir2->pme_order);
774     cmp_real(fp, "inputrec->ewald_rtol", -1, ir1->ewald_rtol, ir2->ewald_rtol, ftol, abstol);
775     cmp_int(fp, "inputrec->ewald_geometry", -1, ir1->ewald_geometry, ir2->ewald_geometry);
776     cmp_real(fp, "inputrec->epsilon_surface", -1, ir1->epsilon_surface, ir2->epsilon_surface, ftol, abstol);
777     cmp_int(fp, "inputrec->bOptFFT", -1, ir1->bOptFFT, ir2->bOptFFT);
778     cmp_int(fp, "inputrec->bContinuation", -1, ir1->bContinuation, ir2->bContinuation);
779     cmp_int(fp, "inputrec->bShakeSOR", -1, ir1->bShakeSOR, ir2->bShakeSOR);
780     cmp_int(fp, "inputrec->etc", -1, ir1->etc, ir2->etc);
781     cmp_int(fp, "inputrec->bPrintNHChains", -1, ir1->bPrintNHChains, ir2->bPrintNHChains);
782     cmp_int(fp, "inputrec->epc", -1, ir1->epc, ir2->epc);
783     cmp_int(fp, "inputrec->epct", -1, ir1->epct, ir2->epct);
784     cmp_real(fp, "inputrec->tau_p", -1, ir1->tau_p, ir2->tau_p, ftol, abstol);
785     cmp_rvec(fp, "inputrec->ref_p(x)", -1, ir1->ref_p[XX], ir2->ref_p[XX], ftol, abstol);
786     cmp_rvec(fp, "inputrec->ref_p(y)", -1, ir1->ref_p[YY], ir2->ref_p[YY], ftol, abstol);
787     cmp_rvec(fp, "inputrec->ref_p(z)", -1, ir1->ref_p[ZZ], ir2->ref_p[ZZ], ftol, abstol);
788     cmp_rvec(fp, "inputrec->compress(x)", -1, ir1->compress[XX], ir2->compress[XX], ftol, abstol);
789     cmp_rvec(fp, "inputrec->compress(y)", -1, ir1->compress[YY], ir2->compress[YY], ftol, abstol);
790     cmp_rvec(fp, "inputrec->compress(z)", -1, ir1->compress[ZZ], ir2->compress[ZZ], ftol, abstol);
791     cmp_int(fp, "refcoord_scaling", -1, ir1->refcoord_scaling, ir2->refcoord_scaling);
792     cmp_rvec(fp, "inputrec->posres_com", -1, ir1->posres_com, ir2->posres_com, ftol, abstol);
793     cmp_rvec(fp, "inputrec->posres_comB", -1, ir1->posres_comB, ir2->posres_comB, ftol, abstol);
794     cmp_real(fp, "inputrec->verletbuf_tol", -1, ir1->verletbuf_tol, ir2->verletbuf_tol, ftol, abstol);
795     cmp_real(fp, "inputrec->rlist", -1, ir1->rlist, ir2->rlist, ftol, abstol);
796     cmp_real(fp, "inputrec->rlistlong", -1, ir1->rlistlong, ir2->rlistlong, ftol, abstol);
797     cmp_int(fp, "inputrec->nstcalclr", -1, ir1->nstcalclr, ir2->nstcalclr);
798     cmp_real(fp, "inputrec->rtpi", -1, ir1->rtpi, ir2->rtpi, ftol, abstol);
799     cmp_int(fp, "inputrec->coulombtype", -1, ir1->coulombtype, ir2->coulombtype);
800     cmp_int(fp, "inputrec->coulomb_modifier", -1, ir1->coulomb_modifier, ir2->coulomb_modifier);
801     cmp_real(fp, "inputrec->rcoulomb_switch", -1, ir1->rcoulomb_switch, ir2->rcoulomb_switch, ftol, abstol);
802     cmp_real(fp, "inputrec->rcoulomb", -1, ir1->rcoulomb, ir2->rcoulomb, ftol, abstol);
803     cmp_int(fp, "inputrec->vdwtype", -1, ir1->vdwtype, ir2->vdwtype);
804     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);
805     cmp_real(fp, "inputrec->rvdw", -1, ir1->rvdw, ir2->rvdw, ftol, abstol);
806     cmp_real(fp, "inputrec->epsilon_r", -1, ir1->epsilon_r, ir2->epsilon_r, ftol, abstol);
807     cmp_real(fp, "inputrec->epsilon_rf", -1, ir1->epsilon_rf, ir2->epsilon_rf, ftol, abstol);
808     cmp_real(fp, "inputrec->tabext", -1, ir1->tabext, ir2->tabext, ftol, abstol);
809     cmp_int(fp, "inputrec->implicit_solvent", -1, ir1->implicit_solvent, ir2->implicit_solvent);
810     cmp_int(fp, "inputrec->gb_algorithm", -1, ir1->gb_algorithm, ir2->gb_algorithm);
811     cmp_int(fp, "inputrec->nstgbradii", -1, ir1->nstgbradii, ir2->nstgbradii);
812     cmp_real(fp, "inputrec->rgbradii", -1, ir1->rgbradii, ir2->rgbradii, ftol, abstol);
813     cmp_real(fp, "inputrec->gb_saltconc", -1, ir1->gb_saltconc, ir2->gb_saltconc, ftol, abstol);
814     cmp_real(fp, "inputrec->gb_epsilon_solvent", -1, ir1->gb_epsilon_solvent, ir2->gb_epsilon_solvent, ftol, abstol);
815     cmp_real(fp, "inputrec->gb_obc_alpha", -1, ir1->gb_obc_alpha, ir2->gb_obc_alpha, ftol, abstol);
816     cmp_real(fp, "inputrec->gb_obc_beta", -1, ir1->gb_obc_beta, ir2->gb_obc_beta, ftol, abstol);
817     cmp_real(fp, "inputrec->gb_obc_gamma", -1, ir1->gb_obc_gamma, ir2->gb_obc_gamma, ftol, abstol);
818     cmp_real(fp, "inputrec->gb_dielectric_offset", -1, ir1->gb_dielectric_offset, ir2->gb_dielectric_offset, ftol, abstol);
819     cmp_int(fp, "inputrec->sa_algorithm", -1, ir1->sa_algorithm, ir2->sa_algorithm);
820     cmp_real(fp, "inputrec->sa_surface_tension", -1, ir1->sa_surface_tension, ir2->sa_surface_tension, ftol, abstol);
821
822     cmp_int(fp, "inputrec->eDispCorr", -1, ir1->eDispCorr, ir2->eDispCorr);
823     cmp_real(fp, "inputrec->shake_tol", -1, ir1->shake_tol, ir2->shake_tol, ftol, abstol);
824     cmp_int(fp, "inputrec->efep", -1, ir1->efep, ir2->efep);
825     cmp_fepvals(fp, ir1->fepvals, ir2->fepvals, ftol, abstol);
826     cmp_int(fp, "inputrec->bSimTemp", -1, ir1->bSimTemp, ir2->bSimTemp);
827     if ((ir1->bSimTemp == ir2->bSimTemp) && (ir1->bSimTemp))
828     {
829         cmp_simtempvals(fp, ir1->simtempvals, ir2->simtempvals, min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda), ftol, abstol);
830     }
831     cmp_int(fp, "inputrec->bExpanded", -1, ir1->bExpanded, ir2->bExpanded);
832     if ((ir1->bExpanded == ir2->bExpanded) && (ir1->bExpanded))
833     {
834         cmp_expandedvals(fp, ir1->expandedvals, ir2->expandedvals, min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda), ftol, abstol);
835     }
836     cmp_int(fp, "inputrec->nwall", -1, ir1->nwall, ir2->nwall);
837     cmp_int(fp, "inputrec->wall_type", -1, ir1->wall_type, ir2->wall_type);
838     cmp_int(fp, "inputrec->wall_atomtype[0]", -1, ir1->wall_atomtype[0], ir2->wall_atomtype[0]);
839     cmp_int(fp, "inputrec->wall_atomtype[1]", -1, ir1->wall_atomtype[1], ir2->wall_atomtype[1]);
840     cmp_real(fp, "inputrec->wall_density[0]", -1, ir1->wall_density[0], ir2->wall_density[0], ftol, abstol);
841     cmp_real(fp, "inputrec->wall_density[1]", -1, ir1->wall_density[1], ir2->wall_density[1], ftol, abstol);
842     cmp_real(fp, "inputrec->wall_ewald_zfac", -1, ir1->wall_ewald_zfac, ir2->wall_ewald_zfac, ftol, abstol);
843
844     cmp_int(fp, "inputrec->ePull", -1, ir1->ePull, ir2->ePull);
845     if (ir1->ePull == ir2->ePull && ir1->ePull != epullNO)
846     {
847         cmp_pull(fp);
848     }
849
850     cmp_int(fp, "inputrec->eDisre", -1, ir1->eDisre, ir2->eDisre);
851     cmp_real(fp, "inputrec->dr_fc", -1, ir1->dr_fc, ir2->dr_fc, ftol, abstol);
852     cmp_int(fp, "inputrec->eDisreWeighting", -1, ir1->eDisreWeighting, ir2->eDisreWeighting);
853     cmp_int(fp, "inputrec->bDisreMixed", -1, ir1->bDisreMixed, ir2->bDisreMixed);
854     cmp_int(fp, "inputrec->nstdisreout", -1, ir1->nstdisreout, ir2->nstdisreout);
855     cmp_real(fp, "inputrec->dr_tau", -1, ir1->dr_tau, ir2->dr_tau, ftol, abstol);
856     cmp_real(fp, "inputrec->orires_fc", -1, ir1->orires_fc, ir2->orires_fc, ftol, abstol);
857     cmp_real(fp, "inputrec->orires_tau", -1, ir1->orires_tau, ir2->orires_tau, ftol, abstol);
858     cmp_int(fp, "inputrec->nstorireout", -1, ir1->nstorireout, ir2->nstorireout);
859     cmp_real(fp, "inputrec->dihre_fc", -1, ir1->dihre_fc, ir2->dihre_fc, ftol, abstol);
860     cmp_real(fp, "inputrec->em_stepsize", -1, ir1->em_stepsize, ir2->em_stepsize, ftol, abstol);
861     cmp_real(fp, "inputrec->em_tol", -1, ir1->em_tol, ir2->em_tol, ftol, abstol);
862     cmp_int(fp, "inputrec->niter", -1, ir1->niter, ir2->niter);
863     cmp_real(fp, "inputrec->fc_stepsize", -1, ir1->fc_stepsize, ir2->fc_stepsize, ftol, abstol);
864     cmp_int(fp, "inputrec->nstcgsteep", -1, ir1->nstcgsteep, ir2->nstcgsteep);
865     cmp_int(fp, "inputrec->nbfgscorr", 0, ir1->nbfgscorr, ir2->nbfgscorr);
866     cmp_int(fp, "inputrec->eConstrAlg", -1, ir1->eConstrAlg, ir2->eConstrAlg);
867     cmp_int(fp, "inputrec->nProjOrder", -1, ir1->nProjOrder, ir2->nProjOrder);
868     cmp_real(fp, "inputrec->LincsWarnAngle", -1, ir1->LincsWarnAngle, ir2->LincsWarnAngle, ftol, abstol);
869     cmp_int(fp, "inputrec->nLincsIter", -1, ir1->nLincsIter, ir2->nLincsIter);
870     cmp_real(fp, "inputrec->bd_fric", -1, ir1->bd_fric, ir2->bd_fric, ftol, abstol);
871     cmp_int64(fp, "inputrec->ld_seed", ir1->ld_seed, ir2->ld_seed);
872     cmp_real(fp, "inputrec->cos_accel", -1, ir1->cos_accel, ir2->cos_accel, ftol, abstol);
873     cmp_rvec(fp, "inputrec->deform(a)", -1, ir1->deform[XX], ir2->deform[XX], ftol, abstol);
874     cmp_rvec(fp, "inputrec->deform(b)", -1, ir1->deform[YY], ir2->deform[YY], ftol, abstol);
875     cmp_rvec(fp, "inputrec->deform(c)", -1, ir1->deform[ZZ], ir2->deform[ZZ], ftol, abstol);
876
877
878     cmp_bool(fp, "ir->bAdress->type", -1, ir1->bAdress, ir2->bAdress);
879     if (ir1->bAdress && ir2->bAdress)
880     {
881         cmp_adress(fp, ir1->adress, ir2->adress, ftol, abstol);
882     }
883
884     cmp_int(fp, "inputrec->userint1", -1, ir1->userint1, ir2->userint1);
885     cmp_int(fp, "inputrec->userint2", -1, ir1->userint2, ir2->userint2);
886     cmp_int(fp, "inputrec->userint3", -1, ir1->userint3, ir2->userint3);
887     cmp_int(fp, "inputrec->userint4", -1, ir1->userint4, ir2->userint4);
888     cmp_real(fp, "inputrec->userreal1", -1, ir1->userreal1, ir2->userreal1, ftol, abstol);
889     cmp_real(fp, "inputrec->userreal2", -1, ir1->userreal2, ir2->userreal2, ftol, abstol);
890     cmp_real(fp, "inputrec->userreal3", -1, ir1->userreal3, ir2->userreal3, ftol, abstol);
891     cmp_real(fp, "inputrec->userreal4", -1, ir1->userreal4, ir2->userreal4, ftol, abstol);
892     cmp_grpopts(fp, &(ir1->opts), &(ir2->opts), ftol, abstol);
893     cmp_cosines(fp, "ex", ir1->ex, ir2->ex, ftol, abstol);
894     cmp_cosines(fp, "et", ir1->et, ir2->et, ftol, abstol);
895 }
896
897 static void comp_pull_AB(FILE *fp, t_pull *pull, real ftol, real abstol)
898 {
899     int i;
900
901     for (i = 0; i < pull->ncoord; i++)
902     {
903         fprintf(fp, "comparing pull coord %d\n", i);
904         cmp_real(fp, "pull-coord->k", -1, pull->coord[i].k, pull->coord[i].kB, ftol, abstol);
905     }
906 }
907
908 static void comp_state(t_state *st1, t_state *st2,
909                        gmx_bool bRMSD, real ftol, real abstol)
910 {
911     int i, j, nc;
912
913     fprintf(stdout, "comparing flags\n");
914     cmp_int(stdout, "flags", -1, st1->flags, st2->flags);
915     fprintf(stdout, "comparing box\n");
916     cmp_rvecs(stdout, "box", DIM, st1->box, st2->box, FALSE, ftol, abstol);
917     fprintf(stdout, "comparing box_rel\n");
918     cmp_rvecs(stdout, "box_rel", DIM, st1->box_rel, st2->box_rel, FALSE, ftol, abstol);
919     fprintf(stdout, "comparing boxv\n");
920     cmp_rvecs(stdout, "boxv", DIM, st1->boxv, st2->boxv, FALSE, ftol, abstol);
921     if (st1->flags & (1<<estSVIR_PREV))
922     {
923         fprintf(stdout, "comparing shake vir_prev\n");
924         cmp_rvecs_rmstol(stdout, "svir_prev", DIM, st1->svir_prev, st2->svir_prev, ftol, abstol);
925     }
926     if (st1->flags & (1<<estFVIR_PREV))
927     {
928         fprintf(stdout, "comparing force vir_prev\n");
929         cmp_rvecs_rmstol(stdout, "fvir_prev", DIM, st1->fvir_prev, st2->fvir_prev, ftol, abstol);
930     }
931     if (st1->flags & (1<<estPRES_PREV))
932     {
933         fprintf(stdout, "comparing prev_pres\n");
934         cmp_rvecs_rmstol(stdout, "pres_prev", DIM, st1->pres_prev, st2->pres_prev, ftol, abstol);
935     }
936     cmp_int(stdout, "ngtc", -1, st1->ngtc, st2->ngtc);
937     cmp_int(stdout, "nhchainlength", -1, st1->nhchainlength, st2->nhchainlength);
938     if (st1->ngtc == st2->ngtc && st1->nhchainlength == st2->nhchainlength)
939     {
940         for (i = 0; i < st1->ngtc; i++)
941         {
942             nc = i*st1->nhchainlength;
943             for (j = 0; j < nc; j++)
944             {
945                 cmp_real(stdout, "nosehoover_xi",
946                          i, st1->nosehoover_xi[nc+j], st2->nosehoover_xi[nc+j], ftol, abstol);
947             }
948         }
949     }
950     cmp_int(stdout, "nnhpres", -1, st1->nnhpres, st2->nnhpres);
951     if (st1->nnhpres == st2->nnhpres && st1->nhchainlength == st2->nhchainlength)
952     {
953         for (i = 0; i < st1->nnhpres; i++)
954         {
955             nc = i*st1->nhchainlength;
956             for (j = 0; j < nc; j++)
957             {
958                 cmp_real(stdout, "nosehoover_xi",
959                          i, st1->nhpres_xi[nc+j], st2->nhpres_xi[nc+j], ftol, abstol);
960             }
961         }
962     }
963
964     cmp_int(stdout, "natoms", -1, st1->natoms, st2->natoms);
965     if (st1->natoms == st2->natoms)
966     {
967         if ((st1->flags & (1<<estX)) && (st2->flags & (1<<estX)))
968         {
969             fprintf(stdout, "comparing x\n");
970             cmp_rvecs(stdout, "x", st1->natoms, st1->x, st2->x, bRMSD, ftol, abstol);
971         }
972         if ((st1->flags & (1<<estV)) && (st2->flags & (1<<estV)))
973         {
974             fprintf(stdout, "comparing v\n");
975             cmp_rvecs(stdout, "v", st1->natoms, st1->v, st2->v, bRMSD, ftol, abstol);
976         }
977     }
978 }
979
980 void comp_tpx(const char *fn1, const char *fn2,
981               gmx_bool bRMSD, real ftol, real abstol)
982 {
983     const char  *ff[2];
984     t_tpxheader  sh[2];
985     t_inputrec   ir[2];
986     t_state      state[2];
987     gmx_mtop_t   mtop[2];
988     t_topology   top[2];
989     int          i;
990
991     ff[0] = fn1;
992     ff[1] = fn2;
993     for (i = 0; i < (fn2 ? 2 : 1); i++)
994     {
995         read_tpx_state(ff[i], &(ir[i]), &state[i], NULL, &(mtop[i]));
996     }
997     if (fn2)
998     {
999         cmp_inputrec(stdout, &ir[0], &ir[1], ftol, abstol);
1000         /* Convert gmx_mtop_t to t_topology.
1001          * We should implement direct mtop comparison,
1002          * but it might be useful to keep t_topology comparison as an option.
1003          */
1004         top[0] = gmx_mtop_t_to_t_topology(&mtop[0]);
1005         top[1] = gmx_mtop_t_to_t_topology(&mtop[1]);
1006         cmp_top(stdout, &top[0], &top[1], ftol, abstol);
1007         cmp_groups(stdout, &mtop[0].groups, &mtop[1].groups,
1008                    mtop[0].natoms, mtop[1].natoms);
1009         comp_state(&state[0], &state[1], bRMSD, ftol, abstol);
1010     }
1011     else
1012     {
1013         if (ir[0].efep == efepNO)
1014         {
1015             fprintf(stdout, "inputrec->efep = %s\n", efep_names[ir[0].efep]);
1016         }
1017         else
1018         {
1019             if (ir[0].ePull != epullNO)
1020             {
1021                 comp_pull_AB(stdout, ir->pull, ftol, abstol);
1022             }
1023             /* Convert gmx_mtop_t to t_topology.
1024              * We should implement direct mtop comparison,
1025              * but it might be useful to keep t_topology comparison as an option.
1026              */
1027             top[0] = gmx_mtop_t_to_t_topology(&mtop[0]);
1028             cmp_top(stdout, &top[0], NULL, ftol, abstol);
1029         }
1030     }
1031 }
1032
1033 void comp_frame(FILE *fp, t_trxframe *fr1, t_trxframe *fr2,
1034                 gmx_bool bRMSD, real ftol, real abstol)
1035 {
1036     fprintf(fp, "\n");
1037     cmp_int(fp, "flags", -1, fr1->flags, fr2->flags);
1038     cmp_int(fp, "not_ok", -1, fr1->not_ok, fr2->not_ok);
1039     cmp_int(fp, "natoms", -1, fr1->natoms, fr2->natoms);
1040     cmp_real(fp, "t0", -1, fr1->t0, fr2->t0, ftol, abstol);
1041     if (cmp_bool(fp, "bTitle", -1, fr1->bTitle, fr2->bTitle))
1042     {
1043         cmp_str(fp, "title", -1, fr1->title, fr2->title);
1044     }
1045     if (cmp_bool(fp, "bStep", -1, fr1->bStep, fr2->bStep))
1046     {
1047         cmp_int(fp, "step", -1, fr1->step, fr2->step);
1048     }
1049     cmp_int(fp, "step", -1, fr1->step, fr2->step);
1050     if (cmp_bool(fp, "bTime", -1, fr1->bTime, fr2->bTime))
1051     {
1052         cmp_real(fp, "time", -1, fr1->time, fr2->time, ftol, abstol);
1053     }
1054     if (cmp_bool(fp, "bLambda", -1, fr1->bLambda, fr2->bLambda))
1055     {
1056         cmp_real(fp, "lambda", -1, fr1->lambda, fr2->lambda, ftol, abstol);
1057     }
1058     if (cmp_bool(fp, "bAtoms", -1, fr1->bAtoms, fr2->bAtoms))
1059     {
1060         cmp_atoms(fp, fr1->atoms, fr2->atoms, ftol, abstol);
1061     }
1062     if (cmp_bool(fp, "bPrec", -1, fr1->bPrec, fr2->bPrec))
1063     {
1064         cmp_real(fp, "prec", -1, fr1->prec, fr2->prec, ftol, abstol);
1065     }
1066     if (cmp_bool(fp, "bX", -1, fr1->bX, fr2->bX))
1067     {
1068         cmp_rvecs(fp, "x", min(fr1->natoms, fr2->natoms), fr1->x, fr2->x, bRMSD, ftol, abstol);
1069     }
1070     if (cmp_bool(fp, "bV", -1, fr1->bV, fr2->bV))
1071     {
1072         cmp_rvecs(fp, "v", min(fr1->natoms, fr2->natoms), fr1->v, fr2->v, bRMSD, ftol, abstol);
1073     }
1074     if (cmp_bool(fp, "bF", -1, fr1->bF, fr2->bF))
1075     {
1076         cmp_rvecs_rmstol(fp, "f", min(fr1->natoms, fr2->natoms), fr1->f, fr2->f, ftol, abstol);
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 }