Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / topology / topology.cpp
blobd43ae865b6597eaa1cadfb1d3d7d19cf7d2b8eb3
1 /*
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,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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "topology.h"
41 #include <cstdio>
43 #include <algorithm>
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"
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 } };
72 const char *shortName(SimulationAtomGroupType type)
74 return c_simulationAtomGroupTypeShortNames[type];
77 void init_top(t_topology *top)
79 top->name = nullptr;
80 init_idef(&top->idef);
81 init_atom(&(top->atoms));
82 init_atomtypes(&(top->atomtypes));
83 init_block(&top->cgs);
84 init_block(&top->mols);
85 init_blocka(&top->excls);
86 open_symtab(&top->symtab);
90 gmx_moltype_t::gmx_moltype_t() :
91 name(nullptr),
92 cgs(),
93 excls()
95 init_t_atoms(&atoms, 0, FALSE);
98 gmx_moltype_t::~gmx_moltype_t()
100 done_atom(&atoms);
101 done_block(&cgs);
102 done_blocka(&excls);
105 gmx_mtop_t::gmx_mtop_t()
107 init_atomtypes(&atomtypes);
108 open_symtab(&symtab);
111 gmx_mtop_t::~gmx_mtop_t()
113 done_symtab(&symtab);
115 moltype.clear();
116 molblock.clear();
117 done_atomtypes(&atomtypes);
120 void done_top(t_topology *top)
122 done_idef(&top->idef);
123 done_atom(&(top->atoms));
125 /* For GB */
126 done_atomtypes(&(top->atomtypes));
128 done_symtab(&(top->symtab));
129 done_block(&(top->cgs));
130 done_block(&(top->mols));
131 done_blocka(&(top->excls));
134 void done_top_mtop(t_topology *top, gmx_mtop_t *mtop)
136 if (mtop != nullptr)
138 if (top != nullptr)
140 done_idef(&top->idef);
141 done_atom(&top->atoms);
142 done_block(&top->cgs);
143 done_blocka(&top->excls);
144 done_block(&top->mols);
145 done_symtab(&top->symtab);
146 open_symtab(&mtop->symtab);
147 done_atomtypes(&(top->atomtypes));
150 // Note that the rest of mtop will be freed by the destructor
154 gmx_localtop_t::gmx_localtop_t()
156 init_blocka_null(&excls);
157 init_idef(&idef);
158 init_atomtypes(&atomtypes);
161 gmx_localtop_t::~gmx_localtop_t()
163 if (!useInDomainDecomp_)
165 done_idef(&idef);
166 done_blocka(&excls);
167 done_atomtypes(&atomtypes);
171 bool gmx_mtop_has_masses(const gmx_mtop_t *mtop)
173 if (mtop == nullptr)
175 return false;
177 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveMass;
180 bool gmx_mtop_has_charges(const gmx_mtop_t *mtop)
182 if (mtop == nullptr)
184 return false;
186 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveCharge;
189 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t &mtop)
191 for (const gmx_moltype_t &moltype : mtop.moltype)
193 const t_atoms &atoms = moltype.atoms;
194 if (atoms.haveBState)
196 for (int a = 0; a < atoms.nr; a++)
198 if (atoms.atom[a].q != atoms.atom[a].qB)
200 return true;
205 return false;
208 bool gmx_mtop_has_atomtypes(const gmx_mtop_t *mtop)
210 if (mtop == nullptr)
212 return false;
214 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveType;
217 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t *mtop)
219 if (mtop == nullptr)
221 return false;
223 return mtop->moltype.empty() || mtop->moltype[0].atoms.havePdbInfo;
226 static void pr_grps(FILE *fp,
227 const char *title,
228 gmx::ArrayRef<const AtomGroupIndices> grps,
229 char ***grpname)
231 int index = 0;
232 for (const auto &group : grps)
234 fprintf(fp, "%s[%-12s] nr=%zu, name=[", title, c_simulationAtomGroupTypeShortNames[index], group.size());
235 for (const auto &entry : group)
237 fprintf(fp, " %s", *(grpname[entry]));
239 fprintf(fp, "]\n");
240 index++;
244 static void pr_groups(FILE *fp, int indent,
245 const SimulationGroups &groups,
246 gmx_bool bShowNumbers)
248 pr_grps(fp,
249 "grp",
250 groups.groups,
251 const_cast<char ***>(groups.groupNames.data()));
252 pr_strings(fp, indent, "grpname", const_cast<char ***>(groups.groupNames.data()), groups.groupNames.size(), bShowNumbers);
254 pr_indent(fp, indent);
255 fprintf(fp, "groups ");
256 for (const auto &group : c_simulationAtomGroupTypeShortNames)
258 printf(" %5.5s", group);
260 printf("\n");
262 pr_indent(fp, indent);
263 fprintf(fp, "allocated ");
264 int nat_max = 0;
265 for (auto group : keysOf(groups.groups))
267 printf(" %5d", groups.numberOfGroupNumbers(group));
268 nat_max = std::max(nat_max, groups.numberOfGroupNumbers(group));
270 printf("\n");
272 if (nat_max == 0)
274 pr_indent(fp, indent);
275 fprintf(fp, "groupnr[%5s] =", "*");
276 for (auto gmx_unused group : keysOf(groups.groups))
278 fprintf(fp, " %3d ", 0);
280 fprintf(fp, "\n");
282 else
284 for (int i = 0; i < nat_max; i++)
286 pr_indent(fp, indent);
287 fprintf(fp, "groupnr[%5d] =", i);
288 for (auto group : keysOf(groups.groups))
290 fprintf(fp, " %3d ",
291 !groups.groupNumbers[group].empty() ?
292 groups.groupNumbers[group][i] : 0);
294 fprintf(fp, "\n");
299 static void pr_moltype(FILE *fp, int indent, const char *title,
300 const gmx_moltype_t *molt, int n,
301 const gmx_ffparams_t *ffparams,
302 gmx_bool bShowNumbers, gmx_bool bShowParameters)
304 int j;
306 indent = pr_title_n(fp, indent, title, n);
307 pr_indent(fp, indent);
308 fprintf(fp, "name=\"%s\"\n", *(molt->name));
309 pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
310 pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
311 pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
312 for (j = 0; (j < F_NRE); j++)
314 pr_ilist(fp, indent, interaction_function[j].longname,
315 ffparams->functype.data(), molt->ilist[j],
316 bShowNumbers, bShowParameters, ffparams->iparams.data());
320 static void pr_molblock(FILE *fp, int indent, const char *title,
321 const gmx_molblock_t *molb, int n,
322 const std::vector<gmx_moltype_t> &molt)
324 indent = pr_title_n(fp, indent, title, n);
325 pr_indent(fp, indent);
326 fprintf(fp, "%-20s = %d \"%s\"\n",
327 "moltype", molb->type, *(molt[molb->type].name));
328 pr_int(fp, indent, "#molecules", molb->nmol);
329 pr_int(fp, indent, "#posres_xA", molb->posres_xA.size());
330 if (!molb->posres_xA.empty())
332 pr_rvecs(fp, indent, "posres_xA", as_rvec_array(molb->posres_xA.data()), molb->posres_xA.size());
334 pr_int(fp, indent, "#posres_xB", molb->posres_xB.size());
335 if (!molb->posres_xB.empty())
337 pr_rvecs(fp, indent, "posres_xB", as_rvec_array(molb->posres_xB.data()), molb->posres_xB.size());
341 void pr_mtop(FILE *fp, int indent, const char *title, const gmx_mtop_t *mtop,
342 gmx_bool bShowNumbers, gmx_bool bShowParameters)
344 if (available(fp, mtop, indent, title))
346 indent = pr_title(fp, indent, title);
347 pr_indent(fp, indent);
348 fprintf(fp, "name=\"%s\"\n", *(mtop->name));
349 pr_int(fp, indent, "#atoms", mtop->natoms);
350 pr_int(fp, indent, "#molblock", mtop->molblock.size());
351 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
353 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
355 pr_str(fp, indent, "bIntermolecularInteractions",
356 gmx::boolToString(mtop->bIntermolecularInteractions));
357 if (mtop->bIntermolecularInteractions)
359 for (int j = 0; j < F_NRE; j++)
361 pr_ilist(fp, indent, interaction_function[j].longname,
362 mtop->ffparams.functype.data(),
363 (*mtop->intermolecular_ilist)[j],
364 bShowNumbers, bShowParameters, mtop->ffparams.iparams.data());
367 pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
368 pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
369 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
371 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
372 &mtop->ffparams, bShowNumbers, bShowParameters);
374 pr_groups(fp, indent, mtop->groups, bShowNumbers);
378 void pr_top(FILE *fp, int indent, const char *title, const t_topology *top,
379 gmx_bool bShowNumbers, gmx_bool bShowParameters)
381 if (available(fp, top, indent, title))
383 indent = pr_title(fp, indent, title);
384 pr_indent(fp, indent);
385 fprintf(fp, "name=\"%s\"\n", *(top->name));
386 pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
387 pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
388 pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
389 pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
390 pr_str(fp, indent, "bIntermolecularInteractions",
391 gmx::boolToString(top->bIntermolecularInteractions));
392 pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
393 pr_idef(fp, indent, "idef", &top->idef, bShowNumbers, bShowParameters);
397 static void cmp_iparm(FILE *fp, const char *s, t_functype ft,
398 const t_iparams &ip1, const t_iparams &ip2, real relativeTolerance, real absoluteTolerance)
400 int i;
401 gmx_bool bDiff;
403 bDiff = FALSE;
404 for (i = 0; i < MAXFORCEPARAM && !bDiff; i++)
406 bDiff = !equal_real(ip1.generic.buf[i], ip2.generic.buf[i], relativeTolerance, absoluteTolerance);
408 if (bDiff)
410 fprintf(fp, "%s1: ", s);
411 pr_iparams(fp, ft, &ip1);
412 fprintf(fp, "%s2: ", s);
413 pr_iparams(fp, ft, &ip2);
417 static void cmp_iparm_AB(FILE *fp, const char *s, t_functype ft,
418 const t_iparams &ip1, real relativeTolerance, real absoluteTolerance)
420 int nrfpA, nrfpB, p0, i;
421 gmx_bool bDiff;
423 /* Normally the first parameter is perturbable */
424 p0 = 0;
425 nrfpA = interaction_function[ft].nrfpA;
426 nrfpB = interaction_function[ft].nrfpB;
427 if (ft == F_PDIHS)
429 nrfpB = 2;
431 else if (interaction_function[ft].flags & IF_TABULATED)
433 /* For tabulated interactions only the second parameter is perturbable */
434 p0 = 1;
435 nrfpB = 1;
437 bDiff = FALSE;
438 for (i = 0; i < nrfpB && !bDiff; i++)
440 bDiff = !equal_real(ip1.generic.buf[p0+i], ip1.generic.buf[nrfpA+i], relativeTolerance, absoluteTolerance);
442 if (bDiff)
444 fprintf(fp, "%s: ", s);
445 pr_iparams(fp, ft, &ip1);
449 static void cmp_cmap(FILE *fp, const gmx_cmap_t *cmap1, const gmx_cmap_t *cmap2, real relativeTolerance, real absoluteTolerance)
451 int cmap1_ngrid = (cmap1 ? cmap1->cmapdata.size() : 0);
452 int cmap2_ngrid = (cmap2 ? cmap2->cmapdata.size() : 0);
454 cmp_int(fp, "cmap ngrid", -1, cmap1_ngrid, cmap2_ngrid);
456 if (cmap1 == nullptr || cmap2 == nullptr)
458 return;
461 cmp_int(fp, "cmap grid_spacing", -1, cmap1->grid_spacing, cmap2->grid_spacing);
462 if (cmap1->cmapdata.size() == cmap2->cmapdata.size() &&
463 cmap1->grid_spacing == cmap2->grid_spacing)
465 for (size_t g = 0; g < cmap1->cmapdata.size(); g++)
467 int i;
469 fprintf(fp, "comparing cmap %zu\n", g);
471 for (i = 0; i < 4*cmap1->grid_spacing*cmap1->grid_spacing; i++)
473 cmp_real(fp, "", i, cmap1->cmapdata[g].cmap[i], cmap2->cmapdata[g].cmap[i], relativeTolerance, absoluteTolerance);
479 static void cmp_block(FILE *fp, const t_block *b1, const t_block *b2, const char *s)
481 char buf[32];
483 fprintf(fp, "comparing block %s\n", s);
484 sprintf(buf, "%s.nr", s);
485 cmp_int(fp, buf, -1, b1->nr, b2->nr);
488 static void cmp_blocka(FILE *fp, const t_blocka *b1, const t_blocka *b2, const char *s)
490 char buf[32];
492 fprintf(fp, "comparing blocka %s\n", s);
493 sprintf(buf, "%s.nr", s);
494 cmp_int(fp, buf, -1, b1->nr, b2->nr);
495 sprintf(buf, "%s.nra", s);
496 cmp_int(fp, buf, -1, b1->nra, b2->nra);
499 static void compareFfparams(FILE *fp, const gmx_ffparams_t &ff1, const gmx_ffparams_t &ff2, real relativeTolerance, real absoluteTolerance)
501 fprintf(fp, "comparing force field parameters\n");
502 cmp_int(fp, "numTypes", -1, ff1.numTypes(), ff2.numTypes());
503 cmp_int(fp, "atnr", -1, ff1.atnr, ff1.atnr);
504 cmp_double(fp, "reppow", -1, ff1.reppow, ff2.reppow, relativeTolerance, absoluteTolerance);
505 cmp_real(fp, "fudgeQQ", -1, ff1.fudgeQQ, ff2.fudgeQQ, relativeTolerance, absoluteTolerance);
506 cmp_cmap(fp, &ff1.cmap_grid, &ff2.cmap_grid, relativeTolerance, absoluteTolerance);
507 for (int i = 0; i < std::min(ff1.numTypes(), ff2.numTypes()); i++)
509 std::string buf = gmx::formatString("ffparams->functype[%d]", i);
510 cmp_int(fp, buf.c_str(), i, ff1.functype[i], ff2.functype[i]);
511 buf = gmx::formatString("ffparams->iparams[%d]", i);
512 cmp_iparm(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], ff2.iparams[i], relativeTolerance, absoluteTolerance);
517 static void compareFfparamAB(FILE *fp, const gmx_ffparams_t &ff1, real relativeTolerance, real absoluteTolerance)
519 fprintf(fp, "comparing free energy parameters\n");
520 for (int i = 0; i < ff1.numTypes(); i++)
522 std::string buf = gmx::formatString("ffparams->iparams[%d]", i);
523 cmp_iparm_AB(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], relativeTolerance, absoluteTolerance);
526 static void compareInteractionLists(FILE *fp, const InteractionLists *il1, const InteractionLists *il2)
528 fprintf(fp, "comparing InteractionLists\n");
529 if ((il1 || il2) && (!il1 || !il2))
531 fprintf(fp, "InteractionLists are present in topology %d but not in the other\n", il1 ? 1 : 2);
533 if (il1 && il2)
535 for (int i = 0; i < F_NRE; i++)
537 cmp_int(fp, "InteractionList size", i, il1->at(i).size(), il2->at(i).size());
538 int nr = std::min(il1->at(i).size(), il2->at(i).size());
539 for (int j = 0; j < nr; j++)
541 cmp_int(fp, "InteractionList entry", j, il1->at(i).iatoms.at(j), il2->at(i).iatoms.at(j));
547 static void compareMoltypes(FILE *fp, gmx::ArrayRef<const gmx_moltype_t> mt1, gmx::ArrayRef<const gmx_moltype_t> mt2, real relativeTolerance, real absoluteTolerance)
549 fprintf(fp, "comparing molecule types\n");
550 cmp_int(fp, "moltype size", -1, mt1.size(), mt2.size());
551 for (int i = 0; i < std::min(mt1.ssize(), mt2.ssize()); i++)
553 cmp_str(fp, "Name", i, *mt1[i].name, *mt2[i].name);
554 compareAtoms(fp, &mt1[i].atoms, &mt2[i].atoms, relativeTolerance, absoluteTolerance);
555 compareInteractionLists(fp, &mt1[i].ilist, &mt2[i].ilist);
556 std::string buf = gmx::formatString("cgs[%d]", i);
557 cmp_block(fp, &mt1[i].cgs, &mt2[i].cgs, buf.c_str());
558 buf = gmx::formatString("excls[%d]", i);
559 cmp_blocka(fp, &mt1[i].excls, &mt2[i].excls, buf.c_str());
563 static void compareMoletypeAB(FILE *fp, gmx::ArrayRef<const gmx_moltype_t> mt1, real relativeTolerance, real absoluteTolerance)
565 fprintf(fp, "comparing free energy molecule types\n");
566 for (int i = 0; i < mt1.ssize(); i++)
568 compareAtoms(fp, &mt1[i].atoms, nullptr, relativeTolerance, absoluteTolerance);
571 static void compareMolblocks(FILE *fp, gmx::ArrayRef<const gmx_molblock_t> mb1, gmx::ArrayRef<const gmx_molblock_t> mb2)
573 fprintf(fp, "comparing molecule blocks\n");
574 cmp_int(fp, "molblock size", -1, mb1.size(), mb2.size());
575 int nr = std::min(mb1.size(), mb2.size());
576 for (int i = 0; i < nr; i++)
578 cmp_int(fp, "type", i, mb1[i].type, mb2[i].type);
579 cmp_int(fp, "nmol", i, mb1[i].nmol, mb2[i].nmol);
580 // Only checking size of restraint vectors for now
581 cmp_int(fp, "posres_xA size", i, mb1[i].posres_xA.size(), mb2[i].posres_xA.size());
582 cmp_int(fp, "posres_xB size", i, mb1[i].posres_xB.size(), mb2[i].posres_xB.size());
587 static void compareAtomtypes(FILE *fp, const t_atomtypes &at1, const t_atomtypes &at2)
589 fprintf(fp, "comparing atomtypes\n");
590 cmp_int(fp, "nr", -1, at1.nr, at2.nr);
591 int nr = std::min(at1.nr, at2.nr);
592 for (int i = 0; i < nr; i++)
594 cmp_int(fp, "atomtype", i, at1.atomnumber[i], at2.atomnumber[i]);
598 static void compareIntermolecularExclusions(FILE *fp, gmx::ArrayRef<const int> ime1, gmx::ArrayRef<const int> ime2)
600 fprintf(fp, "comparing intermolecular exclusions\n");
601 cmp_int(fp, "exclusion number", -1, ime1.size(), ime2.size());
602 int nr = std::min(ime1.size(), ime2.size());
603 for (int i = 0; i < nr; i++)
605 cmp_int(fp, "exclusion", i, ime1[i], ime2[i]);
609 static void compareBlockIndices(FILE *fp, gmx::ArrayRef<const MoleculeBlockIndices> mbi1, gmx::ArrayRef<const MoleculeBlockIndices> mbi2)
611 fprintf(fp, "comparing moleculeBlockIndices\n");
612 cmp_int(fp, "size", -1, mbi1.size(), mbi2.size());
613 int nr = std::min(mbi1.size(), mbi2.size());
614 for (int i = 0; i < nr; i++)
616 cmp_int(fp, "numAtomsPerMolecule", i, mbi1[i].numAtomsPerMolecule, mbi2[i].numAtomsPerMolecule);
617 cmp_int(fp, "globalAtomStart", i, mbi1[i].globalAtomStart, mbi2[i].globalAtomStart);
618 cmp_int(fp, "globalAtomEnd", i, mbi1[i].globalAtomEnd, mbi2[i].globalAtomEnd);
619 cmp_int(fp, "globalResidueStart", i, mbi1[i].globalResidueStart, mbi2[i].globalResidueStart);
620 cmp_int(fp, "moleculeIndexStart", i, mbi1[i].moleculeIndexStart, mbi2[i].moleculeIndexStart);
624 void compareMtop(FILE *fp, const gmx_mtop_t &mtop1, const gmx_mtop_t &mtop2, real relativeTolerance, real absoluteTolerance)
626 fprintf(fp, "comparing mtop topology\n");
627 cmp_str(fp, "Name", -1, *mtop1.name, *mtop2.name);
628 cmp_int(fp, "natoms", -1, mtop1.natoms, mtop2.natoms);
629 cmp_int(fp, "maxres_renum", -1, mtop1.maxres_renum, mtop2.maxres_renum);
630 cmp_int(fp, "maxresnr", -1, mtop1.maxresnr, mtop2.maxresnr);
631 cmp_bool(fp, "bIntermolecularInteractions", -1, mtop1.bIntermolecularInteractions, mtop2.bIntermolecularInteractions);
632 cmp_bool(fp, "haveMoleculeIndices", -1, mtop1.haveMoleculeIndices, mtop2.haveMoleculeIndices);
634 compareFfparams(fp, mtop1.ffparams, mtop2.ffparams, relativeTolerance, absoluteTolerance);
635 compareMoltypes(fp, mtop1.moltype, mtop2.moltype, relativeTolerance, absoluteTolerance);
636 compareMolblocks(fp, mtop1.molblock, mtop2.molblock);
637 compareInteractionLists(fp, mtop1.intermolecular_ilist.get(), mtop2.intermolecular_ilist.get());
638 compareAtomtypes(fp, mtop1.atomtypes, mtop2.atomtypes);
639 compareAtomGroups(fp, mtop1.groups, mtop2.groups, mtop1.natoms, mtop2.natoms);
640 compareIntermolecularExclusions(fp, mtop1.intermolecularExclusionGroup, mtop2.intermolecularExclusionGroup);
641 compareBlockIndices(fp, mtop1.moleculeBlockIndices, mtop2.moleculeBlockIndices);
644 void compareMtopAB(FILE *fp, const gmx_mtop_t &mtop1, real relativeTolerance, real absoluteTolerance)
646 fprintf(fp, "comparing topAB\n");
647 compareFfparamAB(fp, mtop1.ffparams, relativeTolerance, absoluteTolerance);
648 compareMoletypeAB(fp, mtop1.moltype, relativeTolerance, absoluteTolerance);
651 void compareAtomGroups(FILE *fp, const SimulationGroups &g0, const SimulationGroups &g1,
652 int natoms0, int natoms1)
654 fprintf(fp, "comparing groups\n");
656 for (auto group : keysOf(g0.groups))
658 std::string buf = gmx::formatString("grps[%d].nr", static_cast<int>(group));
659 cmp_int(fp, buf.c_str(), -1, g0.groups[group].size(), g1.groups[group].size());
660 if (g0.groups[group].size() == g1.groups[group].size())
662 for (int j = 0; j < gmx::ssize(g0.groups[group]); j++)
664 buf = gmx::formatString("grps[%d].name[%d]", static_cast<int>(group), j);
665 cmp_str(fp, buf.c_str(), -1,
666 *g0.groupNames[g0.groups[group][j]],
667 *g1.groupNames[g1.groups[group][j]]);
670 cmp_int(fp, "ngrpnr", static_cast<int>(group), g0.numberOfGroupNumbers(group), g1.numberOfGroupNumbers(group));
671 if (g0.numberOfGroupNumbers(group) == g1.numberOfGroupNumbers(group) && natoms0 == natoms1 &&
672 (!g0.groupNumbers[group].empty() || !g1.groupNumbers[group].empty()))
674 for (int j = 0; j < natoms0; j++)
676 cmp_int(fp, c_simulationAtomGroupTypeShortNames[group], j, getGroupType(g0, group, j), getGroupType(g1, group, j));
680 /* We have compared the names in the groups lists,
681 * so we can skip the grpname list comparison.
685 int getGroupType(const SimulationGroups &group, SimulationAtomGroupType type, int atom)
687 return (group.groupNumbers[type].empty() ? 0 : group.groupNumbers[type][atom]);
690 void copy_moltype(const gmx_moltype_t *src, gmx_moltype_t *dst)
692 dst->name = src->name;
693 copy_blocka(&src->excls, &dst->excls);
694 copy_block(&src->cgs, &dst->cgs);
695 t_atoms *atomsCopy = copy_t_atoms(&src->atoms);
696 dst->atoms = *atomsCopy;
697 sfree(atomsCopy);
699 for (int i = 0; i < F_NRE; ++i)
701 dst->ilist[i] = src->ilist[i];