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