e7e7d2eccd91e55dd66944bb245d9646baf51102
[alexxy/gromacs.git] / src / gromacs / topology / topology.cpp
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,2015,2016,2017,2018,2019, 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 #include "gmxpre.h"
38
39 #include "topology.h"
40
41 #include <cstdio>
42
43 #include <algorithm>
44
45 #include "gromacs/math/vecdump.h"
46 #include "gromacs/topology/atoms.h"
47 #include "gromacs/topology/idef.h"
48 #include "gromacs/topology/ifunc.h"
49 #include "gromacs/topology/symtab.h"
50 #include "gromacs/utility/arrayref.h"
51 #include "gromacs/utility/compare.h"
52 #include "gromacs/utility/gmxassert.h"
53 #include "gromacs/utility/smalloc.h"
54 #include "gromacs/utility/strconvert.h"
55 #include "gromacs/utility/txtdump.h"
56
57 static gmx::EnumerationArray<SimulationAtomGroupType, const char *>
58 c_simulationAtomGroupTypeShortNames
59     = { {
60             "T-Coupling",
61             "Energy Mon.",
62             "Acceleration",
63             "Freeze",
64             "User1",
65             "User2",
66             "VCM",
67             "Compressed X",
68             "Or. Res. Fit",
69             "QMMM"
70         } };
71
72 const char *shortName(SimulationAtomGroupType type)
73 {
74     return c_simulationAtomGroupTypeShortNames[type];
75 }
76
77 SimulationGroups::SimulationGroups()
78 {
79     for (auto group : keysOf(groups))
80     {
81         groups[group].nr     = 0;
82         groups[group].nm_ind = nullptr;
83     }
84
85 }
86
87 void init_top(t_topology *top)
88 {
89     top->name = nullptr;
90     init_idef(&top->idef);
91     init_atom(&(top->atoms));
92     init_atomtypes(&(top->atomtypes));
93     init_block(&top->cgs);
94     init_block(&top->mols);
95     init_blocka(&top->excls);
96     open_symtab(&top->symtab);
97 }
98
99
100 gmx_moltype_t::gmx_moltype_t() :
101     name(nullptr),
102     cgs(),
103     excls()
104 {
105     init_t_atoms(&atoms, 0, FALSE);
106 }
107
108 gmx_moltype_t::~gmx_moltype_t()
109 {
110     done_atom(&atoms);
111     done_block(&cgs);
112     done_blocka(&excls);
113 }
114
115 SimulationGroups::~SimulationGroups()
116 {
117     for (auto group : keysOf(groups))
118     {
119         if (nullptr != groups[group].nm_ind)
120         {
121             sfree(groups[group].nm_ind);
122             groups[group].nm_ind = nullptr;
123         }
124     }
125 }
126
127 gmx_mtop_t::gmx_mtop_t()
128 {
129     init_atomtypes(&atomtypes);
130     open_symtab(&symtab);
131 }
132
133 gmx_mtop_t::~gmx_mtop_t()
134 {
135     done_symtab(&symtab);
136
137     moltype.clear();
138     molblock.clear();
139     done_atomtypes(&atomtypes);
140 }
141
142 void done_top(t_topology *top)
143 {
144     done_idef(&top->idef);
145     done_atom(&(top->atoms));
146
147     /* For GB */
148     done_atomtypes(&(top->atomtypes));
149
150     done_symtab(&(top->symtab));
151     done_block(&(top->cgs));
152     done_block(&(top->mols));
153     done_blocka(&(top->excls));
154 }
155
156 void done_top_mtop(t_topology *top, gmx_mtop_t *mtop)
157 {
158     if (mtop != nullptr)
159     {
160         if (top != nullptr)
161         {
162             done_idef(&top->idef);
163             done_atom(&top->atoms);
164             done_block(&top->cgs);
165             done_blocka(&top->excls);
166             done_block(&top->mols);
167             done_symtab(&top->symtab);
168             open_symtab(&mtop->symtab);
169             done_atomtypes(&(top->atomtypes));
170         }
171
172         // Note that the rest of mtop will be freed by the destructor
173     }
174 }
175
176 gmx_localtop_t::gmx_localtop_t()
177 {
178     init_block_null(&cgs);
179     init_blocka_null(&excls);
180     init_idef(&idef);
181     init_atomtypes(&atomtypes);
182 }
183
184 gmx_localtop_t::~gmx_localtop_t()
185 {
186     if (!useInDomainDecomp_)
187     {
188         done_idef(&idef);
189         done_block(&cgs);
190         done_blocka(&excls);
191         done_atomtypes(&atomtypes);
192     }
193 }
194
195 bool gmx_mtop_has_masses(const gmx_mtop_t *mtop)
196 {
197     if (mtop == nullptr)
198     {
199         return false;
200     }
201     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveMass;
202 }
203
204 bool gmx_mtop_has_charges(const gmx_mtop_t *mtop)
205 {
206     if (mtop == nullptr)
207     {
208         return false;
209     }
210     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveCharge;
211 }
212
213 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t &mtop)
214 {
215     for (const gmx_moltype_t &moltype : mtop.moltype)
216     {
217         const t_atoms &atoms = moltype.atoms;
218         if (atoms.haveBState)
219         {
220             for (int a = 0; a < atoms.nr; a++)
221             {
222                 if (atoms.atom[a].q != atoms.atom[a].qB)
223                 {
224                     return true;
225                 }
226             }
227         }
228     }
229     return false;
230 }
231
232 bool gmx_mtop_has_atomtypes(const gmx_mtop_t *mtop)
233 {
234     if (mtop == nullptr)
235     {
236         return false;
237     }
238     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveType;
239 }
240
241 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t *mtop)
242 {
243     if (mtop == nullptr)
244     {
245         return false;
246     }
247     return mtop->moltype.empty() || mtop->moltype[0].atoms.havePdbInfo;
248 }
249
250 static void pr_grps(FILE                       *fp,
251                     const char                 *title,
252                     gmx::ArrayRef<const t_grps> grps,
253                     char                     ***grpname)
254 {
255     int index = 0;
256     for (const auto &group : grps)
257     {
258         fprintf(fp, "%s[%-12s] nr=%d, name=[", title, c_simulationAtomGroupTypeShortNames[index], group.nr);
259         for (int j = 0; (j < group.nr); j++)
260         {
261             fprintf(fp, " %s", *(grpname[group.nm_ind[j]]));
262         }
263         fprintf(fp, "]\n");
264         index++;
265     }
266 }
267
268 static void pr_groups(FILE *fp, int indent,
269                       const SimulationGroups &groups,
270                       gmx_bool bShowNumbers)
271 {
272     pr_grps(fp,
273             "grp",
274             groups.groups,
275             const_cast<char ***>(groups.groupNames.data()));
276     pr_strings(fp, indent, "grpname", const_cast<char ***>(groups.groupNames.data()), groups.groupNames.size(), bShowNumbers);
277
278     pr_indent(fp, indent);
279     fprintf(fp, "groups          ");
280     for (const auto &group : c_simulationAtomGroupTypeShortNames)
281     {
282         printf(" %5.5s", group);
283     }
284     printf("\n");
285
286     pr_indent(fp, indent);
287     fprintf(fp, "allocated       ");
288     int nat_max = 0;
289     for (auto group : keysOf(groups.groups))
290     {
291         printf(" %5d", groups.numberOfGroupNumbers(group));
292         nat_max = std::max(nat_max, groups.numberOfGroupNumbers(group));
293     }
294     printf("\n");
295
296     if (nat_max == 0)
297     {
298         pr_indent(fp, indent);
299         fprintf(fp, "groupnr[%5s] =", "*");
300         for (auto gmx_unused group : keysOf(groups.groups))
301         {
302             fprintf(fp, "  %3d ", 0);
303         }
304         fprintf(fp, "\n");
305     }
306     else
307     {
308         for (int i = 0; i < nat_max; i++)
309         {
310             pr_indent(fp, indent);
311             fprintf(fp, "groupnr[%5d] =", i);
312             for (auto group : keysOf(groups.groups))
313             {
314                 fprintf(fp, "  %3d ",
315                         !groups.groupNumbers[group].empty() ?
316                         groups.groupNumbers[group][i] : 0);
317             }
318             fprintf(fp, "\n");
319         }
320     }
321 }
322
323 static void pr_moltype(FILE *fp, int indent, const char *title,
324                        const gmx_moltype_t *molt, int n,
325                        const gmx_ffparams_t *ffparams,
326                        gmx_bool bShowNumbers, gmx_bool bShowParameters)
327 {
328     int j;
329
330     indent = pr_title_n(fp, indent, title, n);
331     pr_indent(fp, indent);
332     fprintf(fp, "name=\"%s\"\n", *(molt->name));
333     pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
334     pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
335     pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
336     for (j = 0; (j < F_NRE); j++)
337     {
338         pr_ilist(fp, indent, interaction_function[j].longname,
339                  ffparams->functype.data(), molt->ilist[j],
340                  bShowNumbers, bShowParameters, ffparams->iparams.data());
341     }
342 }
343
344 static void pr_molblock(FILE *fp, int indent, const char *title,
345                         const gmx_molblock_t *molb, int n,
346                         const std::vector<gmx_moltype_t> &molt)
347 {
348     indent = pr_title_n(fp, indent, title, n);
349     pr_indent(fp, indent);
350     fprintf(fp, "%-20s = %d \"%s\"\n",
351             "moltype", molb->type, *(molt[molb->type].name));
352     pr_int(fp, indent, "#molecules", molb->nmol);
353     pr_int(fp, indent, "#posres_xA", molb->posres_xA.size());
354     if (!molb->posres_xA.empty())
355     {
356         pr_rvecs(fp, indent, "posres_xA", as_rvec_array(molb->posres_xA.data()), molb->posres_xA.size());
357     }
358     pr_int(fp, indent, "#posres_xB", molb->posres_xB.size());
359     if (!molb->posres_xB.empty())
360     {
361         pr_rvecs(fp, indent, "posres_xB", as_rvec_array(molb->posres_xB.data()), molb->posres_xB.size());
362     }
363 }
364
365 void pr_mtop(FILE *fp, int indent, const char *title, const gmx_mtop_t *mtop,
366              gmx_bool bShowNumbers, gmx_bool bShowParameters)
367 {
368     if (available(fp, mtop, indent, title))
369     {
370         indent = pr_title(fp, indent, title);
371         pr_indent(fp, indent);
372         fprintf(fp, "name=\"%s\"\n", *(mtop->name));
373         pr_int(fp, indent, "#atoms", mtop->natoms);
374         pr_int(fp, indent, "#molblock", mtop->molblock.size());
375         for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
376         {
377             pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
378         }
379         pr_str(fp, indent, "bIntermolecularInteractions",
380                gmx::boolToString(mtop->bIntermolecularInteractions));
381         if (mtop->bIntermolecularInteractions)
382         {
383             for (int j = 0; j < F_NRE; j++)
384             {
385                 pr_ilist(fp, indent, interaction_function[j].longname,
386                          mtop->ffparams.functype.data(),
387                          (*mtop->intermolecular_ilist)[j],
388                          bShowNumbers, bShowParameters, mtop->ffparams.iparams.data());
389             }
390         }
391         pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
392         pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
393         for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
394         {
395             pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
396                        &mtop->ffparams, bShowNumbers, bShowParameters);
397         }
398         pr_groups(fp, indent, mtop->groups, bShowNumbers);
399     }
400 }
401
402 void pr_top(FILE *fp, int indent, const char *title, const t_topology *top,
403             gmx_bool bShowNumbers, gmx_bool bShowParameters)
404 {
405     if (available(fp, top, indent, title))
406     {
407         indent = pr_title(fp, indent, title);
408         pr_indent(fp, indent);
409         fprintf(fp, "name=\"%s\"\n", *(top->name));
410         pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
411         pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
412         pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
413         pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
414         pr_str(fp, indent, "bIntermolecularInteractions",
415                gmx::boolToString(top->bIntermolecularInteractions));
416         pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
417         pr_idef(fp, indent, "idef", &top->idef, bShowNumbers, bShowParameters);
418     }
419 }
420
421 static void cmp_iparm(FILE *fp, const char *s, t_functype ft,
422                       const t_iparams &ip1, const t_iparams &ip2, real relativeTolerance, real absoluteTolerance)
423 {
424     int      i;
425     gmx_bool bDiff;
426
427     bDiff = FALSE;
428     for (i = 0; i < MAXFORCEPARAM && !bDiff; i++)
429     {
430         bDiff = !equal_real(ip1.generic.buf[i], ip2.generic.buf[i], relativeTolerance, absoluteTolerance);
431     }
432     if (bDiff)
433     {
434         fprintf(fp, "%s1: ", s);
435         pr_iparams(fp, ft, &ip1);
436         fprintf(fp, "%s2: ", s);
437         pr_iparams(fp, ft, &ip2);
438     }
439 }
440
441 static void cmp_iparm_AB(FILE *fp, const char *s, t_functype ft,
442                          const t_iparams &ip1, real relativeTolerance, real absoluteTolerance)
443 {
444     int      nrfpA, nrfpB, p0, i;
445     gmx_bool bDiff;
446
447     /* Normally the first parameter is perturbable */
448     p0    = 0;
449     nrfpA = interaction_function[ft].nrfpA;
450     nrfpB = interaction_function[ft].nrfpB;
451     if (ft == F_PDIHS)
452     {
453         nrfpB = 2;
454     }
455     else if (interaction_function[ft].flags & IF_TABULATED)
456     {
457         /* For tabulated interactions only the second parameter is perturbable */
458         p0    = 1;
459         nrfpB = 1;
460     }
461     bDiff = FALSE;
462     for (i = 0; i < nrfpB && !bDiff; i++)
463     {
464         bDiff = !equal_real(ip1.generic.buf[p0+i], ip1.generic.buf[nrfpA+i], relativeTolerance, absoluteTolerance);
465     }
466     if (bDiff)
467     {
468         fprintf(fp, "%s: ", s);
469         pr_iparams(fp, ft, &ip1);
470     }
471 }
472
473 static void cmp_cmap(FILE *fp, const gmx_cmap_t *cmap1, const gmx_cmap_t *cmap2, real relativeTolerance, real absoluteTolerance)
474 {
475     int cmap1_ngrid = (cmap1 ? cmap1->cmapdata.size() : 0);
476     int cmap2_ngrid = (cmap2 ? cmap2->cmapdata.size() : 0);
477
478     cmp_int(fp, "cmap ngrid", -1, cmap1_ngrid, cmap2_ngrid);
479
480     if (cmap1 == nullptr || cmap2 == nullptr)
481     {
482         return;
483     }
484
485     cmp_int(fp, "cmap grid_spacing", -1, cmap1->grid_spacing, cmap2->grid_spacing);
486     if (cmap1->cmapdata.size() == cmap2->cmapdata.size() &&
487         cmap1->grid_spacing == cmap2->grid_spacing)
488     {
489         for (size_t g = 0; g < cmap1->cmapdata.size(); g++)
490         {
491             int i;
492
493             fprintf(fp, "comparing cmap %zu\n", g);
494
495             for (i = 0; i < 4*cmap1->grid_spacing*cmap1->grid_spacing; i++)
496             {
497                 cmp_real(fp, "", i, cmap1->cmapdata[g].cmap[i], cmap2->cmapdata[g].cmap[i], relativeTolerance, absoluteTolerance);
498             }
499         }
500     }
501 }
502
503 static void cmp_block(FILE *fp, const t_block *b1, const t_block *b2, const char *s)
504 {
505     char buf[32];
506
507     fprintf(fp, "comparing block %s\n", s);
508     sprintf(buf, "%s.nr", s);
509     cmp_int(fp, buf, -1, b1->nr, b2->nr);
510 }
511
512 static void cmp_blocka(FILE *fp, const t_blocka *b1, const t_blocka *b2, const char *s)
513 {
514     char buf[32];
515
516     fprintf(fp, "comparing blocka %s\n", s);
517     sprintf(buf, "%s.nr", s);
518     cmp_int(fp, buf, -1, b1->nr, b2->nr);
519     sprintf(buf, "%s.nra", s);
520     cmp_int(fp, buf, -1, b1->nra, b2->nra);
521 }
522
523 static void compareFfparams(FILE *fp, const gmx_ffparams_t &ff1, const gmx_ffparams_t &ff2, real relativeTolerance, real absoluteTolerance)
524 {
525     fprintf(fp, "comparing force field parameters\n");
526     cmp_int(fp, "numTypes", -1, ff1.numTypes(), ff2.numTypes());
527     cmp_int(fp, "atnr", -1, ff1.atnr, ff1.atnr);
528     cmp_double(fp, "reppow", -1, ff1.reppow, ff2.reppow, relativeTolerance, absoluteTolerance);
529     cmp_real(fp, "fudgeQQ", -1, ff1.fudgeQQ, ff2.fudgeQQ, relativeTolerance, absoluteTolerance);
530     cmp_cmap(fp, &ff1.cmap_grid, &ff2.cmap_grid, relativeTolerance, absoluteTolerance);
531     for (int i = 0; i < std::min(ff1.numTypes(), ff2.numTypes()); i++)
532     {
533         std::string buf = gmx::formatString("ffparams->functype[%d]", i);
534         cmp_int(fp, buf.c_str(), i, ff1.functype[i], ff2.functype[i]);
535         buf = gmx::formatString("ffparams->iparams[%d]", i);
536         cmp_iparm(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], ff2.iparams[i], relativeTolerance, absoluteTolerance);
537     }
538
539 }
540
541 static void compareFfparamAB(FILE *fp, const gmx_ffparams_t &ff1, real relativeTolerance,   real absoluteTolerance)
542 {
543     fprintf(fp, "comparing free energy parameters\n");
544     for (int i = 0; i < ff1.numTypes(); i++)
545     {
546         std::string buf = gmx::formatString("ffparams->iparams[%d]", i);
547         cmp_iparm_AB(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], relativeTolerance, absoluteTolerance);
548     }
549 }
550 static void compareInteractionLists(FILE *fp, const InteractionLists *il1, const InteractionLists *il2)
551 {
552     fprintf(fp, "comparing InteractionLists\n");
553     if ((il1 || il2) && (!il1 || !il2))
554     {
555         fprintf(fp, "InteractionLists are present in topology %d but not in the other\n", il1 ? 1 : 2);
556     }
557     if (il1 && il2)
558     {
559         for (int i = 0; i < F_NRE; i++)
560         {
561             cmp_int(fp, "InteractionList size", i, il1->at(i).size(), il2->at(i).size());
562             int nr = std::min(il1->at(i).size(), il2->at(i).size());
563             for (int j = 0; j < nr; j++)
564             {
565                 cmp_int(fp, "InteractionList entry", j, il1->at(i).iatoms.at(j), il2->at(i).iatoms.at(j));
566             }
567         }
568     }
569 }
570
571 static void compareMoltypes(FILE *fp, gmx::ArrayRef<const gmx_moltype_t> mt1, gmx::ArrayRef<const gmx_moltype_t> mt2, real relativeTolerance, real absoluteTolerance)
572 {
573     fprintf(fp, "comparing molecule types\n");
574     cmp_int(fp, "moltype size", -1, mt1.size(), mt2.size());
575     for (int i = 0; i < std::min(mt1.ssize(), mt2.ssize()); i++)
576     {
577         cmp_str(fp, "Name", i, *mt1[i].name, *mt2[i].name);
578         compareAtoms(fp, &mt1[i].atoms, &mt2[i].atoms, relativeTolerance, absoluteTolerance);
579         compareInteractionLists(fp, &mt1[i].ilist, &mt2[i].ilist);
580         std::string buf = gmx::formatString("cgs[%d]", i);
581         cmp_block(fp, &mt1[i].cgs, &mt2[i].cgs, buf.c_str());
582         buf = gmx::formatString("excls[%d]", i);
583         cmp_blocka(fp, &mt1[i].excls, &mt2[i].excls, buf.c_str());
584     }
585 }
586
587 static void compareMoletypeAB(FILE *fp, gmx::ArrayRef<const gmx_moltype_t> mt1, real relativeTolerance, real absoluteTolerance)
588 {
589     fprintf(fp, "comparing free energy molecule types\n");
590     for (int i = 0; i < mt1.ssize(); i++)
591     {
592         compareAtoms(fp, &mt1[i].atoms, nullptr, relativeTolerance, absoluteTolerance);
593     }
594 }
595 static void compareMolblocks(FILE *fp, gmx::ArrayRef<const gmx_molblock_t> mb1, gmx::ArrayRef<const gmx_molblock_t> mb2)
596 {
597     fprintf(fp, "comparing molecule blocks\n");
598     cmp_int(fp, "molblock size", -1, mb1.size(), mb2.size());
599     int nr = std::min(mb1.size(), mb2.size());
600     for (int i = 0; i < nr; i++)
601     {
602         cmp_int(fp, "type", i, mb1[i].type, mb2[i].type);
603         cmp_int(fp, "nmol", i, mb1[i].nmol, mb2[i].nmol);
604         // Only checking size of restraint vectors for now
605         cmp_int(fp, "posres_xA size", i, mb1[i].posres_xA.size(), mb2[i].posres_xA.size());
606         cmp_int(fp, "posres_xB size", i, mb1[i].posres_xB.size(), mb2[i].posres_xB.size());
607     }
608
609 }
610
611 static void compareAtomtypes(FILE *fp, const t_atomtypes &at1, const t_atomtypes &at2)
612 {
613     fprintf(fp, "comparing atomtypes\n");
614     cmp_int(fp, "nr", -1, at1.nr, at2.nr);
615     int nr = std::min(at1.nr, at2.nr);
616     for (int i = 0; i < nr; i++)
617     {
618         cmp_int(fp, "atomtype", i, at1.atomnumber[i], at2.atomnumber[i]);
619     }
620 }
621
622 static void compareIntermolecularExclusions(FILE *fp, gmx::ArrayRef<const int> ime1, gmx::ArrayRef<const int> ime2)
623 {
624     fprintf(fp, "comparing intermolecular exclusions\n");
625     cmp_int(fp, "exclusion number", -1, ime1.size(), ime2.size());
626     int nr = std::min(ime1.size(), ime2.size());
627     for (int i = 0; i < nr; i++)
628     {
629         cmp_int(fp, "exclusion", i, ime1[i], ime2[i]);
630     }
631 }
632
633 static void compareBlockIndices(FILE *fp, gmx::ArrayRef<const MoleculeBlockIndices> mbi1, gmx::ArrayRef<const MoleculeBlockIndices> mbi2)
634 {
635     fprintf(fp, "comparing moleculeBlockIndices\n");
636     cmp_int(fp, "size", -1, mbi1.size(), mbi2.size());
637     int nr = std::min(mbi1.size(), mbi2.size());
638     for (int i = 0; i < nr; i++)
639     {
640         cmp_int(fp, "numAtomsPerMolecule", i, mbi1[i].numAtomsPerMolecule, mbi2[i].numAtomsPerMolecule);
641         cmp_int(fp, "globalAtomStart", i, mbi1[i].globalAtomStart, mbi2[i].globalAtomStart);
642         cmp_int(fp, "globalAtomEnd", i, mbi1[i].globalAtomEnd, mbi2[i].globalAtomEnd);
643         cmp_int(fp, "globalResidueStart", i, mbi1[i].globalResidueStart, mbi2[i].globalResidueStart);
644         cmp_int(fp, "moleculeIndexStart", i, mbi1[i].moleculeIndexStart, mbi2[i].moleculeIndexStart);
645     }
646 }
647
648 void compareMtop(FILE *fp, const gmx_mtop_t &mtop1, const gmx_mtop_t &mtop2, real relativeTolerance, real absoluteTolerance)
649 {
650     fprintf(fp, "comparing mtop topology\n");
651     cmp_str(fp, "Name", -1, *mtop1.name, *mtop2.name);
652     cmp_int(fp, "natoms", -1, mtop1.natoms, mtop2.natoms);
653     cmp_int(fp, "maxres_renum", -1, mtop1.maxres_renum, mtop2.maxres_renum);
654     cmp_int(fp, "maxresnr", -1, mtop1.maxresnr, mtop2.maxresnr);
655     cmp_bool(fp, "bIntermolecularInteractions", -1, mtop1.bIntermolecularInteractions, mtop2.bIntermolecularInteractions);
656     cmp_bool(fp, "haveMoleculeIndices", -1, mtop1.haveMoleculeIndices, mtop2.haveMoleculeIndices);
657
658     compareFfparams(fp, mtop1.ffparams, mtop2.ffparams, relativeTolerance, absoluteTolerance);
659     compareMoltypes(fp, mtop1.moltype, mtop2.moltype, relativeTolerance, absoluteTolerance);
660     compareMolblocks(fp, mtop1.molblock, mtop2.molblock);
661     compareInteractionLists(fp, mtop1.intermolecular_ilist.get(), mtop2.intermolecular_ilist.get());
662     compareAtomtypes(fp, mtop1.atomtypes, mtop2.atomtypes);
663     compareAtomGroups(fp, mtop1.groups, mtop2.groups, mtop1.natoms, mtop2.natoms);
664     compareIntermolecularExclusions(fp, mtop1.intermolecularExclusionGroup, mtop2.intermolecularExclusionGroup);
665     compareBlockIndices(fp, mtop1.moleculeBlockIndices, mtop2.moleculeBlockIndices);
666 }
667
668 void compareMtopAB(FILE *fp, const gmx_mtop_t &mtop1, real relativeTolerance, real absoluteTolerance)
669 {
670     fprintf(fp, "comparing topAB\n");
671     compareFfparamAB(fp, mtop1.ffparams, relativeTolerance, absoluteTolerance);
672     compareMoletypeAB(fp, mtop1.moltype, relativeTolerance, absoluteTolerance);
673 }
674
675 void compareAtomGroups(FILE *fp, const SimulationGroups &g0, const SimulationGroups &g1,
676                        int natoms0, int natoms1)
677 {
678     fprintf(fp, "comparing groups\n");
679
680     for (auto group : keysOf(g0.groups))
681     {
682         std::string buf = gmx::formatString("grps[%d].nr", static_cast<int>(group));
683         cmp_int(fp, buf.c_str(), -1, g0.groups[group].nr, g1.groups[group].nr);
684         if (g0.groups[group].nr == g1.groups[group].nr)
685         {
686             for (int j = 0; j < g0.groups[group].nr; j++)
687             {
688                 buf = gmx::formatString("grps[%d].name[%d]", static_cast<int>(group), j);
689                 cmp_str(fp, buf.c_str(), -1,
690                         *g0.groupNames[g0.groups[group].nm_ind[j]],
691                         *g1.groupNames[g1.groups[group].nm_ind[j]]);
692             }
693         }
694         cmp_int(fp, "ngrpnr", static_cast<int>(group), g0.numberOfGroupNumbers(group), g1.numberOfGroupNumbers(group));
695         if (g0.numberOfGroupNumbers(group) == g1.numberOfGroupNumbers(group) && natoms0 == natoms1 &&
696             (!g0.groupNumbers[group].empty() || !g1.groupNumbers[group].empty()))
697         {
698             for (int j = 0; j < natoms0; j++)
699             {
700                 cmp_int(fp, c_simulationAtomGroupTypeShortNames[group], j, getGroupType(g0, group, j), getGroupType(g1, group, j));
701             }
702         }
703     }
704     /* We have compared the names in the groups lists,
705      * so we can skip the grpname list comparison.
706      */
707 }
708
709 int getGroupType(const SimulationGroups &group, SimulationAtomGroupType type, int atom)
710 {
711     return (group.groupNumbers[type].empty() ? 0 : group.groupNumbers[type][atom]);
712 }
713
714 void copy_moltype(const gmx_moltype_t *src, gmx_moltype_t *dst)
715 {
716     dst->name = src->name;
717     copy_blocka(&src->excls, &dst->excls);
718     copy_block(&src->cgs, &dst->cgs);
719     t_atoms *atomsCopy = copy_t_atoms(&src->atoms);
720     dst->atoms = *atomsCopy;
721     sfree(atomsCopy);
722
723     for (int i = 0; i < F_NRE; ++i)
724     {
725         dst->ilist[i] = src->ilist[i];
726     }
727 }