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*> c_simulationAtomGroupTypeShortNames = {
58 { "T-Coupling", "Energy Mon.", "Acceleration", "Freeze", "User1", "User2", "VCM",
59 "Compressed X", "Or. Res. Fit", "QMMM" }
62 const char* shortName(SimulationAtomGroupType type)
64 return c_simulationAtomGroupTypeShortNames[type];
67 void init_top(t_topology* top)
70 init_idef(&top->idef);
71 init_atom(&(top->atoms));
72 init_atomtypes(&(top->atomtypes));
73 init_block(&top->mols);
74 init_blocka(&top->excls);
75 open_symtab(&top->symtab);
79 gmx_moltype_t::gmx_moltype_t() : name(nullptr)
81 init_t_atoms(&atoms, 0, FALSE);
84 gmx_moltype_t::~gmx_moltype_t()
89 gmx_mtop_t::gmx_mtop_t()
91 init_atomtypes(&atomtypes);
95 gmx_mtop_t::~gmx_mtop_t()
101 done_atomtypes(&atomtypes);
104 void done_top(t_topology* top)
106 done_idef(&top->idef);
107 done_atom(&(top->atoms));
110 done_atomtypes(&(top->atomtypes));
112 done_symtab(&(top->symtab));
113 done_block(&(top->mols));
114 done_blocka(&(top->excls));
117 void done_top_mtop(t_topology* top, gmx_mtop_t* mtop)
123 done_idef(&top->idef);
124 done_atom(&top->atoms);
125 done_blocka(&top->excls);
126 done_block(&top->mols);
127 done_symtab(&top->symtab);
128 open_symtab(&mtop->symtab);
129 done_atomtypes(&(top->atomtypes));
132 // Note that the rest of mtop will be freed by the destructor
136 gmx_localtop_t::gmx_localtop_t()
139 init_atomtypes(&atomtypes);
142 gmx_localtop_t::~gmx_localtop_t()
144 if (!useInDomainDecomp_)
147 done_atomtypes(&atomtypes);
151 bool gmx_mtop_has_masses(const gmx_mtop_t* mtop)
157 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveMass;
160 bool gmx_mtop_has_charges(const gmx_mtop_t* mtop)
166 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveCharge;
169 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t& mtop)
171 for (const gmx_moltype_t& moltype : mtop.moltype)
173 const t_atoms& atoms = moltype.atoms;
174 if (atoms.haveBState)
176 for (int a = 0; a < atoms.nr; a++)
178 if (atoms.atom[a].q != atoms.atom[a].qB)
188 bool gmx_mtop_has_atomtypes(const gmx_mtop_t* mtop)
194 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveType;
197 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t* mtop)
203 return mtop->moltype.empty() || mtop->moltype[0].atoms.havePdbInfo;
206 static void pr_grps(FILE* fp, const char* title, gmx::ArrayRef<const AtomGroupIndices> grps, char*** grpname)
209 for (const auto& group : grps)
211 fprintf(fp, "%s[%-12s] nr=%zu, name=[", title, c_simulationAtomGroupTypeShortNames[index],
213 for (const auto& entry : group)
215 fprintf(fp, " %s", *(grpname[entry]));
222 static void pr_groups(FILE* fp, int indent, const SimulationGroups& groups, gmx_bool bShowNumbers)
224 pr_grps(fp, "grp", groups.groups, const_cast<char***>(groups.groupNames.data()));
225 pr_strings(fp, indent, "grpname", const_cast<char***>(groups.groupNames.data()),
226 groups.groupNames.size(), bShowNumbers);
228 pr_indent(fp, indent);
229 fprintf(fp, "groups ");
230 for (const auto& group : c_simulationAtomGroupTypeShortNames)
232 printf(" %5.5s", group);
236 pr_indent(fp, indent);
237 fprintf(fp, "allocated ");
239 for (auto group : keysOf(groups.groups))
241 printf(" %5d", groups.numberOfGroupNumbers(group));
242 nat_max = std::max(nat_max, groups.numberOfGroupNumbers(group));
248 pr_indent(fp, indent);
249 fprintf(fp, "groupnr[%5s] =", "*");
250 for (auto gmx_unused group : keysOf(groups.groups))
252 fprintf(fp, " %3d ", 0);
258 for (int i = 0; i < nat_max; i++)
260 pr_indent(fp, indent);
261 fprintf(fp, "groupnr[%5d] =", i);
262 for (auto group : keysOf(groups.groups))
265 !groups.groupNumbers[group].empty() ? groups.groupNumbers[group][i] : 0);
272 static void pr_moltype(FILE* fp,
275 const gmx_moltype_t* molt,
277 const gmx_ffparams_t* ffparams,
278 gmx_bool bShowNumbers,
279 gmx_bool bShowParameters)
283 indent = pr_title_n(fp, indent, title, n);
284 pr_indent(fp, indent);
285 fprintf(fp, "name=\"%s\"\n", *(molt->name));
286 pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
287 pr_listoflists(fp, indent, "excls", &molt->excls, bShowNumbers);
288 for (j = 0; (j < F_NRE); j++)
290 pr_ilist(fp, indent, interaction_function[j].longname, ffparams->functype.data(),
291 molt->ilist[j], bShowNumbers, bShowParameters, ffparams->iparams.data());
295 static void pr_molblock(FILE* fp,
298 const gmx_molblock_t* molb,
300 const std::vector<gmx_moltype_t>& molt)
302 indent = pr_title_n(fp, indent, title, n);
303 pr_indent(fp, indent);
304 fprintf(fp, "%-20s = %d \"%s\"\n", "moltype", molb->type, *(molt[molb->type].name));
305 pr_int(fp, indent, "#molecules", molb->nmol);
306 pr_int(fp, indent, "#posres_xA", molb->posres_xA.size());
307 if (!molb->posres_xA.empty())
309 pr_rvecs(fp, indent, "posres_xA", as_rvec_array(molb->posres_xA.data()), molb->posres_xA.size());
311 pr_int(fp, indent, "#posres_xB", molb->posres_xB.size());
312 if (!molb->posres_xB.empty())
314 pr_rvecs(fp, indent, "posres_xB", as_rvec_array(molb->posres_xB.data()), molb->posres_xB.size());
318 void pr_mtop(FILE* fp, int indent, const char* title, const gmx_mtop_t* mtop, gmx_bool bShowNumbers, gmx_bool bShowParameters)
320 if (available(fp, mtop, indent, title))
322 indent = pr_title(fp, indent, title);
323 pr_indent(fp, indent);
324 fprintf(fp, "name=\"%s\"\n", *(mtop->name));
325 pr_int(fp, indent, "#atoms", mtop->natoms);
326 pr_int(fp, indent, "#molblock", mtop->molblock.size());
327 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
329 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
331 pr_str(fp, indent, "bIntermolecularInteractions",
332 gmx::boolToString(mtop->bIntermolecularInteractions));
333 if (mtop->bIntermolecularInteractions)
335 for (int j = 0; j < F_NRE; j++)
337 pr_ilist(fp, indent, interaction_function[j].longname,
338 mtop->ffparams.functype.data(), (*mtop->intermolecular_ilist)[j],
339 bShowNumbers, bShowParameters, mtop->ffparams.iparams.data());
342 pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
343 pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
344 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
346 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt, &mtop->ffparams, bShowNumbers,
349 pr_groups(fp, indent, mtop->groups, bShowNumbers);
353 void pr_top(FILE* fp, int indent, const char* title, const t_topology* top, gmx_bool bShowNumbers, gmx_bool bShowParameters)
355 if (available(fp, top, indent, title))
357 indent = pr_title(fp, indent, title);
358 pr_indent(fp, indent);
359 fprintf(fp, "name=\"%s\"\n", *(top->name));
360 pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
361 pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
362 pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
363 pr_str(fp, indent, "bIntermolecularInteractions",
364 gmx::boolToString(top->bIntermolecularInteractions));
365 pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
366 pr_idef(fp, indent, "idef", &top->idef, bShowNumbers, bShowParameters);
370 static void cmp_iparm(FILE* fp,
373 const t_iparams& ip1,
374 const t_iparams& ip2,
375 real relativeTolerance,
376 real absoluteTolerance)
382 for (i = 0; i < MAXFORCEPARAM && !bDiff; i++)
384 bDiff = !equal_real(ip1.generic.buf[i], ip2.generic.buf[i], relativeTolerance, absoluteTolerance);
388 fprintf(fp, "%s1: ", s);
389 pr_iparams(fp, ft, ip1);
390 fprintf(fp, "%s2: ", s);
391 pr_iparams(fp, ft, ip2);
395 static void cmp_iparm_AB(FILE* fp, const char* s, t_functype ft, const t_iparams& ip1, real relativeTolerance, real absoluteTolerance)
397 int nrfpA, nrfpB, p0, i;
400 /* Normally the first parameter is perturbable */
402 nrfpA = interaction_function[ft].nrfpA;
403 nrfpB = interaction_function[ft].nrfpB;
408 else if (interaction_function[ft].flags & IF_TABULATED)
410 /* For tabulated interactions only the second parameter is perturbable */
415 for (i = 0; i < nrfpB && !bDiff; i++)
417 bDiff = !equal_real(ip1.generic.buf[p0 + i], ip1.generic.buf[nrfpA + i], relativeTolerance,
422 fprintf(fp, "%s: ", s);
423 pr_iparams(fp, ft, ip1);
427 static void cmp_cmap(FILE* fp, const gmx_cmap_t* cmap1, const gmx_cmap_t* cmap2, real relativeTolerance, real absoluteTolerance)
429 int cmap1_ngrid = (cmap1 ? cmap1->cmapdata.size() : 0);
430 int cmap2_ngrid = (cmap2 ? cmap2->cmapdata.size() : 0);
432 cmp_int(fp, "cmap ngrid", -1, cmap1_ngrid, cmap2_ngrid);
434 if (cmap1 == nullptr || cmap2 == nullptr)
439 cmp_int(fp, "cmap grid_spacing", -1, cmap1->grid_spacing, cmap2->grid_spacing);
440 if (cmap1->cmapdata.size() == cmap2->cmapdata.size() && cmap1->grid_spacing == cmap2->grid_spacing)
442 for (size_t g = 0; g < cmap1->cmapdata.size(); g++)
446 fprintf(fp, "comparing cmap %zu\n", g);
448 for (i = 0; i < 4 * cmap1->grid_spacing * cmap1->grid_spacing; i++)
450 cmp_real(fp, "", i, cmap1->cmapdata[g].cmap[i], cmap2->cmapdata[g].cmap[i],
451 relativeTolerance, absoluteTolerance);
457 static void cmp_listoflists(FILE* fp,
458 const gmx::ListOfLists<int>& list1,
459 const gmx::ListOfLists<int>& list2,
464 fprintf(fp, "comparing blocka %s\n", s);
465 sprintf(buf, "%s.numLists", s);
466 cmp_int(fp, buf, -1, list1.ssize(), list2.ssize());
467 sprintf(buf, "%s.numElements", s);
468 cmp_int(fp, buf, -1, list1.numElements(), list2.numElements());
471 static void compareFfparams(FILE* fp,
472 const gmx_ffparams_t& ff1,
473 const gmx_ffparams_t& ff2,
474 real relativeTolerance,
475 real absoluteTolerance)
477 fprintf(fp, "comparing force field parameters\n");
478 cmp_int(fp, "numTypes", -1, ff1.numTypes(), ff2.numTypes());
479 cmp_int(fp, "atnr", -1, ff1.atnr, ff1.atnr);
480 cmp_double(fp, "reppow", -1, ff1.reppow, ff2.reppow, relativeTolerance, absoluteTolerance);
481 cmp_real(fp, "fudgeQQ", -1, ff1.fudgeQQ, ff2.fudgeQQ, relativeTolerance, absoluteTolerance);
482 cmp_cmap(fp, &ff1.cmap_grid, &ff2.cmap_grid, relativeTolerance, absoluteTolerance);
483 for (int i = 0; i < std::min(ff1.numTypes(), ff2.numTypes()); i++)
485 std::string buf = gmx::formatString("ffparams->functype[%d]", i);
486 cmp_int(fp, buf.c_str(), i, ff1.functype[i], ff2.functype[i]);
487 buf = gmx::formatString("ffparams->iparams[%d]", i);
488 cmp_iparm(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], ff2.iparams[i],
489 relativeTolerance, absoluteTolerance);
493 static void compareFfparamAB(FILE* fp, const gmx_ffparams_t& ff1, real relativeTolerance, real absoluteTolerance)
495 fprintf(fp, "comparing free energy parameters\n");
496 for (int i = 0; i < ff1.numTypes(); i++)
498 std::string buf = gmx::formatString("ffparams->iparams[%d]", i);
499 cmp_iparm_AB(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], relativeTolerance, absoluteTolerance);
502 static void compareInteractionLists(FILE* fp, const InteractionLists* il1, const InteractionLists* il2)
504 fprintf(fp, "comparing InteractionLists\n");
505 if ((il1 || il2) && (!il1 || !il2))
507 fprintf(fp, "InteractionLists are present in topology %d but not in the other\n", il1 ? 1 : 2);
511 for (int i = 0; i < F_NRE; i++)
513 cmp_int(fp, "InteractionList size", i, il1->at(i).size(), il2->at(i).size());
514 int nr = std::min(il1->at(i).size(), il2->at(i).size());
515 for (int j = 0; j < nr; j++)
517 cmp_int(fp, "InteractionList entry", j, il1->at(i).iatoms.at(j), il2->at(i).iatoms.at(j));
523 static void compareMoltypes(FILE* fp,
524 gmx::ArrayRef<const gmx_moltype_t> mt1,
525 gmx::ArrayRef<const gmx_moltype_t> mt2,
526 real relativeTolerance,
527 real absoluteTolerance)
529 fprintf(fp, "comparing molecule types\n");
530 cmp_int(fp, "moltype size", -1, mt1.size(), mt2.size());
531 for (int i = 0; i < std::min(mt1.ssize(), mt2.ssize()); i++)
533 cmp_str(fp, "Name", i, *mt1[i].name, *mt2[i].name);
534 compareAtoms(fp, &mt1[i].atoms, &mt2[i].atoms, relativeTolerance, absoluteTolerance);
535 compareInteractionLists(fp, &mt1[i].ilist, &mt2[i].ilist);
536 std::string buf = gmx::formatString("excls[%d]", i);
537 cmp_listoflists(fp, mt1[i].excls, mt2[i].excls, buf.c_str());
541 static void compareMoletypeAB(FILE* fp, gmx::ArrayRef<const gmx_moltype_t> mt1, real relativeTolerance, real absoluteTolerance)
543 fprintf(fp, "comparing free energy molecule types\n");
544 for (gmx::index i = 0; i < mt1.ssize(); i++)
546 compareAtoms(fp, &mt1[i].atoms, nullptr, relativeTolerance, absoluteTolerance);
549 static void compareMolblocks(FILE* fp,
550 gmx::ArrayRef<const gmx_molblock_t> mb1,
551 gmx::ArrayRef<const gmx_molblock_t> mb2)
553 fprintf(fp, "comparing molecule blocks\n");
554 cmp_int(fp, "molblock size", -1, mb1.size(), mb2.size());
555 int nr = std::min(mb1.size(), mb2.size());
556 for (int i = 0; i < nr; i++)
558 cmp_int(fp, "type", i, mb1[i].type, mb2[i].type);
559 cmp_int(fp, "nmol", i, mb1[i].nmol, mb2[i].nmol);
560 // Only checking size of restraint vectors for now
561 cmp_int(fp, "posres_xA size", i, mb1[i].posres_xA.size(), mb2[i].posres_xA.size());
562 cmp_int(fp, "posres_xB size", i, mb1[i].posres_xB.size(), mb2[i].posres_xB.size());
566 static void compareAtomtypes(FILE* fp, const t_atomtypes& at1, const t_atomtypes& at2)
568 fprintf(fp, "comparing atomtypes\n");
569 cmp_int(fp, "nr", -1, at1.nr, at2.nr);
570 int nr = std::min(at1.nr, at2.nr);
571 for (int i = 0; i < nr; i++)
573 cmp_int(fp, "atomtype", i, at1.atomnumber[i], at2.atomnumber[i]);
577 static void compareIntermolecularExclusions(FILE* fp,
578 gmx::ArrayRef<const int> ime1,
579 gmx::ArrayRef<const int> ime2)
581 fprintf(fp, "comparing intermolecular exclusions\n");
582 cmp_int(fp, "exclusion number", -1, ime1.size(), ime2.size());
583 int nr = std::min(ime1.size(), ime2.size());
584 for (int i = 0; i < nr; i++)
586 cmp_int(fp, "exclusion", i, ime1[i], ime2[i]);
590 static void compareBlockIndices(FILE* fp,
591 gmx::ArrayRef<const MoleculeBlockIndices> mbi1,
592 gmx::ArrayRef<const MoleculeBlockIndices> mbi2)
594 fprintf(fp, "comparing moleculeBlockIndices\n");
595 cmp_int(fp, "size", -1, mbi1.size(), mbi2.size());
596 int nr = std::min(mbi1.size(), mbi2.size());
597 for (int i = 0; i < nr; i++)
599 cmp_int(fp, "numAtomsPerMolecule", i, mbi1[i].numAtomsPerMolecule, mbi2[i].numAtomsPerMolecule);
600 cmp_int(fp, "globalAtomStart", i, mbi1[i].globalAtomStart, mbi2[i].globalAtomStart);
601 cmp_int(fp, "globalAtomEnd", i, mbi1[i].globalAtomEnd, mbi2[i].globalAtomEnd);
602 cmp_int(fp, "globalResidueStart", i, mbi1[i].globalResidueStart, mbi2[i].globalResidueStart);
603 cmp_int(fp, "moleculeIndexStart", i, mbi1[i].moleculeIndexStart, mbi2[i].moleculeIndexStart);
607 void compareMtop(FILE* fp, const gmx_mtop_t& mtop1, const gmx_mtop_t& mtop2, real relativeTolerance, real absoluteTolerance)
609 fprintf(fp, "comparing mtop topology\n");
610 cmp_str(fp, "Name", -1, *mtop1.name, *mtop2.name);
611 cmp_int(fp, "natoms", -1, mtop1.natoms, mtop2.natoms);
612 cmp_int(fp, "maxres_renum", -1, mtop1.maxres_renum, mtop2.maxres_renum);
613 cmp_int(fp, "maxresnr", -1, mtop1.maxresnr, mtop2.maxresnr);
614 cmp_bool(fp, "bIntermolecularInteractions", -1, mtop1.bIntermolecularInteractions,
615 mtop2.bIntermolecularInteractions);
616 cmp_bool(fp, "haveMoleculeIndices", -1, mtop1.haveMoleculeIndices, mtop2.haveMoleculeIndices);
618 compareFfparams(fp, mtop1.ffparams, mtop2.ffparams, relativeTolerance, absoluteTolerance);
619 compareMoltypes(fp, mtop1.moltype, mtop2.moltype, relativeTolerance, absoluteTolerance);
620 compareMolblocks(fp, mtop1.molblock, mtop2.molblock);
621 compareInteractionLists(fp, mtop1.intermolecular_ilist.get(), mtop2.intermolecular_ilist.get());
622 compareAtomtypes(fp, mtop1.atomtypes, mtop2.atomtypes);
623 compareAtomGroups(fp, mtop1.groups, mtop2.groups, mtop1.natoms, mtop2.natoms);
624 compareIntermolecularExclusions(fp, mtop1.intermolecularExclusionGroup,
625 mtop2.intermolecularExclusionGroup);
626 compareBlockIndices(fp, mtop1.moleculeBlockIndices, mtop2.moleculeBlockIndices);
629 void compareMtopAB(FILE* fp, const gmx_mtop_t& mtop1, real relativeTolerance, real absoluteTolerance)
631 fprintf(fp, "comparing topAB\n");
632 compareFfparamAB(fp, mtop1.ffparams, relativeTolerance, absoluteTolerance);
633 compareMoletypeAB(fp, mtop1.moltype, relativeTolerance, absoluteTolerance);
636 void compareAtomGroups(FILE* fp, const SimulationGroups& g0, const SimulationGroups& g1, int natoms0, int natoms1)
638 fprintf(fp, "comparing groups\n");
640 for (auto group : keysOf(g0.groups))
642 std::string buf = gmx::formatString("grps[%d].nr", static_cast<int>(group));
643 cmp_int(fp, buf.c_str(), -1, g0.groups[group].size(), g1.groups[group].size());
644 if (g0.groups[group].size() == g1.groups[group].size())
646 for (gmx::index j = 0; j < gmx::ssize(g0.groups[group]); j++)
648 buf = gmx::formatString("grps[%d].name[%zd]", static_cast<int>(group), j);
649 cmp_str(fp, buf.c_str(), -1, *g0.groupNames[g0.groups[group][j]],
650 *g1.groupNames[g1.groups[group][j]]);
653 cmp_int(fp, "ngrpnr", static_cast<int>(group), g0.numberOfGroupNumbers(group),
654 g1.numberOfGroupNumbers(group));
655 if (g0.numberOfGroupNumbers(group) == g1.numberOfGroupNumbers(group) && natoms0 == natoms1
656 && (!g0.groupNumbers[group].empty() || !g1.groupNumbers[group].empty()))
658 for (int j = 0; j < natoms0; j++)
660 cmp_int(fp, c_simulationAtomGroupTypeShortNames[group], j,
661 getGroupType(g0, group, j), getGroupType(g1, group, j));
665 /* We have compared the names in the groups lists,
666 * so we can skip the grpname list comparison.
670 int getGroupType(const SimulationGroups& group, SimulationAtomGroupType type, int atom)
672 return (group.groupNumbers[type].empty() ? 0 : group.groupNumbers[type][atom]);
675 void copy_moltype(const gmx_moltype_t* src, gmx_moltype_t* dst)
677 dst->name = src->name;
678 dst->excls = src->excls;
679 t_atoms* atomsCopy = copy_t_atoms(&src->atoms);
680 dst->atoms = *atomsCopy;
683 for (int i = 0; i < F_NRE; ++i)
685 dst->ilist[i] = src->ilist[i];