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.
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
72 const char *shortName(SimulationAtomGroupType type)
74 return c_simulationAtomGroupTypeShortNames[type];
77 SimulationGroups::SimulationGroups()
79 for (auto group : keysOf(groups))
82 groups[group].nm_ind = nullptr;
87 void init_top(t_topology *top)
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);
100 gmx_moltype_t::gmx_moltype_t() :
105 init_t_atoms(&atoms, 0, FALSE);
108 gmx_moltype_t::~gmx_moltype_t()
115 SimulationGroups::~SimulationGroups()
117 for (auto group : keysOf(groups))
119 if (nullptr != groups[group].nm_ind)
121 sfree(groups[group].nm_ind);
122 groups[group].nm_ind = nullptr;
127 gmx_mtop_t::gmx_mtop_t()
129 init_atomtypes(&atomtypes);
130 open_symtab(&symtab);
133 gmx_mtop_t::~gmx_mtop_t()
135 done_symtab(&symtab);
139 done_atomtypes(&atomtypes);
142 void done_top(t_topology *top)
144 done_idef(&top->idef);
145 done_atom(&(top->atoms));
148 done_atomtypes(&(top->atomtypes));
150 done_symtab(&(top->symtab));
151 done_block(&(top->cgs));
152 done_block(&(top->mols));
153 done_blocka(&(top->excls));
156 void done_top_mtop(t_topology *top, gmx_mtop_t *mtop)
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));
172 // Note that the rest of mtop will be freed by the destructor
176 gmx_localtop_t::gmx_localtop_t()
178 init_blocka_null(&excls);
180 init_atomtypes(&atomtypes);
183 gmx_localtop_t::~gmx_localtop_t()
185 if (!useInDomainDecomp_)
189 done_atomtypes(&atomtypes);
193 bool gmx_mtop_has_masses(const gmx_mtop_t *mtop)
199 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveMass;
202 bool gmx_mtop_has_charges(const gmx_mtop_t *mtop)
208 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveCharge;
211 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t &mtop)
213 for (const gmx_moltype_t &moltype : mtop.moltype)
215 const t_atoms &atoms = moltype.atoms;
216 if (atoms.haveBState)
218 for (int a = 0; a < atoms.nr; a++)
220 if (atoms.atom[a].q != atoms.atom[a].qB)
230 bool gmx_mtop_has_atomtypes(const gmx_mtop_t *mtop)
236 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveType;
239 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t *mtop)
245 return mtop->moltype.empty() || mtop->moltype[0].atoms.havePdbInfo;
248 static void pr_grps(FILE *fp,
250 gmx::ArrayRef<const t_grps> grps,
254 for (const auto &group : grps)
256 fprintf(fp, "%s[%-12s] nr=%d, name=[", title, c_simulationAtomGroupTypeShortNames[index], group.nr);
257 for (int j = 0; (j < group.nr); j++)
259 fprintf(fp, " %s", *(grpname[group.nm_ind[j]]));
266 static void pr_groups(FILE *fp, int indent,
267 const SimulationGroups &groups,
268 gmx_bool bShowNumbers)
273 const_cast<char ***>(groups.groupNames.data()));
274 pr_strings(fp, indent, "grpname", const_cast<char ***>(groups.groupNames.data()), groups.groupNames.size(), bShowNumbers);
276 pr_indent(fp, indent);
277 fprintf(fp, "groups ");
278 for (const auto &group : c_simulationAtomGroupTypeShortNames)
280 printf(" %5.5s", group);
284 pr_indent(fp, indent);
285 fprintf(fp, "allocated ");
287 for (auto group : keysOf(groups.groups))
289 printf(" %5d", groups.numberOfGroupNumbers(group));
290 nat_max = std::max(nat_max, groups.numberOfGroupNumbers(group));
296 pr_indent(fp, indent);
297 fprintf(fp, "groupnr[%5s] =", "*");
298 for (auto gmx_unused group : keysOf(groups.groups))
300 fprintf(fp, " %3d ", 0);
306 for (int i = 0; i < nat_max; i++)
308 pr_indent(fp, indent);
309 fprintf(fp, "groupnr[%5d] =", i);
310 for (auto group : keysOf(groups.groups))
313 !groups.groupNumbers[group].empty() ?
314 groups.groupNumbers[group][i] : 0);
321 static void pr_moltype(FILE *fp, int indent, const char *title,
322 const gmx_moltype_t *molt, int n,
323 const gmx_ffparams_t *ffparams,
324 gmx_bool bShowNumbers, gmx_bool bShowParameters)
328 indent = pr_title_n(fp, indent, title, n);
329 pr_indent(fp, indent);
330 fprintf(fp, "name=\"%s\"\n", *(molt->name));
331 pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
332 pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
333 pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
334 for (j = 0; (j < F_NRE); j++)
336 pr_ilist(fp, indent, interaction_function[j].longname,
337 ffparams->functype.data(), molt->ilist[j],
338 bShowNumbers, bShowParameters, ffparams->iparams.data());
342 static void pr_molblock(FILE *fp, int indent, const char *title,
343 const gmx_molblock_t *molb, int n,
344 const std::vector<gmx_moltype_t> &molt)
346 indent = pr_title_n(fp, indent, title, n);
347 pr_indent(fp, indent);
348 fprintf(fp, "%-20s = %d \"%s\"\n",
349 "moltype", molb->type, *(molt[molb->type].name));
350 pr_int(fp, indent, "#molecules", molb->nmol);
351 pr_int(fp, indent, "#posres_xA", molb->posres_xA.size());
352 if (!molb->posres_xA.empty())
354 pr_rvecs(fp, indent, "posres_xA", as_rvec_array(molb->posres_xA.data()), molb->posres_xA.size());
356 pr_int(fp, indent, "#posres_xB", molb->posres_xB.size());
357 if (!molb->posres_xB.empty())
359 pr_rvecs(fp, indent, "posres_xB", as_rvec_array(molb->posres_xB.data()), molb->posres_xB.size());
363 void pr_mtop(FILE *fp, int indent, const char *title, const gmx_mtop_t *mtop,
364 gmx_bool bShowNumbers, gmx_bool bShowParameters)
366 if (available(fp, mtop, indent, title))
368 indent = pr_title(fp, indent, title);
369 pr_indent(fp, indent);
370 fprintf(fp, "name=\"%s\"\n", *(mtop->name));
371 pr_int(fp, indent, "#atoms", mtop->natoms);
372 pr_int(fp, indent, "#molblock", mtop->molblock.size());
373 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
375 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
377 pr_str(fp, indent, "bIntermolecularInteractions",
378 gmx::boolToString(mtop->bIntermolecularInteractions));
379 if (mtop->bIntermolecularInteractions)
381 for (int j = 0; j < F_NRE; j++)
383 pr_ilist(fp, indent, interaction_function[j].longname,
384 mtop->ffparams.functype.data(),
385 (*mtop->intermolecular_ilist)[j],
386 bShowNumbers, bShowParameters, mtop->ffparams.iparams.data());
389 pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
390 pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
391 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
393 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
394 &mtop->ffparams, bShowNumbers, bShowParameters);
396 pr_groups(fp, indent, mtop->groups, bShowNumbers);
400 void pr_top(FILE *fp, int indent, const char *title, const t_topology *top,
401 gmx_bool bShowNumbers, gmx_bool bShowParameters)
403 if (available(fp, top, indent, title))
405 indent = pr_title(fp, indent, title);
406 pr_indent(fp, indent);
407 fprintf(fp, "name=\"%s\"\n", *(top->name));
408 pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
409 pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
410 pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
411 pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
412 pr_str(fp, indent, "bIntermolecularInteractions",
413 gmx::boolToString(top->bIntermolecularInteractions));
414 pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
415 pr_idef(fp, indent, "idef", &top->idef, bShowNumbers, bShowParameters);
419 static void cmp_iparm(FILE *fp, const char *s, t_functype ft,
420 const t_iparams &ip1, const t_iparams &ip2, real relativeTolerance, real absoluteTolerance)
426 for (i = 0; i < MAXFORCEPARAM && !bDiff; i++)
428 bDiff = !equal_real(ip1.generic.buf[i], ip2.generic.buf[i], relativeTolerance, absoluteTolerance);
432 fprintf(fp, "%s1: ", s);
433 pr_iparams(fp, ft, &ip1);
434 fprintf(fp, "%s2: ", s);
435 pr_iparams(fp, ft, &ip2);
439 static void cmp_iparm_AB(FILE *fp, const char *s, t_functype ft,
440 const t_iparams &ip1, real relativeTolerance, real absoluteTolerance)
442 int nrfpA, nrfpB, p0, i;
445 /* Normally the first parameter is perturbable */
447 nrfpA = interaction_function[ft].nrfpA;
448 nrfpB = interaction_function[ft].nrfpB;
453 else if (interaction_function[ft].flags & IF_TABULATED)
455 /* For tabulated interactions only the second parameter is perturbable */
460 for (i = 0; i < nrfpB && !bDiff; i++)
462 bDiff = !equal_real(ip1.generic.buf[p0+i], ip1.generic.buf[nrfpA+i], relativeTolerance, absoluteTolerance);
466 fprintf(fp, "%s: ", s);
467 pr_iparams(fp, ft, &ip1);
471 static void cmp_cmap(FILE *fp, const gmx_cmap_t *cmap1, const gmx_cmap_t *cmap2, real relativeTolerance, real absoluteTolerance)
473 int cmap1_ngrid = (cmap1 ? cmap1->cmapdata.size() : 0);
474 int cmap2_ngrid = (cmap2 ? cmap2->cmapdata.size() : 0);
476 cmp_int(fp, "cmap ngrid", -1, cmap1_ngrid, cmap2_ngrid);
478 if (cmap1 == nullptr || cmap2 == nullptr)
483 cmp_int(fp, "cmap grid_spacing", -1, cmap1->grid_spacing, cmap2->grid_spacing);
484 if (cmap1->cmapdata.size() == cmap2->cmapdata.size() &&
485 cmap1->grid_spacing == cmap2->grid_spacing)
487 for (size_t g = 0; g < cmap1->cmapdata.size(); g++)
491 fprintf(fp, "comparing cmap %zu\n", g);
493 for (i = 0; i < 4*cmap1->grid_spacing*cmap1->grid_spacing; i++)
495 cmp_real(fp, "", i, cmap1->cmapdata[g].cmap[i], cmap2->cmapdata[g].cmap[i], relativeTolerance, absoluteTolerance);
501 static void cmp_block(FILE *fp, const t_block *b1, const t_block *b2, const char *s)
505 fprintf(fp, "comparing block %s\n", s);
506 sprintf(buf, "%s.nr", s);
507 cmp_int(fp, buf, -1, b1->nr, b2->nr);
510 static void cmp_blocka(FILE *fp, const t_blocka *b1, const t_blocka *b2, const char *s)
514 fprintf(fp, "comparing blocka %s\n", s);
515 sprintf(buf, "%s.nr", s);
516 cmp_int(fp, buf, -1, b1->nr, b2->nr);
517 sprintf(buf, "%s.nra", s);
518 cmp_int(fp, buf, -1, b1->nra, b2->nra);
521 static void compareFfparams(FILE *fp, const gmx_ffparams_t &ff1, const gmx_ffparams_t &ff2, real relativeTolerance, real absoluteTolerance)
523 fprintf(fp, "comparing force field parameters\n");
524 cmp_int(fp, "numTypes", -1, ff1.numTypes(), ff2.numTypes());
525 cmp_int(fp, "atnr", -1, ff1.atnr, ff1.atnr);
526 cmp_double(fp, "reppow", -1, ff1.reppow, ff2.reppow, relativeTolerance, absoluteTolerance);
527 cmp_real(fp, "fudgeQQ", -1, ff1.fudgeQQ, ff2.fudgeQQ, relativeTolerance, absoluteTolerance);
528 cmp_cmap(fp, &ff1.cmap_grid, &ff2.cmap_grid, relativeTolerance, absoluteTolerance);
529 for (int i = 0; i < std::min(ff1.numTypes(), ff2.numTypes()); i++)
531 std::string buf = gmx::formatString("ffparams->functype[%d]", i);
532 cmp_int(fp, buf.c_str(), i, ff1.functype[i], ff2.functype[i]);
533 buf = gmx::formatString("ffparams->iparams[%d]", i);
534 cmp_iparm(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], ff2.iparams[i], relativeTolerance, absoluteTolerance);
539 static void compareFfparamAB(FILE *fp, const gmx_ffparams_t &ff1, real relativeTolerance, real absoluteTolerance)
541 fprintf(fp, "comparing free energy parameters\n");
542 for (int i = 0; i < ff1.numTypes(); i++)
544 std::string buf = gmx::formatString("ffparams->iparams[%d]", i);
545 cmp_iparm_AB(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], relativeTolerance, absoluteTolerance);
548 static void compareInteractionLists(FILE *fp, const InteractionLists *il1, const InteractionLists *il2)
550 fprintf(fp, "comparing InteractionLists\n");
551 if ((il1 || il2) && (!il1 || !il2))
553 fprintf(fp, "InteractionLists are present in topology %d but not in the other\n", il1 ? 1 : 2);
557 for (int i = 0; i < F_NRE; i++)
559 cmp_int(fp, "InteractionList size", i, il1->at(i).size(), il2->at(i).size());
560 int nr = std::min(il1->at(i).size(), il2->at(i).size());
561 for (int j = 0; j < nr; j++)
563 cmp_int(fp, "InteractionList entry", j, il1->at(i).iatoms.at(j), il2->at(i).iatoms.at(j));
569 static void compareMoltypes(FILE *fp, gmx::ArrayRef<const gmx_moltype_t> mt1, gmx::ArrayRef<const gmx_moltype_t> mt2, real relativeTolerance, real absoluteTolerance)
571 fprintf(fp, "comparing molecule types\n");
572 cmp_int(fp, "moltype size", -1, mt1.size(), mt2.size());
573 for (int i = 0; i < std::min(mt1.ssize(), mt2.ssize()); i++)
575 cmp_str(fp, "Name", i, *mt1[i].name, *mt2[i].name);
576 compareAtoms(fp, &mt1[i].atoms, &mt2[i].atoms, relativeTolerance, absoluteTolerance);
577 compareInteractionLists(fp, &mt1[i].ilist, &mt2[i].ilist);
578 std::string buf = gmx::formatString("cgs[%d]", i);
579 cmp_block(fp, &mt1[i].cgs, &mt2[i].cgs, buf.c_str());
580 buf = gmx::formatString("excls[%d]", i);
581 cmp_blocka(fp, &mt1[i].excls, &mt2[i].excls, buf.c_str());
585 static void compareMoletypeAB(FILE *fp, gmx::ArrayRef<const gmx_moltype_t> mt1, real relativeTolerance, real absoluteTolerance)
587 fprintf(fp, "comparing free energy molecule types\n");
588 for (int i = 0; i < mt1.ssize(); i++)
590 compareAtoms(fp, &mt1[i].atoms, nullptr, relativeTolerance, absoluteTolerance);
593 static void compareMolblocks(FILE *fp, gmx::ArrayRef<const gmx_molblock_t> mb1, gmx::ArrayRef<const gmx_molblock_t> mb2)
595 fprintf(fp, "comparing molecule blocks\n");
596 cmp_int(fp, "molblock size", -1, mb1.size(), mb2.size());
597 int nr = std::min(mb1.size(), mb2.size());
598 for (int i = 0; i < nr; i++)
600 cmp_int(fp, "type", i, mb1[i].type, mb2[i].type);
601 cmp_int(fp, "nmol", i, mb1[i].nmol, mb2[i].nmol);
602 // Only checking size of restraint vectors for now
603 cmp_int(fp, "posres_xA size", i, mb1[i].posres_xA.size(), mb2[i].posres_xA.size());
604 cmp_int(fp, "posres_xB size", i, mb1[i].posres_xB.size(), mb2[i].posres_xB.size());
609 static void compareAtomtypes(FILE *fp, const t_atomtypes &at1, const t_atomtypes &at2)
611 fprintf(fp, "comparing atomtypes\n");
612 cmp_int(fp, "nr", -1, at1.nr, at2.nr);
613 int nr = std::min(at1.nr, at2.nr);
614 for (int i = 0; i < nr; i++)
616 cmp_int(fp, "atomtype", i, at1.atomnumber[i], at2.atomnumber[i]);
620 static void compareIntermolecularExclusions(FILE *fp, gmx::ArrayRef<const int> ime1, gmx::ArrayRef<const int> ime2)
622 fprintf(fp, "comparing intermolecular exclusions\n");
623 cmp_int(fp, "exclusion number", -1, ime1.size(), ime2.size());
624 int nr = std::min(ime1.size(), ime2.size());
625 for (int i = 0; i < nr; i++)
627 cmp_int(fp, "exclusion", i, ime1[i], ime2[i]);
631 static void compareBlockIndices(FILE *fp, gmx::ArrayRef<const MoleculeBlockIndices> mbi1, gmx::ArrayRef<const MoleculeBlockIndices> mbi2)
633 fprintf(fp, "comparing moleculeBlockIndices\n");
634 cmp_int(fp, "size", -1, mbi1.size(), mbi2.size());
635 int nr = std::min(mbi1.size(), mbi2.size());
636 for (int i = 0; i < nr; i++)
638 cmp_int(fp, "numAtomsPerMolecule", i, mbi1[i].numAtomsPerMolecule, mbi2[i].numAtomsPerMolecule);
639 cmp_int(fp, "globalAtomStart", i, mbi1[i].globalAtomStart, mbi2[i].globalAtomStart);
640 cmp_int(fp, "globalAtomEnd", i, mbi1[i].globalAtomEnd, mbi2[i].globalAtomEnd);
641 cmp_int(fp, "globalResidueStart", i, mbi1[i].globalResidueStart, mbi2[i].globalResidueStart);
642 cmp_int(fp, "moleculeIndexStart", i, mbi1[i].moleculeIndexStart, mbi2[i].moleculeIndexStart);
646 void compareMtop(FILE *fp, const gmx_mtop_t &mtop1, const gmx_mtop_t &mtop2, real relativeTolerance, real absoluteTolerance)
648 fprintf(fp, "comparing mtop topology\n");
649 cmp_str(fp, "Name", -1, *mtop1.name, *mtop2.name);
650 cmp_int(fp, "natoms", -1, mtop1.natoms, mtop2.natoms);
651 cmp_int(fp, "maxres_renum", -1, mtop1.maxres_renum, mtop2.maxres_renum);
652 cmp_int(fp, "maxresnr", -1, mtop1.maxresnr, mtop2.maxresnr);
653 cmp_bool(fp, "bIntermolecularInteractions", -1, mtop1.bIntermolecularInteractions, mtop2.bIntermolecularInteractions);
654 cmp_bool(fp, "haveMoleculeIndices", -1, mtop1.haveMoleculeIndices, mtop2.haveMoleculeIndices);
656 compareFfparams(fp, mtop1.ffparams, mtop2.ffparams, relativeTolerance, absoluteTolerance);
657 compareMoltypes(fp, mtop1.moltype, mtop2.moltype, relativeTolerance, absoluteTolerance);
658 compareMolblocks(fp, mtop1.molblock, mtop2.molblock);
659 compareInteractionLists(fp, mtop1.intermolecular_ilist.get(), mtop2.intermolecular_ilist.get());
660 compareAtomtypes(fp, mtop1.atomtypes, mtop2.atomtypes);
661 compareAtomGroups(fp, mtop1.groups, mtop2.groups, mtop1.natoms, mtop2.natoms);
662 compareIntermolecularExclusions(fp, mtop1.intermolecularExclusionGroup, mtop2.intermolecularExclusionGroup);
663 compareBlockIndices(fp, mtop1.moleculeBlockIndices, mtop2.moleculeBlockIndices);
666 void compareMtopAB(FILE *fp, const gmx_mtop_t &mtop1, real relativeTolerance, real absoluteTolerance)
668 fprintf(fp, "comparing topAB\n");
669 compareFfparamAB(fp, mtop1.ffparams, relativeTolerance, absoluteTolerance);
670 compareMoletypeAB(fp, mtop1.moltype, relativeTolerance, absoluteTolerance);
673 void compareAtomGroups(FILE *fp, const SimulationGroups &g0, const SimulationGroups &g1,
674 int natoms0, int natoms1)
676 fprintf(fp, "comparing groups\n");
678 for (auto group : keysOf(g0.groups))
680 std::string buf = gmx::formatString("grps[%d].nr", static_cast<int>(group));
681 cmp_int(fp, buf.c_str(), -1, g0.groups[group].nr, g1.groups[group].nr);
682 if (g0.groups[group].nr == g1.groups[group].nr)
684 for (int j = 0; j < g0.groups[group].nr; j++)
686 buf = gmx::formatString("grps[%d].name[%d]", static_cast<int>(group), j);
687 cmp_str(fp, buf.c_str(), -1,
688 *g0.groupNames[g0.groups[group].nm_ind[j]],
689 *g1.groupNames[g1.groups[group].nm_ind[j]]);
692 cmp_int(fp, "ngrpnr", static_cast<int>(group), g0.numberOfGroupNumbers(group), g1.numberOfGroupNumbers(group));
693 if (g0.numberOfGroupNumbers(group) == g1.numberOfGroupNumbers(group) && natoms0 == natoms1 &&
694 (!g0.groupNumbers[group].empty() || !g1.groupNumbers[group].empty()))
696 for (int j = 0; j < natoms0; j++)
698 cmp_int(fp, c_simulationAtomGroupTypeShortNames[group], j, getGroupType(g0, group, j), getGroupType(g1, group, j));
702 /* We have compared the names in the groups lists,
703 * so we can skip the grpname list comparison.
707 int getGroupType(const SimulationGroups &group, SimulationAtomGroupType type, int atom)
709 return (group.groupNumbers[type].empty() ? 0 : group.groupNumbers[type][atom]);
712 void copy_moltype(const gmx_moltype_t *src, gmx_moltype_t *dst)
714 dst->name = src->name;
715 copy_blocka(&src->excls, &dst->excls);
716 copy_block(&src->cgs, &dst->cgs);
717 t_atoms *atomsCopy = copy_t_atoms(&src->atoms);
718 dst->atoms = *atomsCopy;
721 for (int i = 0; i < F_NRE; ++i)
723 dst->ilist[i] = src->ilist[i];