2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
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"
58 static gmx::EnumerationArray<SimulationAtomGroupType, const char*> c_simulationAtomGroupTypeShortNames = {
71 const char* shortName(SimulationAtomGroupType type)
73 return c_simulationAtomGroupTypeShortNames[type];
76 void init_top(t_topology* top)
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);
87 gmx_moltype_t::gmx_moltype_t() : name(nullptr)
89 init_t_atoms(&atoms, 0, FALSE);
92 gmx_moltype_t::~gmx_moltype_t()
97 static int gmx_mtop_maxresnr(const gmx::ArrayRef<const gmx_moltype_t> moltypes, int maxres_renum)
101 for (const gmx_moltype_t& moltype : moltypes)
103 const t_atoms& atoms = moltype.atoms;
104 if (atoms.nres > maxres_renum)
106 for (int r = 0; r < atoms.nres; r++)
108 if (atoms.resinfo[r].nr > maxresnr)
110 maxresnr = atoms.resinfo[r].nr;
119 gmx_mtop_t::gmx_mtop_t()
121 init_atomtypes(&atomtypes);
122 open_symtab(&symtab);
125 gmx_mtop_t::~gmx_mtop_t()
127 done_symtab(&symtab);
131 done_atomtypes(&atomtypes);
134 void gmx_mtop_t::finalize()
136 if (molblock.size() == 1 && molblock[0].nmol == 1)
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.
142 maxResiduesPerMoleculeToTriggerRenumber_ = 0;
146 /* We only renumber single residue molecules. Their intra-molecular
147 * residue numbering is anyhow irrelevant.
149 maxResiduesPerMoleculeToTriggerRenumber_ = 1;
152 const char* env = getenv("GMX_MAXRESRENUM");
155 sscanf(env, "%d", &maxResiduesPerMoleculeToTriggerRenumber_);
157 if (maxResiduesPerMoleculeToTriggerRenumber_ == -1)
159 /* -1 signals renumber residues in all molecules */
160 maxResiduesPerMoleculeToTriggerRenumber_ = std::numeric_limits<int>::max();
163 maxResNumberNotRenumbered_ = gmx_mtop_maxresnr(moltype, maxResiduesPerMoleculeToTriggerRenumber_);
165 buildMolblockIndices();
168 void gmx_mtop_t::buildMolblockIndices()
170 moleculeBlockIndices.resize(molblock.size());
173 int residueIndex = 0;
174 int residueNumberStart = maxResNumberNotRenumbered_ + 1;
175 int moleculeIndexStart = 0;
176 for (size_t mb = 0; mb < molblock.size(); mb++)
178 const gmx_molblock_t& molb = molblock[mb];
179 MoleculeBlockIndices& indices = moleculeBlockIndices[mb];
180 const int numResPerMol = moltype[molb.type].atoms.nres;
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_)
191 residueNumberStart += molb.nmol * numResPerMol;
193 indices.moleculeIndexStart = moleculeIndexStart;
194 moleculeIndexStart += molb.nmol;
198 void done_top(t_topology* top)
200 done_idef(&top->idef);
201 done_atom(&(top->atoms));
204 done_atomtypes(&(top->atomtypes));
206 done_symtab(&(top->symtab));
207 done_block(&(top->mols));
210 void done_top_mtop(t_topology* top, gmx_mtop_t* mtop)
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));
224 // Note that the rest of mtop will be freed by the destructor
228 gmx_localtop_t::gmx_localtop_t(const gmx_ffparams_t& ffparams) : idef(ffparams) {}
230 bool gmx_mtop_has_masses(const gmx_mtop_t* mtop)
236 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveMass;
239 bool gmx_mtop_has_charges(const gmx_mtop_t* mtop)
245 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveCharge;
248 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t& mtop)
250 for (const gmx_moltype_t& moltype : mtop.moltype)
252 const t_atoms& atoms = moltype.atoms;
253 if (atoms.haveBState)
255 for (int a = 0; a < atoms.nr; a++)
257 if (atoms.atom[a].q != atoms.atom[a].qB)
267 bool gmx_mtop_has_atomtypes(const gmx_mtop_t* mtop)
273 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveType;
276 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t* mtop)
282 return mtop->moltype.empty() || mtop->moltype[0].atoms.havePdbInfo;
285 static void pr_grps(FILE* fp, const char* title, gmx::ArrayRef<const AtomGroupIndices> grps, char*** grpname)
288 for (const auto& group : grps)
290 fprintf(fp, "%s[%-12s] nr=%zu, name=[", title, c_simulationAtomGroupTypeShortNames[index], group.size());
291 for (const auto& entry : group)
293 fprintf(fp, " %s", *(grpname[entry]));
300 static void pr_groups(FILE* fp, int indent, const SimulationGroups& groups, gmx_bool bShowNumbers)
302 pr_grps(fp, "grp", groups.groups, const_cast<char***>(groups.groupNames.data()));
306 const_cast<char***>(groups.groupNames.data()),
307 groups.groupNames.size(),
310 pr_indent(fp, indent);
311 fprintf(fp, "groups ");
312 for (const auto& group : c_simulationAtomGroupTypeShortNames)
314 printf(" %5.5s", group);
318 pr_indent(fp, indent);
319 fprintf(fp, "allocated ");
321 for (auto group : keysOf(groups.groups))
323 printf(" %5d", groups.numberOfGroupNumbers(group));
324 nat_max = std::max(nat_max, groups.numberOfGroupNumbers(group));
330 pr_indent(fp, indent);
331 fprintf(fp, "groupnr[%5s] =", "*");
332 for (auto gmx_unused group : keysOf(groups.groups))
334 fprintf(fp, " %3d ", 0);
340 for (int i = 0; i < nat_max; i++)
342 pr_indent(fp, indent);
343 fprintf(fp, "groupnr[%5d] =", i);
344 for (auto group : keysOf(groups.groups))
348 !groups.groupNumbers[group].empty() ? groups.groupNumbers[group][i] : 0);
355 static void pr_moltype(FILE* fp,
358 const gmx_moltype_t* molt,
360 const gmx_ffparams_t* ffparams,
361 gmx_bool bShowNumbers,
362 gmx_bool bShowParameters)
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++)
375 interaction_function[j].longname,
376 ffparams->functype.data(),
380 ffparams->iparams.data());
384 static void pr_molblock(FILE* fp,
387 const gmx_molblock_t* molb,
389 const std::vector<gmx_moltype_t>& molt)
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())
398 pr_rvecs(fp, indent, "posres_xA", as_rvec_array(molb->posres_xA.data()), molb->posres_xA.size());
400 pr_int(fp, indent, "#posres_xB", molb->posres_xB.size());
401 if (!molb->posres_xB.empty())
403 pr_rvecs(fp, indent, "posres_xB", as_rvec_array(molb->posres_xB.data()), molb->posres_xB.size());
407 void pr_mtop(FILE* fp, int indent, const char* title, const gmx_mtop_t* mtop, gmx_bool bShowNumbers, gmx_bool bShowParameters)
409 if (available(fp, mtop, indent, title))
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++)
418 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
420 pr_str(fp, indent, "bIntermolecularInteractions", gmx::boolToString(mtop->bIntermolecularInteractions));
421 if (mtop->bIntermolecularInteractions)
423 for (int j = 0; j < F_NRE; j++)
427 interaction_function[j].longname,
428 mtop->ffparams.functype.data(),
429 (*mtop->intermolecular_ilist)[j],
432 mtop->ffparams.iparams.data());
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++)
439 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt, &mtop->ffparams, bShowNumbers, bShowParameters);
441 pr_groups(fp, indent, mtop->groups, bShowNumbers);
445 void pr_top(FILE* fp, int indent, const char* title, const t_topology* top, gmx_bool bShowNumbers, gmx_bool bShowParameters)
447 if (available(fp, top, indent, title))
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);
460 static void cmp_iparm(FILE* fp,
463 const t_iparams& ip1,
464 const t_iparams& ip2,
465 real relativeTolerance,
466 real absoluteTolerance)
472 for (i = 0; i < MAXFORCEPARAM && !bDiff; i++)
474 bDiff = !equal_real(ip1.generic.buf[i], ip2.generic.buf[i], relativeTolerance, absoluteTolerance);
478 fprintf(fp, "%s1: ", s);
479 pr_iparams(fp, ft, ip1);
480 fprintf(fp, "%s2: ", s);
481 pr_iparams(fp, ft, ip2);
485 static void cmp_iparm_AB(FILE* fp, const char* s, t_functype ft, const t_iparams& ip1, real relativeTolerance, real absoluteTolerance)
487 int nrfpA, nrfpB, p0, i;
490 /* Normally the first parameter is perturbable */
492 nrfpA = interaction_function[ft].nrfpA;
493 nrfpB = interaction_function[ft].nrfpB;
498 else if (interaction_function[ft].flags & IF_TABULATED)
500 /* For tabulated interactions only the second parameter is perturbable */
505 for (i = 0; i < nrfpB && !bDiff; i++)
508 ip1.generic.buf[p0 + i], ip1.generic.buf[nrfpA + i], relativeTolerance, absoluteTolerance);
512 fprintf(fp, "%s: ", s);
513 pr_iparams(fp, ft, ip1);
517 static void cmp_cmap(FILE* fp, const gmx_cmap_t* cmap1, const gmx_cmap_t* cmap2, real relativeTolerance, real absoluteTolerance)
519 int cmap1_ngrid = (cmap1 ? cmap1->cmapdata.size() : 0);
520 int cmap2_ngrid = (cmap2 ? cmap2->cmapdata.size() : 0);
522 cmp_int(fp, "cmap ngrid", -1, cmap1_ngrid, cmap2_ngrid);
524 if (cmap1 == nullptr || cmap2 == nullptr)
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)
532 for (size_t g = 0; g < cmap1->cmapdata.size(); g++)
536 fprintf(fp, "comparing cmap %zu\n", g);
538 for (i = 0; i < 4 * cmap1->grid_spacing * cmap1->grid_spacing; i++)
540 cmp_real(fp, "", i, cmap1->cmapdata[g].cmap[i], cmap2->cmapdata[g].cmap[i], relativeTolerance, absoluteTolerance);
546 static void cmp_listoflists(FILE* fp,
547 const gmx::ListOfLists<int>& list1,
548 const gmx::ListOfLists<int>& list2,
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());
560 static void compareFfparams(FILE* fp,
561 const gmx_ffparams_t& ff1,
562 const gmx_ffparams_t& ff2,
563 real relativeTolerance,
564 real absoluteTolerance)
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++)
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);
581 static void compareFfparamAB(FILE* fp, const gmx_ffparams_t& ff1, real relativeTolerance, real absoluteTolerance)
583 fprintf(fp, "comparing free energy parameters\n");
584 for (int i = 0; i < ff1.numTypes(); i++)
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);
590 static void compareInteractionLists(FILE* fp, const InteractionLists* il1, const InteractionLists* il2)
592 fprintf(fp, "comparing InteractionLists\n");
593 if ((il1 || il2) && (!il1 || !il2))
595 fprintf(fp, "InteractionLists are present in topology %d but not in the other\n", il1 ? 1 : 2);
599 for (int i = 0; i < F_NRE; i++)
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++)
605 cmp_int(fp, "InteractionList entry", j, il1->at(i).iatoms.at(j), il2->at(i).iatoms.at(j));
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)
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++)
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());
629 static void compareMoletypeAB(FILE* fp, gmx::ArrayRef<const gmx_moltype_t> mt1, real relativeTolerance, real absoluteTolerance)
631 fprintf(fp, "comparing free energy molecule types\n");
632 for (gmx::index i = 0; i < mt1.ssize(); i++)
634 compareAtoms(fp, &mt1[i].atoms, nullptr, relativeTolerance, absoluteTolerance);
637 static void compareMolblocks(FILE* fp,
638 gmx::ArrayRef<const gmx_molblock_t> mb1,
639 gmx::ArrayRef<const gmx_molblock_t> mb2)
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++)
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());
654 static void compareAtomtypes(FILE* fp, const t_atomtypes& at1, const t_atomtypes& at2)
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++)
661 cmp_int(fp, "atomtype", i, at1.atomnumber[i], at2.atomnumber[i]);
665 static void compareIntermolecularExclusions(FILE* fp,
666 gmx::ArrayRef<const int> ime1,
667 gmx::ArrayRef<const int> ime2)
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++)
674 cmp_int(fp, "exclusion", i, ime1[i], ime2[i]);
678 static void compareBlockIndices(FILE* fp,
679 gmx::ArrayRef<const MoleculeBlockIndices> mbi1,
680 gmx::ArrayRef<const MoleculeBlockIndices> mbi2)
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++)
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);
695 void compareMtop(FILE* fp, const gmx_mtop_t& mtop1, const gmx_mtop_t& mtop2, real relativeTolerance, real absoluteTolerance)
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);
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);
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);
720 void compareMtopAB(FILE* fp, const gmx_mtop_t& mtop1, real relativeTolerance, real absoluteTolerance)
722 fprintf(fp, "comparing topAB\n");
723 compareFfparamAB(fp, mtop1.ffparams, relativeTolerance, absoluteTolerance);
724 compareMoletypeAB(fp, mtop1.moltype, relativeTolerance, absoluteTolerance);
727 void compareAtomGroups(FILE* fp, const SimulationGroups& g0, const SimulationGroups& g1, int natoms0, int natoms1)
729 fprintf(fp, "comparing groups\n");
731 for (auto group : keysOf(g0.groups))
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())
737 for (gmx::index j = 0; j < gmx::ssize(g0.groups[group]); j++)
739 buf = gmx::formatString("grps[%d].name[%zd]", static_cast<int>(group), j);
743 *g0.groupNames[g0.groups[group][j]],
744 *g1.groupNames[g1.groups[group][j]]);
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()))
755 for (int j = 0; j < natoms0; j++)
758 c_simulationAtomGroupTypeShortNames[group],
760 getGroupType(g0, group, j),
761 getGroupType(g1, group, j));
765 /* We have compared the names in the groups lists,
766 * so we can skip the grpname list comparison.
770 int getGroupType(const SimulationGroups& group, SimulationAtomGroupType type, int atom)
772 return (group.groupNumbers[type].empty() ? 0 : group.groupNumbers[type][atom]);
775 void copy_moltype(const gmx_moltype_t* src, gmx_moltype_t* dst)
777 dst->name = src->name;
778 dst->excls = src->excls;
779 t_atoms* atomsCopy = copy_t_atoms(&src->atoms);
780 dst->atoms = *atomsCopy;
783 for (int i = 0; i < F_NRE; ++i)
785 dst->ilist[i] = src->ilist[i];