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 void init_top(t_topology *top)
80 init_idef(&top->idef);
81 init_atom(&(top->atoms));
82 init_atomtypes(&(top->atomtypes));
83 init_block(&top->mols);
84 init_blocka(&top->excls);
85 open_symtab(&top->symtab);
89 gmx_moltype_t::gmx_moltype_t() :
93 init_t_atoms(&atoms, 0, FALSE);
96 gmx_moltype_t::~gmx_moltype_t()
102 gmx_mtop_t::gmx_mtop_t()
104 init_atomtypes(&atomtypes);
105 open_symtab(&symtab);
108 gmx_mtop_t::~gmx_mtop_t()
110 done_symtab(&symtab);
114 done_atomtypes(&atomtypes);
117 void done_top(t_topology *top)
119 done_idef(&top->idef);
120 done_atom(&(top->atoms));
123 done_atomtypes(&(top->atomtypes));
125 done_symtab(&(top->symtab));
126 done_block(&(top->mols));
127 done_blocka(&(top->excls));
130 void done_top_mtop(t_topology *top, gmx_mtop_t *mtop)
136 done_idef(&top->idef);
137 done_atom(&top->atoms);
138 done_blocka(&top->excls);
139 done_block(&top->mols);
140 done_symtab(&top->symtab);
141 open_symtab(&mtop->symtab);
142 done_atomtypes(&(top->atomtypes));
145 // Note that the rest of mtop will be freed by the destructor
149 gmx_localtop_t::gmx_localtop_t()
151 init_blocka_null(&excls);
153 init_atomtypes(&atomtypes);
156 gmx_localtop_t::~gmx_localtop_t()
158 if (!useInDomainDecomp_)
162 done_atomtypes(&atomtypes);
166 bool gmx_mtop_has_masses(const gmx_mtop_t *mtop)
172 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveMass;
175 bool gmx_mtop_has_charges(const gmx_mtop_t *mtop)
181 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveCharge;
184 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t &mtop)
186 for (const gmx_moltype_t &moltype : mtop.moltype)
188 const t_atoms &atoms = moltype.atoms;
189 if (atoms.haveBState)
191 for (int a = 0; a < atoms.nr; a++)
193 if (atoms.atom[a].q != atoms.atom[a].qB)
203 bool gmx_mtop_has_atomtypes(const gmx_mtop_t *mtop)
209 return mtop->moltype.empty() || mtop->moltype[0].atoms.haveType;
212 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t *mtop)
218 return mtop->moltype.empty() || mtop->moltype[0].atoms.havePdbInfo;
221 static void pr_grps(FILE *fp,
223 gmx::ArrayRef<const AtomGroupIndices> grps,
227 for (const auto &group : grps)
229 fprintf(fp, "%s[%-12s] nr=%zu, name=[", title, c_simulationAtomGroupTypeShortNames[index], group.size());
230 for (const auto &entry : group)
232 fprintf(fp, " %s", *(grpname[entry]));
239 static void pr_groups(FILE *fp, int indent,
240 const SimulationGroups &groups,
241 gmx_bool bShowNumbers)
246 const_cast<char ***>(groups.groupNames.data()));
247 pr_strings(fp, indent, "grpname", const_cast<char ***>(groups.groupNames.data()), groups.groupNames.size(), bShowNumbers);
249 pr_indent(fp, indent);
250 fprintf(fp, "groups ");
251 for (const auto &group : c_simulationAtomGroupTypeShortNames)
253 printf(" %5.5s", group);
257 pr_indent(fp, indent);
258 fprintf(fp, "allocated ");
260 for (auto group : keysOf(groups.groups))
262 printf(" %5d", groups.numberOfGroupNumbers(group));
263 nat_max = std::max(nat_max, groups.numberOfGroupNumbers(group));
269 pr_indent(fp, indent);
270 fprintf(fp, "groupnr[%5s] =", "*");
271 for (auto gmx_unused group : keysOf(groups.groups))
273 fprintf(fp, " %3d ", 0);
279 for (int i = 0; i < nat_max; i++)
281 pr_indent(fp, indent);
282 fprintf(fp, "groupnr[%5d] =", i);
283 for (auto group : keysOf(groups.groups))
286 !groups.groupNumbers[group].empty() ?
287 groups.groupNumbers[group][i] : 0);
294 static void pr_moltype(FILE *fp, int indent, const char *title,
295 const gmx_moltype_t *molt, int n,
296 const gmx_ffparams_t *ffparams,
297 gmx_bool bShowNumbers, gmx_bool bShowParameters)
301 indent = pr_title_n(fp, indent, title, n);
302 pr_indent(fp, indent);
303 fprintf(fp, "name=\"%s\"\n", *(molt->name));
304 pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
305 pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
306 for (j = 0; (j < F_NRE); j++)
308 pr_ilist(fp, indent, interaction_function[j].longname,
309 ffparams->functype.data(), molt->ilist[j],
310 bShowNumbers, bShowParameters, ffparams->iparams.data());
314 static void pr_molblock(FILE *fp, int indent, const char *title,
315 const gmx_molblock_t *molb, int n,
316 const std::vector<gmx_moltype_t> &molt)
318 indent = pr_title_n(fp, indent, title, n);
319 pr_indent(fp, indent);
320 fprintf(fp, "%-20s = %d \"%s\"\n",
321 "moltype", molb->type, *(molt[molb->type].name));
322 pr_int(fp, indent, "#molecules", molb->nmol);
323 pr_int(fp, indent, "#posres_xA", molb->posres_xA.size());
324 if (!molb->posres_xA.empty())
326 pr_rvecs(fp, indent, "posres_xA", as_rvec_array(molb->posres_xA.data()), molb->posres_xA.size());
328 pr_int(fp, indent, "#posres_xB", molb->posres_xB.size());
329 if (!molb->posres_xB.empty())
331 pr_rvecs(fp, indent, "posres_xB", as_rvec_array(molb->posres_xB.data()), molb->posres_xB.size());
335 void pr_mtop(FILE *fp, int indent, const char *title, const gmx_mtop_t *mtop,
336 gmx_bool bShowNumbers, gmx_bool bShowParameters)
338 if (available(fp, mtop, indent, title))
340 indent = pr_title(fp, indent, title);
341 pr_indent(fp, indent);
342 fprintf(fp, "name=\"%s\"\n", *(mtop->name));
343 pr_int(fp, indent, "#atoms", mtop->natoms);
344 pr_int(fp, indent, "#molblock", mtop->molblock.size());
345 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
347 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
349 pr_str(fp, indent, "bIntermolecularInteractions",
350 gmx::boolToString(mtop->bIntermolecularInteractions));
351 if (mtop->bIntermolecularInteractions)
353 for (int j = 0; j < F_NRE; j++)
355 pr_ilist(fp, indent, interaction_function[j].longname,
356 mtop->ffparams.functype.data(),
357 (*mtop->intermolecular_ilist)[j],
358 bShowNumbers, bShowParameters, mtop->ffparams.iparams.data());
361 pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
362 pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
363 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
365 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
366 &mtop->ffparams, bShowNumbers, bShowParameters);
368 pr_groups(fp, indent, mtop->groups, bShowNumbers);
372 void pr_top(FILE *fp, int indent, const char *title, const t_topology *top,
373 gmx_bool bShowNumbers, gmx_bool bShowParameters)
375 if (available(fp, top, indent, title))
377 indent = pr_title(fp, indent, title);
378 pr_indent(fp, indent);
379 fprintf(fp, "name=\"%s\"\n", *(top->name));
380 pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
381 pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
382 pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
383 pr_str(fp, indent, "bIntermolecularInteractions",
384 gmx::boolToString(top->bIntermolecularInteractions));
385 pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
386 pr_idef(fp, indent, "idef", &top->idef, bShowNumbers, bShowParameters);
390 static void cmp_iparm(FILE *fp, const char *s, t_functype ft,
391 const t_iparams &ip1, const t_iparams &ip2, real relativeTolerance, real absoluteTolerance)
397 for (i = 0; i < MAXFORCEPARAM && !bDiff; i++)
399 bDiff = !equal_real(ip1.generic.buf[i], ip2.generic.buf[i], relativeTolerance, absoluteTolerance);
403 fprintf(fp, "%s1: ", s);
404 pr_iparams(fp, ft, &ip1);
405 fprintf(fp, "%s2: ", s);
406 pr_iparams(fp, ft, &ip2);
410 static void cmp_iparm_AB(FILE *fp, const char *s, t_functype ft,
411 const t_iparams &ip1, real relativeTolerance, real absoluteTolerance)
413 int nrfpA, nrfpB, p0, i;
416 /* Normally the first parameter is perturbable */
418 nrfpA = interaction_function[ft].nrfpA;
419 nrfpB = interaction_function[ft].nrfpB;
424 else if (interaction_function[ft].flags & IF_TABULATED)
426 /* For tabulated interactions only the second parameter is perturbable */
431 for (i = 0; i < nrfpB && !bDiff; i++)
433 bDiff = !equal_real(ip1.generic.buf[p0+i], ip1.generic.buf[nrfpA+i], relativeTolerance, absoluteTolerance);
437 fprintf(fp, "%s: ", s);
438 pr_iparams(fp, ft, &ip1);
442 static void cmp_cmap(FILE *fp, const gmx_cmap_t *cmap1, const gmx_cmap_t *cmap2, real relativeTolerance, real absoluteTolerance)
444 int cmap1_ngrid = (cmap1 ? cmap1->cmapdata.size() : 0);
445 int cmap2_ngrid = (cmap2 ? cmap2->cmapdata.size() : 0);
447 cmp_int(fp, "cmap ngrid", -1, cmap1_ngrid, cmap2_ngrid);
449 if (cmap1 == nullptr || cmap2 == nullptr)
454 cmp_int(fp, "cmap grid_spacing", -1, cmap1->grid_spacing, cmap2->grid_spacing);
455 if (cmap1->cmapdata.size() == cmap2->cmapdata.size() &&
456 cmap1->grid_spacing == cmap2->grid_spacing)
458 for (size_t g = 0; g < cmap1->cmapdata.size(); g++)
462 fprintf(fp, "comparing cmap %zu\n", g);
464 for (i = 0; i < 4*cmap1->grid_spacing*cmap1->grid_spacing; i++)
466 cmp_real(fp, "", i, cmap1->cmapdata[g].cmap[i], cmap2->cmapdata[g].cmap[i], relativeTolerance, absoluteTolerance);
472 static void cmp_blocka(FILE *fp, const t_blocka *b1, const t_blocka *b2, const char *s)
476 fprintf(fp, "comparing blocka %s\n", s);
477 sprintf(buf, "%s.nr", s);
478 cmp_int(fp, buf, -1, b1->nr, b2->nr);
479 sprintf(buf, "%s.nra", s);
480 cmp_int(fp, buf, -1, b1->nra, b2->nra);
483 static void compareFfparams(FILE *fp, const gmx_ffparams_t &ff1, const gmx_ffparams_t &ff2, real relativeTolerance, real absoluteTolerance)
485 fprintf(fp, "comparing force field parameters\n");
486 cmp_int(fp, "numTypes", -1, ff1.numTypes(), ff2.numTypes());
487 cmp_int(fp, "atnr", -1, ff1.atnr, ff1.atnr);
488 cmp_double(fp, "reppow", -1, ff1.reppow, ff2.reppow, relativeTolerance, absoluteTolerance);
489 cmp_real(fp, "fudgeQQ", -1, ff1.fudgeQQ, ff2.fudgeQQ, relativeTolerance, absoluteTolerance);
490 cmp_cmap(fp, &ff1.cmap_grid, &ff2.cmap_grid, relativeTolerance, absoluteTolerance);
491 for (int i = 0; i < std::min(ff1.numTypes(), ff2.numTypes()); i++)
493 std::string buf = gmx::formatString("ffparams->functype[%d]", i);
494 cmp_int(fp, buf.c_str(), i, ff1.functype[i], ff2.functype[i]);
495 buf = gmx::formatString("ffparams->iparams[%d]", i);
496 cmp_iparm(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], ff2.iparams[i], relativeTolerance, absoluteTolerance);
501 static void compareFfparamAB(FILE *fp, const gmx_ffparams_t &ff1, real relativeTolerance, real absoluteTolerance)
503 fprintf(fp, "comparing free energy parameters\n");
504 for (int i = 0; i < ff1.numTypes(); i++)
506 std::string buf = gmx::formatString("ffparams->iparams[%d]", i);
507 cmp_iparm_AB(fp, buf.c_str(), ff1.functype[i], ff1.iparams[i], relativeTolerance, absoluteTolerance);
510 static void compareInteractionLists(FILE *fp, const InteractionLists *il1, const InteractionLists *il2)
512 fprintf(fp, "comparing InteractionLists\n");
513 if ((il1 || il2) && (!il1 || !il2))
515 fprintf(fp, "InteractionLists are present in topology %d but not in the other\n", il1 ? 1 : 2);
519 for (int i = 0; i < F_NRE; i++)
521 cmp_int(fp, "InteractionList size", i, il1->at(i).size(), il2->at(i).size());
522 int nr = std::min(il1->at(i).size(), il2->at(i).size());
523 for (int j = 0; j < nr; j++)
525 cmp_int(fp, "InteractionList entry", j, il1->at(i).iatoms.at(j), il2->at(i).iatoms.at(j));
531 static void compareMoltypes(FILE *fp, gmx::ArrayRef<const gmx_moltype_t> mt1, gmx::ArrayRef<const gmx_moltype_t> mt2, real relativeTolerance, real absoluteTolerance)
533 fprintf(fp, "comparing molecule types\n");
534 cmp_int(fp, "moltype size", -1, mt1.size(), mt2.size());
535 for (int i = 0; i < std::min(mt1.ssize(), mt2.ssize()); i++)
537 cmp_str(fp, "Name", i, *mt1[i].name, *mt2[i].name);
538 compareAtoms(fp, &mt1[i].atoms, &mt2[i].atoms, relativeTolerance, absoluteTolerance);
539 compareInteractionLists(fp, &mt1[i].ilist, &mt2[i].ilist);
540 std::string buf = gmx::formatString("excls[%d]", i);
541 cmp_blocka(fp, &mt1[i].excls, &mt2[i].excls, buf.c_str());
545 static void compareMoletypeAB(FILE *fp, gmx::ArrayRef<const gmx_moltype_t> mt1, real relativeTolerance, real absoluteTolerance)
547 fprintf(fp, "comparing free energy molecule types\n");
548 for (gmx::index i = 0; i < mt1.ssize(); i++)
550 compareAtoms(fp, &mt1[i].atoms, nullptr, relativeTolerance, absoluteTolerance);
553 static void compareMolblocks(FILE *fp, gmx::ArrayRef<const gmx_molblock_t> mb1, gmx::ArrayRef<const gmx_molblock_t> mb2)
555 fprintf(fp, "comparing molecule blocks\n");
556 cmp_int(fp, "molblock size", -1, mb1.size(), mb2.size());
557 int nr = std::min(mb1.size(), mb2.size());
558 for (int i = 0; i < nr; i++)
560 cmp_int(fp, "type", i, mb1[i].type, mb2[i].type);
561 cmp_int(fp, "nmol", i, mb1[i].nmol, mb2[i].nmol);
562 // Only checking size of restraint vectors for now
563 cmp_int(fp, "posres_xA size", i, mb1[i].posres_xA.size(), mb2[i].posres_xA.size());
564 cmp_int(fp, "posres_xB size", i, mb1[i].posres_xB.size(), mb2[i].posres_xB.size());
569 static void compareAtomtypes(FILE *fp, const t_atomtypes &at1, const t_atomtypes &at2)
571 fprintf(fp, "comparing atomtypes\n");
572 cmp_int(fp, "nr", -1, at1.nr, at2.nr);
573 int nr = std::min(at1.nr, at2.nr);
574 for (int i = 0; i < nr; i++)
576 cmp_int(fp, "atomtype", i, at1.atomnumber[i], at2.atomnumber[i]);
580 static void compareIntermolecularExclusions(FILE *fp, gmx::ArrayRef<const int> ime1, gmx::ArrayRef<const int> ime2)
582 fprintf(fp, "comparing intermolecular exclusions\n");
583 cmp_int(fp, "exclusion number", -1, ime1.size(), ime2.size());
584 int nr = std::min(ime1.size(), ime2.size());
585 for (int i = 0; i < nr; i++)
587 cmp_int(fp, "exclusion", i, ime1[i], ime2[i]);
591 static void compareBlockIndices(FILE *fp, gmx::ArrayRef<const MoleculeBlockIndices> mbi1, gmx::ArrayRef<const MoleculeBlockIndices> mbi2)
593 fprintf(fp, "comparing moleculeBlockIndices\n");
594 cmp_int(fp, "size", -1, mbi1.size(), mbi2.size());
595 int nr = std::min(mbi1.size(), mbi2.size());
596 for (int i = 0; i < nr; i++)
598 cmp_int(fp, "numAtomsPerMolecule", i, mbi1[i].numAtomsPerMolecule, mbi2[i].numAtomsPerMolecule);
599 cmp_int(fp, "globalAtomStart", i, mbi1[i].globalAtomStart, mbi2[i].globalAtomStart);
600 cmp_int(fp, "globalAtomEnd", i, mbi1[i].globalAtomEnd, mbi2[i].globalAtomEnd);
601 cmp_int(fp, "globalResidueStart", i, mbi1[i].globalResidueStart, mbi2[i].globalResidueStart);
602 cmp_int(fp, "moleculeIndexStart", i, mbi1[i].moleculeIndexStart, mbi2[i].moleculeIndexStart);
606 void compareMtop(FILE *fp, const gmx_mtop_t &mtop1, const gmx_mtop_t &mtop2, real relativeTolerance, real absoluteTolerance)
608 fprintf(fp, "comparing mtop topology\n");
609 cmp_str(fp, "Name", -1, *mtop1.name, *mtop2.name);
610 cmp_int(fp, "natoms", -1, mtop1.natoms, mtop2.natoms);
611 cmp_int(fp, "maxres_renum", -1, mtop1.maxres_renum, mtop2.maxres_renum);
612 cmp_int(fp, "maxresnr", -1, mtop1.maxresnr, mtop2.maxresnr);
613 cmp_bool(fp, "bIntermolecularInteractions", -1, mtop1.bIntermolecularInteractions, mtop2.bIntermolecularInteractions);
614 cmp_bool(fp, "haveMoleculeIndices", -1, mtop1.haveMoleculeIndices, mtop2.haveMoleculeIndices);
616 compareFfparams(fp, mtop1.ffparams, mtop2.ffparams, relativeTolerance, absoluteTolerance);
617 compareMoltypes(fp, mtop1.moltype, mtop2.moltype, relativeTolerance, absoluteTolerance);
618 compareMolblocks(fp, mtop1.molblock, mtop2.molblock);
619 compareInteractionLists(fp, mtop1.intermolecular_ilist.get(), mtop2.intermolecular_ilist.get());
620 compareAtomtypes(fp, mtop1.atomtypes, mtop2.atomtypes);
621 compareAtomGroups(fp, mtop1.groups, mtop2.groups, mtop1.natoms, mtop2.natoms);
622 compareIntermolecularExclusions(fp, mtop1.intermolecularExclusionGroup, mtop2.intermolecularExclusionGroup);
623 compareBlockIndices(fp, mtop1.moleculeBlockIndices, mtop2.moleculeBlockIndices);
626 void compareMtopAB(FILE *fp, const gmx_mtop_t &mtop1, real relativeTolerance, real absoluteTolerance)
628 fprintf(fp, "comparing topAB\n");
629 compareFfparamAB(fp, mtop1.ffparams, relativeTolerance, absoluteTolerance);
630 compareMoletypeAB(fp, mtop1.moltype, relativeTolerance, absoluteTolerance);
633 void compareAtomGroups(FILE *fp, const SimulationGroups &g0, const SimulationGroups &g1,
634 int natoms0, int natoms1)
636 fprintf(fp, "comparing groups\n");
638 for (auto group : keysOf(g0.groups))
640 std::string buf = gmx::formatString("grps[%d].nr", static_cast<int>(group));
641 cmp_int(fp, buf.c_str(), -1, g0.groups[group].size(), g1.groups[group].size());
642 if (g0.groups[group].size() == g1.groups[group].size())
644 for (gmx::index j = 0; j < gmx::ssize(g0.groups[group]); j++)
646 buf = gmx::formatString("grps[%d].name[%zd]", static_cast<int>(group), j);
647 cmp_str(fp, buf.c_str(), -1,
648 *g0.groupNames[g0.groups[group][j]],
649 *g1.groupNames[g1.groups[group][j]]);
652 cmp_int(fp, "ngrpnr", static_cast<int>(group), g0.numberOfGroupNumbers(group), g1.numberOfGroupNumbers(group));
653 if (g0.numberOfGroupNumbers(group) == g1.numberOfGroupNumbers(group) && natoms0 == natoms1 &&
654 (!g0.groupNumbers[group].empty() || !g1.groupNumbers[group].empty()))
656 for (int j = 0; j < natoms0; j++)
658 cmp_int(fp, c_simulationAtomGroupTypeShortNames[group], j, getGroupType(g0, group, j), getGroupType(g1, group, j));
662 /* We have compared the names in the groups lists,
663 * so we can skip the grpname list comparison.
667 int getGroupType(const SimulationGroups &group, SimulationAtomGroupType type, int atom)
669 return (group.groupNumbers[type].empty() ? 0 : group.groupNumbers[type][atom]);
672 void copy_moltype(const gmx_moltype_t *src, gmx_moltype_t *dst)
674 dst->name = src->name;
675 copy_blocka(&src->excls, &dst->excls);
676 t_atoms *atomsCopy = copy_t_atoms(&src->atoms);
677 dst->atoms = *atomsCopy;
680 for (int i = 0; i < F_NRE; ++i)
682 dst->ilist[i] = src->ilist[i];