807eb7de2312e3ce3f6387700fe3794878606242
[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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "topology.h"
41
42 #include <cstdio>
43
44 #include <algorithm>
45
46 #include "gromacs/math/vecdump.h"
47 #include "gromacs/topology/atoms.h"
48 #include "gromacs/topology/idef.h"
49 #include "gromacs/topology/ifunc.h"
50 #include "gromacs/topology/symtab.h"
51 #include "gromacs/utility/arrayref.h"
52 #include "gromacs/utility/compare.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/smalloc.h"
55 #include "gromacs/utility/strconvert.h"
56 #include "gromacs/utility/txtdump.h"
57
58 static gmx::EnumerationArray<SimulationAtomGroupType, const char*> c_simulationAtomGroupTypeShortNames = {
59     { "T-Coupling",
60       "Energy Mon.",
61       "Acc. not used",
62       "Freeze",
63       "User1",
64       "User2",
65       "VCM",
66       "Compressed X",
67       "Or. Res. Fit",
68       "QMMM" }
69 };
70
71 const char* shortName(SimulationAtomGroupType type)
72 {
73     return c_simulationAtomGroupTypeShortNames[type];
74 }
75
76 void init_top(t_topology* top)
77 {
78     top->name = nullptr;
79     init_idef(&top->idef);
80     init_atom(&(top->atoms));
81     init_atomtypes(&(top->atomtypes));
82     init_block(&top->mols);
83     open_symtab(&top->symtab);
84 }
85
86
87 gmx_moltype_t::gmx_moltype_t() : name(nullptr)
88 {
89     init_t_atoms(&atoms, 0, FALSE);
90 }
91
92 gmx_moltype_t::~gmx_moltype_t()
93 {
94     done_atom(&atoms);
95 }
96
97 static int gmx_mtop_maxresnr(const gmx::ArrayRef<const gmx_moltype_t> moltypes, int maxres_renum)
98 {
99     int maxresnr = 0;
100
101     for (const gmx_moltype_t& moltype : moltypes)
102     {
103         const t_atoms& atoms = moltype.atoms;
104         if (atoms.nres > maxres_renum)
105         {
106             for (int r = 0; r < atoms.nres; r++)
107             {
108                 if (atoms.resinfo[r].nr > maxresnr)
109                 {
110                     maxresnr = atoms.resinfo[r].nr;
111                 }
112             }
113         }
114     }
115
116     return maxresnr;
117 }
118
119 gmx_mtop_t::gmx_mtop_t()
120 {
121     init_atomtypes(&atomtypes);
122     open_symtab(&symtab);
123 }
124
125 gmx_mtop_t::~gmx_mtop_t()
126 {
127     done_symtab(&symtab);
128
129     moltype.clear();
130     molblock.clear();
131     done_atomtypes(&atomtypes);
132 }
133
134 void gmx_mtop_t::finalize()
135 {
136     if (molblock.size() == 1 && molblock[0].nmol == 1)
137     {
138         /* We have a single molecule only, no renumbering needed.
139          * This case also covers an mtop converted from pdb/gro/... input,
140          * so we retain the original residue numbering.
141          */
142         maxResiduesPerMoleculeToTriggerRenumber_ = 0;
143     }
144     else
145     {
146         /* We only renumber single residue molecules. Their intra-molecular
147          * residue numbering is anyhow irrelevant.
148          */
149         maxResiduesPerMoleculeToTriggerRenumber_ = 1;
150     }
151
152     const char* env = getenv("GMX_MAXRESRENUM");
153     if (env != nullptr)
154     {
155         sscanf(env, "%d", &maxResiduesPerMoleculeToTriggerRenumber_);
156     }
157     if (maxResiduesPerMoleculeToTriggerRenumber_ == -1)
158     {
159         /* -1 signals renumber residues in all molecules */
160         maxResiduesPerMoleculeToTriggerRenumber_ = std::numeric_limits<int>::max();
161     }
162
163     maxResNumberNotRenumbered_ = gmx_mtop_maxresnr(moltype, maxResiduesPerMoleculeToTriggerRenumber_);
164
165     buildMolblockIndices();
166 }
167
168 void gmx_mtop_t::buildMolblockIndices()
169 {
170     moleculeBlockIndices.resize(molblock.size());
171
172     int atomIndex          = 0;
173     int residueIndex       = 0;
174     int residueNumberStart = maxResNumberNotRenumbered_ + 1;
175     int moleculeIndexStart = 0;
176     for (size_t mb = 0; mb < molblock.size(); mb++)
177     {
178         const gmx_molblock_t& molb         = molblock[mb];
179         MoleculeBlockIndices& indices      = moleculeBlockIndices[mb];
180         const int             numResPerMol = moltype[molb.type].atoms.nres;
181
182         indices.numAtomsPerMolecule = moltype[molb.type].atoms.nr;
183         indices.globalAtomStart     = atomIndex;
184         indices.globalResidueStart  = residueIndex;
185         atomIndex += molb.nmol * indices.numAtomsPerMolecule;
186         residueIndex += molb.nmol * numResPerMol;
187         indices.globalAtomEnd      = atomIndex;
188         indices.residueNumberStart = residueNumberStart;
189         if (numResPerMol <= maxResiduesPerMoleculeToTriggerRenumber_)
190         {
191             residueNumberStart += molb.nmol * numResPerMol;
192         }
193         indices.moleculeIndexStart = moleculeIndexStart;
194         moleculeIndexStart += molb.nmol;
195     }
196 }
197
198 void done_top(t_topology* top)
199 {
200     done_idef(&top->idef);
201     done_atom(&(top->atoms));
202
203     /* For GB */
204     done_atomtypes(&(top->atomtypes));
205
206     done_symtab(&(top->symtab));
207     done_block(&(top->mols));
208 }
209
210 void done_top_mtop(t_topology* top, gmx_mtop_t* mtop)
211 {
212     if (mtop != nullptr)
213     {
214         if (top != nullptr)
215         {
216             done_idef(&top->idef);
217             done_atom(&top->atoms);
218             done_block(&top->mols);
219             done_symtab(&top->symtab);
220             open_symtab(&mtop->symtab);
221             done_atomtypes(&(top->atomtypes));
222         }
223
224         // Note that the rest of mtop will be freed by the destructor
225     }
226 }
227
228 gmx_localtop_t::gmx_localtop_t(const gmx_ffparams_t& ffparams) : idef(ffparams) {}
229
230 bool gmx_mtop_has_masses(const gmx_mtop_t* mtop)
231 {
232     if (mtop == nullptr)
233     {
234         return false;
235     }
236     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveMass;
237 }
238
239 bool gmx_mtop_has_charges(const gmx_mtop_t* mtop)
240 {
241     if (mtop == nullptr)
242     {
243         return false;
244     }
245     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveCharge;
246 }
247
248 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t& mtop)
249 {
250     for (const gmx_moltype_t& moltype : mtop.moltype)
251     {
252         const t_atoms& atoms = moltype.atoms;
253         if (atoms.haveBState)
254         {
255             for (int a = 0; a < atoms.nr; a++)
256             {
257                 if (atoms.atom[a].q != atoms.atom[a].qB)
258                 {
259                     return true;
260                 }
261             }
262         }
263     }
264     return false;
265 }
266
267 bool gmx_mtop_has_atomtypes(const gmx_mtop_t* mtop)
268 {
269     if (mtop == nullptr)
270     {
271         return false;
272     }
273     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveType;
274 }
275
276 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t* mtop)
277 {
278     if (mtop == nullptr)
279     {
280         return false;
281     }
282     return mtop->moltype.empty() || mtop->moltype[0].atoms.havePdbInfo;
283 }
284
285 static void pr_grps(FILE* fp, const char* title, gmx::ArrayRef<const AtomGroupIndices> grps, char*** grpname)
286 {
287     int index = 0;
288     for (const auto& group : grps)
289     {
290         fprintf(fp, "%s[%-12s] nr=%zu, name=[", title, c_simulationAtomGroupTypeShortNames[index], group.size());
291         for (const auto& entry : group)
292         {
293             fprintf(fp, " %s", *(grpname[entry]));
294         }
295         fprintf(fp, "]\n");
296         index++;
297     }
298 }
299
300 static void pr_groups(FILE* fp, int indent, const SimulationGroups& groups, gmx_bool bShowNumbers)
301 {
302     pr_grps(fp, "grp", groups.groups, const_cast<char***>(groups.groupNames.data()));
303     pr_strings(fp,
304                indent,
305                "grpname",
306                const_cast<char***>(groups.groupNames.data()),
307                groups.groupNames.size(),
308                bShowNumbers);
309
310     pr_indent(fp, indent);
311     fprintf(fp, "groups          ");
312     for (const auto& group : c_simulationAtomGroupTypeShortNames)
313     {
314         printf(" %5.5s", group);
315     }
316     printf("\n");
317
318     pr_indent(fp, indent);
319     fprintf(fp, "allocated       ");
320     int nat_max = 0;
321     for (auto group : keysOf(groups.groups))
322     {
323         printf(" %5d", groups.numberOfGroupNumbers(group));
324         nat_max = std::max(nat_max, groups.numberOfGroupNumbers(group));
325     }
326     printf("\n");
327
328     if (nat_max == 0)
329     {
330         pr_indent(fp, indent);
331         fprintf(fp, "groupnr[%5s] =", "*");
332         for (auto gmx_unused group : keysOf(groups.groups))
333         {
334             fprintf(fp, "  %3d ", 0);
335         }
336         fprintf(fp, "\n");
337     }
338     else
339     {
340         for (int i = 0; i < nat_max; i++)
341         {
342             pr_indent(fp, indent);
343             fprintf(fp, "groupnr[%5d] =", i);
344             for (auto group : keysOf(groups.groups))
345             {
346                 fprintf(fp,
347                         "  %3d ",
348                         !groups.groupNumbers[group].empty() ? groups.groupNumbers[group][i] : 0);
349             }
350             fprintf(fp, "\n");
351         }
352     }
353 }
354
355 static void pr_moltype(FILE*                 fp,
356                        int                   indent,
357                        const char*           title,
358                        const gmx_moltype_t*  molt,
359                        int                   n,
360                        const gmx_ffparams_t* ffparams,
361                        gmx_bool              bShowNumbers,
362                        gmx_bool              bShowParameters)
363 {
364     int j;
365
366     indent = pr_title_n(fp, indent, title, n);
367     pr_indent(fp, indent);
368     fprintf(fp, "name=\"%s\"\n", *(molt->name));
369     pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
370     pr_listoflists(fp, indent, "excls", &molt->excls, bShowNumbers);
371     for (j = 0; (j < F_NRE); j++)
372     {
373         pr_ilist(fp,
374                  indent,
375                  interaction_function[j].longname,
376                  ffparams->functype.data(),
377                  molt->ilist[j],
378                  bShowNumbers,
379                  bShowParameters,
380                  ffparams->iparams.data());
381     }
382 }
383
384 static void pr_molblock(FILE*                             fp,
385                         int                               indent,
386                         const char*                       title,
387                         const gmx_molblock_t*             molb,
388                         int                               n,
389                         const std::vector<gmx_moltype_t>& molt)
390 {
391     indent = pr_title_n(fp, indent, title, n);
392     pr_indent(fp, indent);
393     fprintf(fp, "%-20s = %d \"%s\"\n", "moltype", molb->type, *(molt[molb->type].name));
394     pr_int(fp, indent, "#molecules", molb->nmol);
395     pr_int(fp, indent, "#posres_xA", molb->posres_xA.size());
396     if (!molb->posres_xA.empty())
397     {
398         pr_rvecs(fp, indent, "posres_xA", as_rvec_array(molb->posres_xA.data()), molb->posres_xA.size());
399     }
400     pr_int(fp, indent, "#posres_xB", molb->posres_xB.size());
401     if (!molb->posres_xB.empty())
402     {
403         pr_rvecs(fp, indent, "posres_xB", as_rvec_array(molb->posres_xB.data()), molb->posres_xB.size());
404     }
405 }
406
407 void pr_mtop(FILE* fp, int indent, const char* title, const gmx_mtop_t* mtop, gmx_bool bShowNumbers, gmx_bool bShowParameters)
408 {
409     if (available(fp, mtop, indent, title))
410     {
411         indent = pr_title(fp, indent, title);
412         pr_indent(fp, indent);
413         fprintf(fp, "name=\"%s\"\n", *(mtop->name));
414         pr_int(fp, indent, "#atoms", mtop->natoms);
415         pr_int(fp, indent, "#molblock", mtop->molblock.size());
416         for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
417         {
418             pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
419         }
420         pr_str(fp, indent, "bIntermolecularInteractions", gmx::boolToString(mtop->bIntermolecularInteractions));
421         if (mtop->bIntermolecularInteractions)
422         {
423             for (int j = 0; j < F_NRE; j++)
424             {
425                 pr_ilist(fp,
426                          indent,
427                          interaction_function[j].longname,
428                          mtop->ffparams.functype.data(),
429                          (*mtop->intermolecular_ilist)[j],
430                          bShowNumbers,
431                          bShowParameters,
432                          mtop->ffparams.iparams.data());
433             }
434         }
435         pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
436         pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
437         for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
438         {
439             pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt, &mtop->ffparams, bShowNumbers, bShowParameters);
440         }
441         pr_groups(fp, indent, mtop->groups, bShowNumbers);
442     }
443 }
444
445 void pr_top(FILE* fp, int indent, const char* title, const t_topology* top, gmx_bool bShowNumbers, gmx_bool bShowParameters)
446 {
447     if (available(fp, top, indent, title))
448     {
449         indent = pr_title(fp, indent, title);
450         pr_indent(fp, indent);
451         fprintf(fp, "name=\"%s\"\n", *(top->name));
452         pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
453         pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
454         pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
455         pr_str(fp, indent, "bIntermolecularInteractions", gmx::boolToString(top->bIntermolecularInteractions));
456         pr_idef(fp, indent, "idef", &top->idef, bShowNumbers, bShowParameters);
457     }
458 }
459
460 static void cmp_iparm(FILE*            fp,
461                       const char*      s,
462                       t_functype       ft,
463                       const t_iparams& ip1,
464                       const t_iparams& ip2,
465                       real             relativeTolerance,
466                       real             absoluteTolerance)
467 {
468     int      i;
469     gmx_bool bDiff;
470
471     bDiff = FALSE;
472     for (i = 0; i < MAXFORCEPARAM && !bDiff; i++)
473     {
474         bDiff = !equal_real(ip1.generic.buf[i], ip2.generic.buf[i], relativeTolerance, absoluteTolerance);
475     }
476     if (bDiff)
477     {
478         fprintf(fp, "%s1: ", s);
479         pr_iparams(fp, ft, ip1);
480         fprintf(fp, "%s2: ", s);
481         pr_iparams(fp, ft, ip2);
482     }
483 }
484
485 static void cmp_iparm_AB(FILE* fp, const char* s, t_functype ft, const t_iparams& ip1, real relativeTolerance, real absoluteTolerance)
486 {
487     int      nrfpA, nrfpB, p0, i;
488     gmx_bool bDiff;
489
490     /* Normally the first parameter is perturbable */
491     p0    = 0;
492     nrfpA = interaction_function[ft].nrfpA;
493     nrfpB = interaction_function[ft].nrfpB;
494     if (ft == F_PDIHS)
495     {
496         nrfpB = 2;
497     }
498     else if (interaction_function[ft].flags & IF_TABULATED)
499     {
500         /* For tabulated interactions only the second parameter is perturbable */
501         p0    = 1;
502         nrfpB = 1;
503     }
504     bDiff = FALSE;
505     for (i = 0; i < nrfpB && !bDiff; i++)
506     {
507         bDiff = !equal_real(
508                 ip1.generic.buf[p0 + i], ip1.generic.buf[nrfpA + i], relativeTolerance, absoluteTolerance);
509     }
510     if (bDiff)
511     {
512         fprintf(fp, "%s: ", s);
513         pr_iparams(fp, ft, ip1);
514     }
515 }
516
517 static void cmp_cmap(FILE* fp, const gmx_cmap_t* cmap1, const gmx_cmap_t* cmap2, real relativeTolerance, real absoluteTolerance)
518 {
519     int cmap1_ngrid = (cmap1 ? cmap1->cmapdata.size() : 0);
520     int cmap2_ngrid = (cmap2 ? cmap2->cmapdata.size() : 0);
521
522     cmp_int(fp, "cmap ngrid", -1, cmap1_ngrid, cmap2_ngrid);
523
524     if (cmap1 == nullptr || cmap2 == nullptr)
525     {
526         return;
527     }
528
529     cmp_int(fp, "cmap grid_spacing", -1, cmap1->grid_spacing, cmap2->grid_spacing);
530     if (cmap1->cmapdata.size() == cmap2->cmapdata.size() && cmap1->grid_spacing == cmap2->grid_spacing)
531     {
532         for (size_t g = 0; g < cmap1->cmapdata.size(); g++)
533         {
534             int i;
535
536             fprintf(fp, "comparing cmap %zu\n", g);
537
538             for (i = 0; i < 4 * cmap1->grid_spacing * cmap1->grid_spacing; i++)
539             {
540                 cmp_real(fp, "", i, cmap1->cmapdata[g].cmap[i], cmap2->cmapdata[g].cmap[i], relativeTolerance, absoluteTolerance);
541             }
542         }
543     }
544 }
545
546 static void cmp_listoflists(FILE*                        fp,
547                             const gmx::ListOfLists<int>& list1,
548                             const gmx::ListOfLists<int>& list2,
549                             const char*                  s)
550 {
551     char buf[32];
552
553     fprintf(fp, "comparing blocka %s\n", s);
554     sprintf(buf, "%s.numLists", s);
555     cmp_int(fp, buf, -1, list1.ssize(), list2.ssize());
556     sprintf(buf, "%s.numElements", s);
557     cmp_int(fp, buf, -1, list1.numElements(), list2.numElements());
558 }
559
560 static void compareFfparams(FILE*                 fp,
561                             const gmx_ffparams_t& ff1,
562                             const gmx_ffparams_t& ff2,
563                             real                  relativeTolerance,
564                             real                  absoluteTolerance)
565 {
566     fprintf(fp, "comparing force field parameters\n");
567     cmp_int(fp, "numTypes", -1, ff1.numTypes(), ff2.numTypes());
568     cmp_int(fp, "atnr", -1, ff1.atnr, ff1.atnr);
569     cmp_double(fp, "reppow", -1, ff1.reppow, ff2.reppow, relativeTolerance, absoluteTolerance);
570     cmp_real(fp, "fudgeQQ", -1, ff1.fudgeQQ, ff2.fudgeQQ, relativeTolerance, absoluteTolerance);
571     cmp_cmap(fp, &ff1.cmap_grid, &ff2.cmap_grid, relativeTolerance, absoluteTolerance);
572     for (int i = 0; i < std::min(ff1.numTypes(), ff2.numTypes()); i++)
573     {
574         std::string buf = gmx::formatString("ffparams->functype[%d]", i);
575         cmp_int(fp, buf.c_str(), i, ff1.functype[i], ff2.functype[i]);
576         buf = gmx::formatString("ffparams->iparams[%d]", i);
577         cmp_iparm(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], ff2.iparams[i], relativeTolerance, absoluteTolerance);
578     }
579 }
580
581 static void compareFfparamAB(FILE* fp, const gmx_ffparams_t& ff1, real relativeTolerance, real absoluteTolerance)
582 {
583     fprintf(fp, "comparing free energy parameters\n");
584     for (int i = 0; i < ff1.numTypes(); i++)
585     {
586         std::string buf = gmx::formatString("ffparams->iparams[%d]", i);
587         cmp_iparm_AB(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], relativeTolerance, absoluteTolerance);
588     }
589 }
590 static void compareInteractionLists(FILE* fp, const InteractionLists* il1, const InteractionLists* il2)
591 {
592     fprintf(fp, "comparing InteractionLists\n");
593     if ((il1 || il2) && (!il1 || !il2))
594     {
595         fprintf(fp, "InteractionLists are present in topology %d but not in the other\n", il1 ? 1 : 2);
596     }
597     if (il1 && il2)
598     {
599         for (int i = 0; i < F_NRE; i++)
600         {
601             cmp_int(fp, "InteractionList size", i, il1->at(i).size(), il2->at(i).size());
602             int nr = std::min(il1->at(i).size(), il2->at(i).size());
603             for (int j = 0; j < nr; j++)
604             {
605                 cmp_int(fp, "InteractionList entry", j, il1->at(i).iatoms.at(j), il2->at(i).iatoms.at(j));
606             }
607         }
608     }
609 }
610
611 static void compareMoltypes(FILE*                              fp,
612                             gmx::ArrayRef<const gmx_moltype_t> mt1,
613                             gmx::ArrayRef<const gmx_moltype_t> mt2,
614                             real                               relativeTolerance,
615                             real                               absoluteTolerance)
616 {
617     fprintf(fp, "comparing molecule types\n");
618     cmp_int(fp, "moltype size", -1, mt1.size(), mt2.size());
619     for (int i = 0; i < std::min(mt1.ssize(), mt2.ssize()); i++)
620     {
621         cmp_str(fp, "Name", i, *mt1[i].name, *mt2[i].name);
622         compareAtoms(fp, &mt1[i].atoms, &mt2[i].atoms, relativeTolerance, absoluteTolerance);
623         compareInteractionLists(fp, &mt1[i].ilist, &mt2[i].ilist);
624         std::string buf = gmx::formatString("excls[%d]", i);
625         cmp_listoflists(fp, mt1[i].excls, mt2[i].excls, buf.c_str());
626     }
627 }
628
629 static void compareMoletypeAB(FILE* fp, gmx::ArrayRef<const gmx_moltype_t> mt1, real relativeTolerance, real absoluteTolerance)
630 {
631     fprintf(fp, "comparing free energy molecule types\n");
632     for (gmx::index i = 0; i < mt1.ssize(); i++)
633     {
634         compareAtoms(fp, &mt1[i].atoms, nullptr, relativeTolerance, absoluteTolerance);
635     }
636 }
637 static void compareMolblocks(FILE*                               fp,
638                              gmx::ArrayRef<const gmx_molblock_t> mb1,
639                              gmx::ArrayRef<const gmx_molblock_t> mb2)
640 {
641     fprintf(fp, "comparing molecule blocks\n");
642     cmp_int(fp, "molblock size", -1, mb1.size(), mb2.size());
643     int nr = std::min(mb1.size(), mb2.size());
644     for (int i = 0; i < nr; i++)
645     {
646         cmp_int(fp, "type", i, mb1[i].type, mb2[i].type);
647         cmp_int(fp, "nmol", i, mb1[i].nmol, mb2[i].nmol);
648         // Only checking size of restraint vectors for now
649         cmp_int(fp, "posres_xA size", i, mb1[i].posres_xA.size(), mb2[i].posres_xA.size());
650         cmp_int(fp, "posres_xB size", i, mb1[i].posres_xB.size(), mb2[i].posres_xB.size());
651     }
652 }
653
654 static void compareAtomtypes(FILE* fp, const t_atomtypes& at1, const t_atomtypes& at2)
655 {
656     fprintf(fp, "comparing atomtypes\n");
657     cmp_int(fp, "nr", -1, at1.nr, at2.nr);
658     int nr = std::min(at1.nr, at2.nr);
659     for (int i = 0; i < nr; i++)
660     {
661         cmp_int(fp, "atomtype", i, at1.atomnumber[i], at2.atomnumber[i]);
662     }
663 }
664
665 static void compareIntermolecularExclusions(FILE*                    fp,
666                                             gmx::ArrayRef<const int> ime1,
667                                             gmx::ArrayRef<const int> ime2)
668 {
669     fprintf(fp, "comparing intermolecular exclusions\n");
670     cmp_int(fp, "exclusion number", -1, ime1.size(), ime2.size());
671     int nr = std::min(ime1.size(), ime2.size());
672     for (int i = 0; i < nr; i++)
673     {
674         cmp_int(fp, "exclusion", i, ime1[i], ime2[i]);
675     }
676 }
677
678 static void compareBlockIndices(FILE*                                     fp,
679                                 gmx::ArrayRef<const MoleculeBlockIndices> mbi1,
680                                 gmx::ArrayRef<const MoleculeBlockIndices> mbi2)
681 {
682     fprintf(fp, "comparing moleculeBlockIndices\n");
683     cmp_int(fp, "size", -1, mbi1.size(), mbi2.size());
684     int nr = std::min(mbi1.size(), mbi2.size());
685     for (int i = 0; i < nr; i++)
686     {
687         cmp_int(fp, "numAtomsPerMolecule", i, mbi1[i].numAtomsPerMolecule, mbi2[i].numAtomsPerMolecule);
688         cmp_int(fp, "globalAtomStart", i, mbi1[i].globalAtomStart, mbi2[i].globalAtomStart);
689         cmp_int(fp, "globalAtomEnd", i, mbi1[i].globalAtomEnd, mbi2[i].globalAtomEnd);
690         cmp_int(fp, "globalResidueStart", i, mbi1[i].globalResidueStart, mbi2[i].globalResidueStart);
691         cmp_int(fp, "moleculeIndexStart", i, mbi1[i].moleculeIndexStart, mbi2[i].moleculeIndexStart);
692     }
693 }
694
695 void compareMtop(FILE* fp, const gmx_mtop_t& mtop1, const gmx_mtop_t& mtop2, real relativeTolerance, real absoluteTolerance)
696 {
697     fprintf(fp, "comparing mtop topology\n");
698     cmp_str(fp, "Name", -1, *mtop1.name, *mtop2.name);
699     cmp_int(fp, "natoms", -1, mtop1.natoms, mtop2.natoms);
700     cmp_int(fp,
701             "maxres_renum",
702             -1,
703             mtop1.maxResiduesPerMoleculeToTriggerRenumber(),
704             mtop2.maxResiduesPerMoleculeToTriggerRenumber());
705     cmp_int(fp, "maxresnr", -1, mtop1.maxResNumberNotRenumbered(), mtop2.maxResNumberNotRenumbered());
706     cmp_bool(fp, "bIntermolecularInteractions", -1, mtop1.bIntermolecularInteractions, mtop2.bIntermolecularInteractions);
707     cmp_bool(fp, "haveMoleculeIndices", -1, mtop1.haveMoleculeIndices, mtop2.haveMoleculeIndices);
708
709     compareFfparams(fp, mtop1.ffparams, mtop2.ffparams, relativeTolerance, absoluteTolerance);
710     compareMoltypes(fp, mtop1.moltype, mtop2.moltype, relativeTolerance, absoluteTolerance);
711     compareMolblocks(fp, mtop1.molblock, mtop2.molblock);
712     compareInteractionLists(fp, mtop1.intermolecular_ilist.get(), mtop2.intermolecular_ilist.get());
713     compareAtomtypes(fp, mtop1.atomtypes, mtop2.atomtypes);
714     compareAtomGroups(fp, mtop1.groups, mtop2.groups, mtop1.natoms, mtop2.natoms);
715     compareIntermolecularExclusions(
716             fp, mtop1.intermolecularExclusionGroup, mtop2.intermolecularExclusionGroup);
717     compareBlockIndices(fp, mtop1.moleculeBlockIndices, mtop2.moleculeBlockIndices);
718 }
719
720 void compareMtopAB(FILE* fp, const gmx_mtop_t& mtop1, real relativeTolerance, real absoluteTolerance)
721 {
722     fprintf(fp, "comparing topAB\n");
723     compareFfparamAB(fp, mtop1.ffparams, relativeTolerance, absoluteTolerance);
724     compareMoletypeAB(fp, mtop1.moltype, relativeTolerance, absoluteTolerance);
725 }
726
727 void compareAtomGroups(FILE* fp, const SimulationGroups& g0, const SimulationGroups& g1, int natoms0, int natoms1)
728 {
729     fprintf(fp, "comparing groups\n");
730
731     for (auto group : keysOf(g0.groups))
732     {
733         std::string buf = gmx::formatString("grps[%d].nr", static_cast<int>(group));
734         cmp_int(fp, buf.c_str(), -1, g0.groups[group].size(), g1.groups[group].size());
735         if (g0.groups[group].size() == g1.groups[group].size())
736         {
737             for (gmx::index j = 0; j < gmx::ssize(g0.groups[group]); j++)
738             {
739                 buf = gmx::formatString("grps[%d].name[%zd]", static_cast<int>(group), j);
740                 cmp_str(fp,
741                         buf.c_str(),
742                         -1,
743                         *g0.groupNames[g0.groups[group][j]],
744                         *g1.groupNames[g1.groups[group][j]]);
745             }
746         }
747         cmp_int(fp,
748                 "ngrpnr",
749                 static_cast<int>(group),
750                 g0.numberOfGroupNumbers(group),
751                 g1.numberOfGroupNumbers(group));
752         if (g0.numberOfGroupNumbers(group) == g1.numberOfGroupNumbers(group) && natoms0 == natoms1
753             && (!g0.groupNumbers[group].empty() || !g1.groupNumbers[group].empty()))
754         {
755             for (int j = 0; j < natoms0; j++)
756             {
757                 cmp_int(fp,
758                         c_simulationAtomGroupTypeShortNames[group],
759                         j,
760                         getGroupType(g0, group, j),
761                         getGroupType(g1, group, j));
762             }
763         }
764     }
765     /* We have compared the names in the groups lists,
766      * so we can skip the grpname list comparison.
767      */
768 }
769
770 int getGroupType(const SimulationGroups& group, SimulationAtomGroupType type, int atom)
771 {
772     return (group.groupNumbers[type].empty() ? 0 : group.groupNumbers[type][atom]);
773 }
774
775 void copy_moltype(const gmx_moltype_t* src, gmx_moltype_t* dst)
776 {
777     dst->name          = src->name;
778     dst->excls         = src->excls;
779     t_atoms* atomsCopy = copy_t_atoms(&src->atoms);
780     dst->atoms         = *atomsCopy;
781     sfree(atomsCopy);
782
783     for (int i = 0; i < F_NRE; ++i)
784     {
785         dst->ilist[i] = src->ilist[i];
786     }
787 }