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