2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2008,2009,2010, The GROMACS development team.
5 * Copyright (c) 2012,2013,2014,2015,2016 The GROMACS development team.
6 * Copyright (c) 2017,2018,2019,2020, 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.
39 #include "mtop_util.h"
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/topology/atoms.h"
49 #include "gromacs/topology/block.h"
50 #include "gromacs/topology/exclusionblocks.h"
51 #include "gromacs/topology/idef.h"
52 #include "gromacs/topology/ifunc.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/topology/topsort.h"
55 #include "gromacs/utility/arrayref.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/real.h"
58 #include "gromacs/utility/smalloc.h"
60 static int gmx_mtop_maxresnr(const gmx_mtop_t* mtop, int maxres_renum)
64 for (const gmx_moltype_t& moltype : mtop->moltype)
66 const t_atoms& atoms = moltype.atoms;
67 if (atoms.nres > maxres_renum)
69 for (int r = 0; r < atoms.nres; r++)
71 if (atoms.resinfo[r].nr > maxresnr)
73 maxresnr = atoms.resinfo[r].nr;
82 static void buildMolblockIndices(gmx_mtop_t* mtop)
84 mtop->moleculeBlockIndices.resize(mtop->molblock.size());
88 int residueNumberStart = mtop->maxresnr + 1;
89 int moleculeIndexStart = 0;
90 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
92 const gmx_molblock_t& molb = mtop->molblock[mb];
93 MoleculeBlockIndices& indices = mtop->moleculeBlockIndices[mb];
94 const int numResPerMol = mtop->moltype[molb.type].atoms.nres;
96 indices.numAtomsPerMolecule = mtop->moltype[molb.type].atoms.nr;
97 indices.globalAtomStart = atomIndex;
98 indices.globalResidueStart = residueIndex;
99 atomIndex += molb.nmol * indices.numAtomsPerMolecule;
100 residueIndex += molb.nmol * numResPerMol;
101 indices.globalAtomEnd = atomIndex;
102 indices.residueNumberStart = residueNumberStart;
103 if (numResPerMol <= mtop->maxres_renum)
105 residueNumberStart += molb.nmol * numResPerMol;
107 indices.moleculeIndexStart = moleculeIndexStart;
108 moleculeIndexStart += molb.nmol;
112 void gmx_mtop_finalize(gmx_mtop_t* mtop)
116 if (mtop->molblock.size() == 1 && mtop->molblock[0].nmol == 1)
118 /* We have a single molecule only, no renumbering needed.
119 * This case also covers an mtop converted from pdb/gro/... input,
120 * so we retain the original residue numbering.
122 mtop->maxres_renum = 0;
126 /* We only renumber single residue molecules. Their intra-molecular
127 * residue numbering is anyhow irrelevant.
129 mtop->maxres_renum = 1;
132 env = getenv("GMX_MAXRESRENUM");
135 sscanf(env, "%d", &mtop->maxres_renum);
137 if (mtop->maxres_renum == -1)
139 /* -1 signals renumber residues in all molecules */
140 mtop->maxres_renum = INT_MAX;
143 mtop->maxresnr = gmx_mtop_maxresnr(mtop, mtop->maxres_renum);
145 buildMolblockIndices(mtop);
148 void gmx_mtop_count_atomtypes(const gmx_mtop_t* mtop, int state, int typecount[])
150 for (int i = 0; i < mtop->ffparams.atnr; ++i)
154 for (const gmx_molblock_t& molb : mtop->molblock)
156 const t_atoms& atoms = mtop->moltype[molb.type].atoms;
157 for (int i = 0; i < atoms.nr; ++i)
162 tpi = atoms.atom[i].type;
166 tpi = atoms.atom[i].typeB;
168 typecount[tpi] += molb.nmol;
173 int gmx_mtop_num_molecules(const gmx_mtop_t& mtop)
175 int numMolecules = 0;
176 for (const gmx_molblock_t& molb : mtop.molblock)
178 numMolecules += molb.nmol;
183 int gmx_mtop_nres(const gmx_mtop_t* mtop)
186 for (const gmx_molblock_t& molb : mtop->molblock)
188 nres += molb.nmol * mtop->moltype[molb.type].atoms.nres;
193 AtomIterator::AtomIterator(const gmx_mtop_t& mtop, int globalAtomNumber) :
196 atoms_(&mtop.moltype[mtop.molblock[0].type].atoms),
198 highestResidueNumber_(mtop.maxresnr),
200 globalAtomNumber_(globalAtomNumber)
202 GMX_ASSERT(globalAtomNumber == 0 || globalAtomNumber == mtop.natoms,
203 "Starting at other atoms not implemented yet");
206 AtomIterator& AtomIterator::operator++()
211 if (localAtomNumber_ >= atoms_->nr)
213 if (atoms_->nres <= mtop_->maxresnr)
215 /* Single residue molecule, increase the count with one */
216 highestResidueNumber_ += atoms_->nres;
219 localAtomNumber_ = 0;
220 if (currentMolecule_ >= mtop_->molblock[mblock_].nmol)
223 if (mblock_ >= mtop_->molblock.size())
227 atoms_ = &mtop_->moltype[mtop_->molblock[mblock_].type].atoms;
228 currentMolecule_ = 0;
234 AtomIterator AtomIterator::operator++(int)
236 AtomIterator temp = *this;
241 bool AtomIterator::operator==(const AtomIterator& o) const
243 return mtop_ == o.mtop_ && globalAtomNumber_ == o.globalAtomNumber_;
246 bool AtomIterator::operator!=(const AtomIterator& o) const
248 return !(*this == o);
251 const t_atom& AtomProxy::atom() const
253 return it_->atoms_->atom[it_->localAtomNumber_];
256 int AtomProxy::globalAtomNumber() const
258 return it_->globalAtomNumber_;
261 const char* AtomProxy::atomName() const
263 return *(it_->atoms_->atomname[it_->localAtomNumber_]);
266 const char* AtomProxy::residueName() const
268 int residueIndexInMolecule = it_->atoms_->atom[it_->localAtomNumber_].resind;
269 return *(it_->atoms_->resinfo[residueIndexInMolecule].name);
272 int AtomProxy::residueNumber() const
274 int residueIndexInMolecule = it_->atoms_->atom[it_->localAtomNumber_].resind;
275 if (it_->atoms_->nres <= it_->mtop_->maxres_renum)
277 return it_->highestResidueNumber_ + 1 + residueIndexInMolecule;
281 return it_->atoms_->resinfo[residueIndexInMolecule].nr;
285 const gmx_moltype_t& AtomProxy::moleculeType() const
287 return it_->mtop_->moltype[it_->mtop_->molblock[it_->mblock_].type];
290 int AtomProxy::atomNumberInMol() const
292 return it_->localAtomNumber_;
295 typedef struct gmx_mtop_atomloop_block
297 const gmx_mtop_t* mtop;
299 const t_atoms* atoms;
301 } t_gmx_mtop_atomloop_block;
303 gmx_mtop_atomloop_block_t gmx_mtop_atomloop_block_init(const gmx_mtop_t* mtop)
305 struct gmx_mtop_atomloop_block* aloop;
311 aloop->atoms = &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
312 aloop->at_local = -1;
317 static void gmx_mtop_atomloop_block_destroy(gmx_mtop_atomloop_block_t aloop)
322 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop, const t_atom** atom, int* nmol)
324 if (aloop == nullptr)
326 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
331 if (aloop->at_local >= aloop->atoms->nr)
334 if (aloop->mblock >= aloop->mtop->molblock.size())
336 gmx_mtop_atomloop_block_destroy(aloop);
339 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
343 *atom = &aloop->atoms->atom[aloop->at_local];
344 *nmol = aloop->mtop->molblock[aloop->mblock].nmol;
349 typedef struct gmx_mtop_ilistloop
351 const gmx_mtop_t* mtop;
355 gmx_mtop_ilistloop_t gmx_mtop_ilistloop_init(const gmx_mtop_t* mtop)
357 struct gmx_mtop_ilistloop* iloop;
367 gmx_mtop_ilistloop_t gmx_mtop_ilistloop_init(const gmx_mtop_t& mtop)
369 return gmx_mtop_ilistloop_init(&mtop);
372 static void gmx_mtop_ilistloop_destroy(gmx_mtop_ilistloop_t iloop)
377 const InteractionLists* gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop, int* nmol)
379 if (iloop == nullptr)
381 gmx_incons("gmx_mtop_ilistloop_next called without calling gmx_mtop_ilistloop_init");
385 if (iloop->mblock >= gmx::ssize(iloop->mtop->molblock))
387 if (iloop->mblock == gmx::ssize(iloop->mtop->molblock) && iloop->mtop->bIntermolecularInteractions)
390 return iloop->mtop->intermolecular_ilist.get();
393 gmx_mtop_ilistloop_destroy(iloop);
397 *nmol = iloop->mtop->molblock[iloop->mblock].nmol;
399 return &iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
401 typedef struct gmx_mtop_ilistloop_all
403 const gmx_mtop_t* mtop;
407 } t_gmx_mtop_ilist_all;
409 int gmx_mtop_ftype_count(const gmx_mtop_t* mtop, int ftype)
411 gmx_mtop_ilistloop_t iloop;
416 iloop = gmx_mtop_ilistloop_init(mtop);
417 while (const InteractionLists* il = gmx_mtop_ilistloop_next(iloop, &nmol))
419 n += nmol * (*il)[ftype].size() / (1 + NRAL(ftype));
422 if (mtop->bIntermolecularInteractions)
424 n += (*mtop->intermolecular_ilist)[ftype].size() / (1 + NRAL(ftype));
430 int gmx_mtop_ftype_count(const gmx_mtop_t& mtop, int ftype)
432 return gmx_mtop_ftype_count(&mtop, ftype);
435 int gmx_mtop_interaction_count(const gmx_mtop_t& mtop, const int unsigned if_flags)
439 gmx_mtop_ilistloop_t iloop = gmx_mtop_ilistloop_init(mtop);
441 while (const InteractionLists* il = gmx_mtop_ilistloop_next(iloop, &nmol))
443 for (int ftype = 0; ftype < F_NRE; ftype++)
445 if ((interaction_function[ftype].flags & if_flags) == if_flags)
447 n += nmol * (*il)[ftype].size() / (1 + NRAL(ftype));
452 if (mtop.bIntermolecularInteractions)
454 for (int ftype = 0; ftype < F_NRE; ftype++)
456 if ((interaction_function[ftype].flags & if_flags) == if_flags)
458 n += (*mtop.intermolecular_ilist)[ftype].size() / (1 + NRAL(ftype));
466 std::array<int, eptNR> gmx_mtop_particletype_count(const gmx_mtop_t& mtop)
468 std::array<int, eptNR> count = { { 0 } };
470 for (const auto& molblock : mtop.molblock)
472 const t_atoms& atoms = mtop.moltype[molblock.type].atoms;
473 for (int a = 0; a < atoms.nr; a++)
475 count[atoms.atom[a].ptype] += molblock.nmol;
482 static void atomcat(t_atoms* dest, const t_atoms* src, int copies, int maxres_renum, int* maxresnr)
486 int destnr = dest->nr;
490 dest->haveMass = src->haveMass;
491 dest->haveType = src->haveType;
492 dest->haveCharge = src->haveCharge;
493 dest->haveBState = src->haveBState;
494 dest->havePdbInfo = src->havePdbInfo;
498 dest->haveMass = dest->haveMass && src->haveMass;
499 dest->haveType = dest->haveType && src->haveType;
500 dest->haveCharge = dest->haveCharge && src->haveCharge;
501 dest->haveBState = dest->haveBState && src->haveBState;
502 dest->havePdbInfo = dest->havePdbInfo && src->havePdbInfo;
507 size = destnr + copies * srcnr;
508 srenew(dest->atom, size);
509 srenew(dest->atomname, size);
512 srenew(dest->atomtype, size);
513 if (dest->haveBState)
515 srenew(dest->atomtypeB, size);
518 if (dest->havePdbInfo)
520 srenew(dest->pdbinfo, size);
525 size = dest->nres + copies * src->nres;
526 srenew(dest->resinfo, size);
529 /* residue information */
530 for (l = dest->nres, j = 0; (j < copies); j++, l += src->nres)
532 memcpy(reinterpret_cast<char*>(&(dest->resinfo[l])), reinterpret_cast<char*>(&(src->resinfo[0])),
533 static_cast<size_t>(src->nres * sizeof(src->resinfo[0])));
536 for (l = destnr, j = 0; (j < copies); j++, l += srcnr)
538 memcpy(reinterpret_cast<char*>(&(dest->atom[l])), reinterpret_cast<char*>(&(src->atom[0])),
539 static_cast<size_t>(srcnr * sizeof(src->atom[0])));
540 memcpy(reinterpret_cast<char*>(&(dest->atomname[l])),
541 reinterpret_cast<char*>(&(src->atomname[0])),
542 static_cast<size_t>(srcnr * sizeof(src->atomname[0])));
545 memcpy(reinterpret_cast<char*>(&(dest->atomtype[l])),
546 reinterpret_cast<char*>(&(src->atomtype[0])),
547 static_cast<size_t>(srcnr * sizeof(src->atomtype[0])));
548 if (dest->haveBState)
550 memcpy(reinterpret_cast<char*>(&(dest->atomtypeB[l])),
551 reinterpret_cast<char*>(&(src->atomtypeB[0])),
552 static_cast<size_t>(srcnr * sizeof(src->atomtypeB[0])));
555 if (dest->havePdbInfo)
557 memcpy(reinterpret_cast<char*>(&(dest->pdbinfo[l])),
558 reinterpret_cast<char*>(&(src->pdbinfo[0])),
559 static_cast<size_t>(srcnr * sizeof(src->pdbinfo[0])));
563 /* Increment residue indices */
564 for (l = destnr, j = 0; (j < copies); j++)
566 for (i = 0; (i < srcnr); i++, l++)
568 dest->atom[l].resind = dest->nres + j * src->nres + src->atom[i].resind;
572 if (src->nres <= maxres_renum)
574 /* Single residue molecule, continue counting residues */
575 for (j = 0; (j < copies); j++)
577 for (l = 0; l < src->nres; l++)
580 dest->resinfo[dest->nres + j * src->nres + l].nr = *maxresnr;
585 dest->nres += copies * src->nres;
586 dest->nr += copies * src->nr;
589 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t* mtop)
593 init_t_atoms(&atoms, 0, FALSE);
595 int maxresnr = mtop->maxresnr;
596 for (const gmx_molblock_t& molb : mtop->molblock)
598 atomcat(&atoms, &mtop->moltype[molb.type].atoms, molb.nmol, mtop->maxres_renum, &maxresnr);
605 * The cat routines below are old code from src/kernel/topcat.c
608 static void ilistcat(int ftype, InteractionList* dest, const InteractionList& src, int copies, int dnum, int snum)
614 size_t destIndex = dest->iatoms.size();
615 dest->iatoms.resize(dest->iatoms.size() + copies * src.size());
617 for (c = 0; c < copies; c++)
619 for (i = 0; i < src.size();)
621 dest->iatoms[destIndex++] = src.iatoms[i++];
622 for (a = 0; a < nral; a++)
624 dest->iatoms[destIndex++] = dnum + src.iatoms[i++];
631 static void ilistcat(int ftype, t_ilist* dest, const InteractionList& src, int copies, int dnum, int snum)
637 dest->nalloc = dest->nr + copies * src.size();
638 srenew(dest->iatoms, dest->nalloc);
640 for (c = 0; c < copies; c++)
642 for (i = 0; i < src.size();)
644 dest->iatoms[dest->nr++] = src.iatoms[i++];
645 for (a = 0; a < nral; a++)
647 dest->iatoms[dest->nr++] = dnum + src.iatoms[i++];
654 static const t_iparams& getIparams(const InteractionDefinitions& idef, const int index)
656 return idef.iparams[index];
659 static const t_iparams& getIparams(const t_idef& idef, const int index)
661 return idef.iparams[index];
664 static void resizeIParams(std::vector<t_iparams>* iparams, const int newSize)
666 iparams->resize(newSize);
669 static void resizeIParams(t_iparams** iparams, const int newSize)
671 srenew(*iparams, newSize);
674 template<typename IdefType>
675 static void set_posres_params(IdefType* idef, const gmx_molblock_t* molb, int i0, int a_offset)
680 auto* il = &idef->il[F_POSRES];
682 resizeIParams(&idef->iparams_posres, i1);
683 for (i = i0; i < i1; i++)
685 ip = &idef->iparams_posres[i];
686 /* Copy the force constants */
687 *ip = getIparams(*idef, il->iatoms[i * 2]);
688 a_molb = il->iatoms[i * 2 + 1] - a_offset;
689 if (molb->posres_xA.empty())
691 gmx_incons("Position restraint coordinates are missing");
693 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
694 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
695 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
696 if (!molb->posres_xB.empty())
698 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
699 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
700 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
704 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
705 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
706 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
708 /* Set the parameter index for idef->iparams_posre */
709 il->iatoms[i * 2] = i;
713 template<typename IdefType>
714 static void set_fbposres_params(IdefType* idef, const gmx_molblock_t* molb, int i0, int a_offset)
719 auto* il = &idef->il[F_FBPOSRES];
721 resizeIParams(&idef->iparams_fbposres, i1);
722 for (i = i0; i < i1; i++)
724 ip = &idef->iparams_fbposres[i];
725 /* Copy the force constants */
726 *ip = getIparams(*idef, il->iatoms[i * 2]);
727 a_molb = il->iatoms[i * 2 + 1] - a_offset;
728 if (molb->posres_xA.empty())
730 gmx_incons("Position restraint coordinates are missing");
732 /* Take flat-bottom posres reference from normal position restraints */
733 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
734 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
735 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
736 /* Note: no B-type for flat-bottom posres */
738 /* Set the parameter index for idef->iparams_posre */
739 il->iatoms[i * 2] = i;
743 /*! \brief Copy parameters to idef structure from mtop.
745 * Makes a deep copy of the force field parameters data structure from a gmx_mtop_t.
746 * Used to initialize legacy topology types.
748 * \param[in] mtop Reference to input mtop.
749 * \param[in] idef Pointer to idef to populate.
751 static void copyFFParametersFromMtop(const gmx_mtop_t& mtop, t_idef* idef)
753 const gmx_ffparams_t* ffp = &mtop.ffparams;
755 idef->ntypes = ffp->numTypes();
756 idef->atnr = ffp->atnr;
757 /* we can no longer copy the pointers to the mtop members,
758 * because they will become invalid as soon as mtop gets free'd.
759 * We also need to make sure to only operate on valid data!
762 if (!ffp->functype.empty())
764 snew(idef->functype, ffp->functype.size());
765 std::copy(ffp->functype.data(), ffp->functype.data() + ffp->functype.size(), idef->functype);
769 idef->functype = nullptr;
771 if (!ffp->iparams.empty())
773 snew(idef->iparams, ffp->iparams.size());
774 std::copy(ffp->iparams.data(), ffp->iparams.data() + ffp->iparams.size(), idef->iparams);
778 idef->iparams = nullptr;
780 idef->iparams_posres = nullptr;
781 idef->iparams_fbposres = nullptr;
782 idef->fudgeQQ = ffp->fudgeQQ;
783 idef->ilsort = ilsortUNKNOWN;
786 /*! \brief Copy idef structure from mtop.
788 * Makes a deep copy of an idef data structure from a gmx_mtop_t.
789 * Used to initialize legacy topology types.
791 * \param[in] mtop Reference to input mtop.
792 * \param[in] idef Pointer to idef to populate.
793 * \param[in] mergeConstr Decide if constraints will be merged.
795 template<typename IdefType>
796 static void copyIListsFromMtop(const gmx_mtop_t& mtop, IdefType* idef, bool mergeConstr)
799 for (const gmx_molblock_t& molb : mtop.molblock)
801 const gmx_moltype_t& molt = mtop.moltype[molb.type];
803 int srcnr = molt.atoms.nr;
806 int nposre_old = idef->il[F_POSRES].size();
807 int nfbposre_old = idef->il[F_FBPOSRES].size();
808 for (int ftype = 0; ftype < F_NRE; ftype++)
810 if (mergeConstr && ftype == F_CONSTR && !molt.ilist[F_CONSTRNC].empty())
812 /* Merge all constrains into one ilist.
813 * This simplifies the constraint code.
815 for (int mol = 0; mol < molb.nmol; mol++)
817 ilistcat(ftype, &idef->il[F_CONSTR], molt.ilist[F_CONSTR], 1,
818 destnr + mol * srcnr, srcnr);
819 ilistcat(ftype, &idef->il[F_CONSTR], molt.ilist[F_CONSTRNC], 1,
820 destnr + mol * srcnr, srcnr);
823 else if (!(mergeConstr && ftype == F_CONSTRNC))
825 ilistcat(ftype, &idef->il[ftype], molt.ilist[ftype], molb.nmol, destnr, srcnr);
828 if (idef->il[F_POSRES].size() > nposre_old)
830 /* Executing this line line stops gmxdump -sys working
831 * correctly. I'm not aware there's an elegant fix. */
832 set_posres_params(idef, &molb, nposre_old / 2, natoms);
834 if (idef->il[F_FBPOSRES].size() > nfbposre_old)
836 set_fbposres_params(idef, &molb, nfbposre_old / 2, natoms);
839 natoms += molb.nmol * srcnr;
842 if (mtop.bIntermolecularInteractions)
844 for (int ftype = 0; ftype < F_NRE; ftype++)
846 ilistcat(ftype, &idef->il[ftype], (*mtop.intermolecular_ilist)[ftype], 1, 0, mtop.natoms);
850 // We have not (yet) sorted free-energy interactions to the end of the ilists
851 idef->ilsort = ilsortNO_FE;
854 /*! \brief Copy atomtypes from mtop
856 * Makes a deep copy of t_atomtypes from gmx_mtop_t.
857 * Used to initialize legacy topology types.
859 * \param[in] mtop Reference to input mtop.
860 * \param[in] atomtypes Pointer to atomtypes to populate.
862 static void copyAtomtypesFromMtop(const gmx_mtop_t& mtop, t_atomtypes* atomtypes)
864 atomtypes->nr = mtop.atomtypes.nr;
865 if (mtop.atomtypes.atomnumber)
867 snew(atomtypes->atomnumber, mtop.atomtypes.nr);
868 std::copy(mtop.atomtypes.atomnumber, mtop.atomtypes.atomnumber + mtop.atomtypes.nr,
869 atomtypes->atomnumber);
873 atomtypes->atomnumber = nullptr;
877 /*! \brief Generate a single list of lists of exclusions for the whole system
879 * \param[in] mtop Reference to input mtop.
881 static gmx::ListOfLists<int> globalExclusionLists(const gmx_mtop_t& mtop)
883 gmx::ListOfLists<int> excls;
886 for (const gmx_molblock_t& molb : mtop.molblock)
888 const gmx_moltype_t& molt = mtop.moltype[molb.type];
890 for (int mol = 0; mol < molb.nmol; mol++)
892 excls.appendListOfLists(molt.excls, atomIndex);
894 atomIndex += molt.atoms.nr;
901 /*! \brief Updates inter-molecular exclusion lists
903 * This function updates inter-molecular exclusions to exclude all
904 * non-bonded interactions between a given list of atoms
906 * \param[inout] excls existing exclusions in local topology
907 * \param[in] ids list of global IDs of atoms
909 static void addMimicExclusions(gmx::ListOfLists<int>* excls, const gmx::ArrayRef<const int> ids)
911 t_blocka inter_excl{};
912 init_blocka(&inter_excl);
913 size_t n_q = ids.size();
915 inter_excl.nr = excls->ssize();
916 inter_excl.nra = n_q * n_q;
918 size_t total_nra = n_q * n_q;
920 snew(inter_excl.index, excls->ssize() + 1);
921 snew(inter_excl.a, total_nra);
923 for (int i = 0; i < inter_excl.nr; ++i)
925 inter_excl.index[i] = 0;
928 /* Here we loop over the list of QM atom ids
929 * and create exclusions between all of them resulting in
930 * n_q * n_q sized exclusion list
933 for (int k = 0; k < inter_excl.nr; ++k)
935 inter_excl.index[k] = prev_index;
936 for (long i = 0; i < ids.ssize(); i++)
942 size_t index = n_q * i;
943 inter_excl.index[ids[i]] = index;
944 prev_index = index + n_q;
945 for (size_t j = 0; j < n_q; ++j)
947 inter_excl.a[n_q * i + j] = ids[j];
951 inter_excl.index[ids[n_q - 1] + 1] = n_q * n_q;
953 inter_excl.index[inter_excl.nr] = n_q * n_q;
955 std::vector<gmx::ExclusionBlock> qmexcl2(excls->size());
956 gmx::blockaToExclusionBlocks(&inter_excl, qmexcl2);
958 // Merge the created exclusion list with the existing one
959 gmx::mergeExclusions(excls, qmexcl2);
962 static void sortFreeEnergyInteractionsAtEnd(const gmx_mtop_t& mtop, InteractionDefinitions* idef)
964 std::vector<real> qA(mtop.natoms);
965 std::vector<real> qB(mtop.natoms);
966 for (const AtomProxy atomP : AtomRange(mtop))
968 const t_atom& local = atomP.atom();
969 int index = atomP.globalAtomNumber();
971 qB[index] = local.qB;
973 gmx_sort_ilist_fe(idef, qA.data(), qB.data());
976 static void gen_local_top(const gmx_mtop_t& mtop,
977 bool freeEnergyInteractionsAtEnd,
981 copyIListsFromMtop(mtop, &top->idef, bMergeConstr);
982 if (freeEnergyInteractionsAtEnd)
984 sortFreeEnergyInteractionsAtEnd(mtop, &top->idef);
986 top->excls = globalExclusionLists(mtop);
987 if (!mtop.intermolecularExclusionGroup.empty())
989 addMimicExclusions(&top->excls, mtop.intermolecularExclusionGroup);
993 void gmx_mtop_generate_local_top(const gmx_mtop_t& mtop, gmx_localtop_t* top, bool freeEnergyInteractionsAtEnd)
995 gen_local_top(mtop, freeEnergyInteractionsAtEnd, true, top);
998 /*! \brief Fills an array with molecule begin/end atom indices
1000 * \param[in] mtop The global topology
1001 * \param[out] index Array of size nr. of molecules + 1 to be filled with molecule begin/end indices
1003 static void fillMoleculeIndices(const gmx_mtop_t& mtop, gmx::ArrayRef<int> index)
1005 int globalAtomIndex = 0;
1006 int globalMolIndex = 0;
1007 index[globalMolIndex] = globalAtomIndex;
1008 for (const gmx_molblock_t& molb : mtop.molblock)
1010 int numAtomsPerMolecule = mtop.moltype[molb.type].atoms.nr;
1011 for (int mol = 0; mol < molb.nmol; mol++)
1013 globalAtomIndex += numAtomsPerMolecule;
1014 globalMolIndex += 1;
1015 index[globalMolIndex] = globalAtomIndex;
1020 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t& mtop)
1022 gmx::RangePartitioning mols;
1024 for (const gmx_molblock_t& molb : mtop.molblock)
1026 int numAtomsPerMolecule = mtop.moltype[molb.type].atoms.nr;
1027 for (int mol = 0; mol < molb.nmol; mol++)
1029 mols.appendBlock(numAtomsPerMolecule);
1036 /*! \brief Creates and returns a deprecated t_block struct with molecule indices
1038 * \param[in] mtop The global topology
1040 static t_block gmx_mtop_molecules_t_block(const gmx_mtop_t& mtop)
1044 mols.nr = gmx_mtop_num_molecules(mtop);
1045 mols.nalloc_index = mols.nr + 1;
1046 snew(mols.index, mols.nalloc_index);
1048 fillMoleculeIndices(mtop, gmx::arrayRefFromArray(mols.index, mols.nr + 1));
1053 static void gen_t_topology(const gmx_mtop_t& mtop, bool bMergeConstr, t_topology* top)
1055 copyAtomtypesFromMtop(mtop, &top->atomtypes);
1056 for (int ftype = 0; ftype < F_NRE; ftype++)
1058 top->idef.il[ftype].nr = 0;
1059 top->idef.il[ftype].nalloc = 0;
1060 top->idef.il[ftype].iatoms = nullptr;
1062 copyFFParametersFromMtop(mtop, &top->idef);
1063 copyIListsFromMtop(mtop, &top->idef, bMergeConstr);
1065 top->name = mtop.name;
1066 top->atoms = gmx_mtop_global_atoms(&mtop);
1067 top->mols = gmx_mtop_molecules_t_block(mtop);
1068 top->bIntermolecularInteractions = mtop.bIntermolecularInteractions;
1069 top->symtab = mtop.symtab;
1072 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t* mtop, bool freeMTop)
1076 gen_t_topology(*mtop, false, &top);
1080 // Clear pointers and counts, such that the pointers copied to top
1081 // keep pointing to valid data after destroying mtop.
1082 mtop->symtab.symbuf = nullptr;
1083 mtop->symtab.nr = 0;
1088 std::vector<int> get_atom_index(const gmx_mtop_t* mtop)
1091 std::vector<int> atom_index;
1092 for (const AtomProxy atomP : AtomRange(*mtop))
1094 const t_atom& local = atomP.atom();
1095 int index = atomP.globalAtomNumber();
1096 if (local.ptype == eptAtom)
1098 atom_index.push_back(index);
1104 void convertAtomsToMtop(t_symtab* symtab, char** name, t_atoms* atoms, gmx_mtop_t* mtop)
1106 mtop->symtab = *symtab;
1110 mtop->moltype.clear();
1111 mtop->moltype.resize(1);
1112 mtop->moltype.back().atoms = *atoms;
1114 mtop->molblock.resize(1);
1115 mtop->molblock[0].type = 0;
1116 mtop->molblock[0].nmol = 1;
1118 mtop->bIntermolecularInteractions = FALSE;
1120 mtop->natoms = atoms->nr;
1122 mtop->haveMoleculeIndices = false;
1124 gmx_mtop_finalize(mtop);
1127 bool haveFepPerturbedNBInteractions(const gmx_mtop_t* mtop)
1129 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
1131 const gmx_molblock_t& molb = mtop->molblock[mb];
1132 const gmx_moltype_t& molt = mtop->moltype[molb.type];
1133 for (int m = 0; m < molb.nmol; m++)
1135 for (int a = 0; a < molt.atoms.nr; a++)
1137 const t_atom& atom = molt.atoms.atom[a];
1138 if (PERTURBED(atom))