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,2021, 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 void gmx_mtop_count_atomtypes(const gmx_mtop_t& mtop, int state, int typecount[])
62 for (int i = 0; i < mtop.ffparams.atnr; ++i)
66 for (const gmx_molblock_t& molb : mtop.molblock)
68 const t_atoms& atoms = mtop.moltype[molb.type].atoms;
69 for (int i = 0; i < atoms.nr; ++i)
71 const int tpi = (state == 0) ? atoms.atom[i].type : atoms.atom[i].typeB;
72 typecount[tpi] += molb.nmol;
77 int gmx_mtop_num_molecules(const gmx_mtop_t& mtop)
80 for (const gmx_molblock_t& molb : mtop.molblock)
82 numMolecules += molb.nmol;
87 int gmx_mtop_nres(const gmx_mtop_t& mtop)
90 for (const gmx_molblock_t& molb : mtop.molblock)
92 nres += molb.nmol * mtop.moltype[molb.type].atoms.nres;
97 AtomIterator::AtomIterator(const gmx_mtop_t& mtop, int globalAtomNumber) :
100 atoms_(&mtop.moltype[mtop.molblock[0].type].atoms),
102 highestResidueNumber_(mtop.maxResNumberNotRenumbered()),
104 globalAtomNumber_(globalAtomNumber)
106 GMX_ASSERT(globalAtomNumber == 0 || globalAtomNumber == mtop.natoms,
107 "Starting at other atoms not implemented yet");
110 AtomIterator& AtomIterator::operator++()
115 if (localAtomNumber_ >= atoms_->nr)
117 if (atoms_->nres <= mtop_->maxResiduesPerMoleculeToTriggerRenumber())
119 /* Single residue molecule, increase the count with one */
120 highestResidueNumber_ += atoms_->nres;
123 localAtomNumber_ = 0;
124 if (currentMolecule_ >= mtop_->molblock[mblock_].nmol)
127 if (mblock_ >= mtop_->molblock.size())
131 atoms_ = &mtop_->moltype[mtop_->molblock[mblock_].type].atoms;
132 currentMolecule_ = 0;
138 bool AtomIterator::operator==(const AtomIterator& o) const
140 return mtop_ == o.mtop_ && globalAtomNumber_ == o.globalAtomNumber_;
143 const t_atom& AtomProxy::atom() const
145 return it_->atoms_->atom[it_->localAtomNumber_];
148 int AtomProxy::globalAtomNumber() const
150 return it_->globalAtomNumber_;
153 const char* AtomProxy::atomName() const
155 return *(it_->atoms_->atomname[it_->localAtomNumber_]);
158 const char* AtomProxy::residueName() const
160 int residueIndexInMolecule = it_->atoms_->atom[it_->localAtomNumber_].resind;
161 return *(it_->atoms_->resinfo[residueIndexInMolecule].name);
164 int AtomProxy::residueNumber() const
166 int residueIndexInMolecule = it_->atoms_->atom[it_->localAtomNumber_].resind;
167 if (it_->atoms_->nres <= it_->mtop_->maxResiduesPerMoleculeToTriggerRenumber())
169 return it_->highestResidueNumber_ + 1 + residueIndexInMolecule;
173 return it_->atoms_->resinfo[residueIndexInMolecule].nr;
177 const gmx_moltype_t& AtomProxy::moleculeType() const
179 return it_->mtop_->moltype[it_->mtop_->molblock[it_->mblock_].type];
182 int AtomProxy::atomNumberInMol() const
184 return it_->localAtomNumber_;
187 struct gmx_mtop_atomloop_block
189 const gmx_mtop_t* mtop;
191 const t_atoms* atoms;
195 gmx_mtop_atomloop_block_t gmx_mtop_atomloop_block_init(const gmx_mtop_t& mtop)
197 struct gmx_mtop_atomloop_block* aloop = nullptr;
203 aloop->atoms = &mtop.moltype[mtop.molblock[aloop->mblock].type].atoms;
204 aloop->at_local = -1;
209 static void gmx_mtop_atomloop_block_destroy(gmx_mtop_atomloop_block_t aloop)
214 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop, const t_atom** atom, int* nmol)
216 if (aloop == nullptr)
218 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
223 if (aloop->at_local >= aloop->atoms->nr)
226 if (aloop->mblock >= aloop->mtop->molblock.size())
228 gmx_mtop_atomloop_block_destroy(aloop);
231 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
235 *atom = &aloop->atoms->atom[aloop->at_local];
236 *nmol = aloop->mtop->molblock[aloop->mblock].nmol;
241 IListIterator::IListIterator(const gmx_mtop_t& mtop, size_t molblock) :
247 IListIterator& IListIterator::operator++()
253 bool IListIterator::operator==(const IListIterator& o) const
255 return mtop_ == o.mtop_ && mblock_ == o.mblock_;
258 const InteractionLists& IListProxy::list() const
260 // one past the end means we want to take the
261 // intermolecular list instead.
262 if (it_->mblock_ == it_->mtop_->molblock.size())
264 return *it_->mtop_->intermolecular_ilist;
268 return it_->mtop_->moltype[it_->mtop_->molblock[it_->mblock_].type].ilist;
272 int IListProxy::nmol() const
274 // one past the end means we want to take the
275 // intermolecular list instead.
276 if (it_->mblock_ == it_->mtop_->molblock.size())
282 return it_->mtop_->molblock[it_->mblock_].nmol;
286 IListRange::IListRange(const gmx_mtop_t& mtop) : begin_(mtop), end_(mtop, mtop.molblock.size())
288 if (mtop.bIntermolecularInteractions)
290 end_ = IListIterator(mtop, mtop.molblock.size() + 1);
294 int gmx_mtop_ftype_count(const gmx_mtop_t& mtop, int ftype)
298 for (const IListProxy il : IListRange(mtop))
300 n += il.nmol() * il.list()[ftype].size() / (1 + NRAL(ftype));
303 if (mtop.bIntermolecularInteractions)
305 n += (*mtop.intermolecular_ilist)[ftype].size() / (1 + NRAL(ftype));
311 int gmx_mtop_interaction_count(const gmx_mtop_t& mtop, const int unsigned if_flags)
315 for (const IListProxy il : IListRange(mtop))
317 for (int ftype = 0; ftype < F_NRE; ftype++)
319 if ((interaction_function[ftype].flags & if_flags) == if_flags)
321 n += il.nmol() * il.list()[ftype].size() / (1 + NRAL(ftype));
326 if (mtop.bIntermolecularInteractions)
328 for (int ftype = 0; ftype < F_NRE; ftype++)
330 if ((interaction_function[ftype].flags & if_flags) == if_flags)
332 n += (*mtop.intermolecular_ilist)[ftype].size() / (1 + NRAL(ftype));
340 gmx::EnumerationArray<ParticleType, int> gmx_mtop_particletype_count(const gmx_mtop_t& mtop)
342 gmx::EnumerationArray<ParticleType, int> count = { { 0 } };
344 for (const auto& molblock : mtop.molblock)
346 const t_atoms& atoms = mtop.moltype[molblock.type].atoms;
347 for (int a = 0; a < atoms.nr; a++)
349 count[atoms.atom[a].ptype] += molblock.nmol;
356 static void atomcat(t_atoms* dest, const t_atoms* src, int copies, int maxres_renum, int* maxresnr)
358 int i = 0, j = 0, l = 0, size = 0;
360 int destnr = dest->nr;
364 dest->haveMass = src->haveMass;
365 dest->haveType = src->haveType;
366 dest->haveCharge = src->haveCharge;
367 dest->haveBState = src->haveBState;
368 dest->havePdbInfo = src->havePdbInfo;
372 dest->haveMass = dest->haveMass && src->haveMass;
373 dest->haveType = dest->haveType && src->haveType;
374 dest->haveCharge = dest->haveCharge && src->haveCharge;
375 dest->haveBState = dest->haveBState && src->haveBState;
376 dest->havePdbInfo = dest->havePdbInfo && src->havePdbInfo;
381 size = destnr + copies * srcnr;
382 srenew(dest->atom, size);
383 srenew(dest->atomname, size);
386 srenew(dest->atomtype, size);
387 if (dest->haveBState)
389 srenew(dest->atomtypeB, size);
392 if (dest->havePdbInfo)
394 srenew(dest->pdbinfo, size);
399 size = dest->nres + copies * src->nres;
400 srenew(dest->resinfo, size);
403 /* residue information */
404 for (l = dest->nres, j = 0; (j < copies); j++, l += src->nres)
406 memcpy(reinterpret_cast<char*>(&(dest->resinfo[l])),
407 reinterpret_cast<char*>(&(src->resinfo[0])),
408 static_cast<size_t>(src->nres * sizeof(src->resinfo[0])));
411 for (l = destnr, j = 0; (j < copies); j++, l += srcnr)
413 memcpy(reinterpret_cast<char*>(&(dest->atom[l])),
414 reinterpret_cast<char*>(&(src->atom[0])),
415 static_cast<size_t>(srcnr * sizeof(src->atom[0])));
416 memcpy(reinterpret_cast<char*>(&(dest->atomname[l])),
417 reinterpret_cast<char*>(&(src->atomname[0])),
418 static_cast<size_t>(srcnr * sizeof(src->atomname[0])));
421 memcpy(reinterpret_cast<char*>(&(dest->atomtype[l])),
422 reinterpret_cast<char*>(&(src->atomtype[0])),
423 static_cast<size_t>(srcnr * sizeof(src->atomtype[0])));
424 if (dest->haveBState)
426 memcpy(reinterpret_cast<char*>(&(dest->atomtypeB[l])),
427 reinterpret_cast<char*>(&(src->atomtypeB[0])),
428 static_cast<size_t>(srcnr * sizeof(src->atomtypeB[0])));
431 if (dest->havePdbInfo)
433 memcpy(reinterpret_cast<char*>(&(dest->pdbinfo[l])),
434 reinterpret_cast<char*>(&(src->pdbinfo[0])),
435 static_cast<size_t>(srcnr * sizeof(src->pdbinfo[0])));
439 /* Increment residue indices */
440 for (l = destnr, j = 0; (j < copies); j++)
442 for (i = 0; (i < srcnr); i++, l++)
444 dest->atom[l].resind = dest->nres + j * src->nres + src->atom[i].resind;
448 if (src->nres <= maxres_renum)
450 /* Single residue molecule, continue counting residues */
451 for (j = 0; (j < copies); j++)
453 for (l = 0; l < src->nres; l++)
456 dest->resinfo[dest->nres + j * src->nres + l].nr = *maxresnr;
461 dest->nres += copies * src->nres;
462 dest->nr += copies * src->nr;
465 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t& mtop)
469 init_t_atoms(&atoms, 0, FALSE);
471 int maxresnr = mtop.maxResNumberNotRenumbered();
472 for (const gmx_molblock_t& molb : mtop.molblock)
475 &mtop.moltype[molb.type].atoms,
477 mtop.maxResiduesPerMoleculeToTriggerRenumber(),
485 * The cat routines below are old code from src/kernel/topcat.c
488 static void ilistcat(int ftype, InteractionList* dest, const InteractionList& src, int copies, int dnum, int snum)
490 const int nral = NRAL(ftype);
492 size_t destIndex = dest->iatoms.size();
493 dest->iatoms.resize(dest->iatoms.size() + copies * src.size());
495 for (int c = 0; c < copies; c++)
497 for (int i = 0; i < src.size();)
499 dest->iatoms[destIndex++] = src.iatoms[i++];
500 for (int a = 0; a < nral; a++)
502 dest->iatoms[destIndex++] = dnum + src.iatoms[i++];
509 static void ilistcat(int ftype, t_ilist* dest, const InteractionList& src, int copies, int dnum, int snum)
511 const int nral = NRAL(ftype);
513 dest->nalloc = dest->nr + copies * src.size();
514 srenew(dest->iatoms, dest->nalloc);
516 for (int c = 0; c < copies; c++)
518 for (int i = 0; i < src.size();)
520 dest->iatoms[dest->nr++] = src.iatoms[i++];
521 for (int a = 0; a < nral; a++)
523 dest->iatoms[dest->nr++] = dnum + src.iatoms[i++];
530 static const t_iparams& getIparams(const InteractionDefinitions& idef, const int index)
532 return idef.iparams[index];
535 static const t_iparams& getIparams(const t_idef& idef, const int index)
537 return idef.iparams[index];
540 static void resizeIParams(std::vector<t_iparams>* iparams, const int newSize)
542 iparams->resize(newSize);
545 static void resizeIParams(t_iparams** iparams, const int newSize)
547 srenew(*iparams, newSize);
550 template<typename IdefType>
551 static void set_posres_params(IdefType* idef, const gmx_molblock_t* molb, int i0, int a_offset)
553 auto* il = &idef->il[F_POSRES];
554 int i1 = il->size() / 2;
555 resizeIParams(&idef->iparams_posres, i1);
556 for (int i = i0; i < i1; i++)
558 t_iparams* ip = &idef->iparams_posres[i];
559 /* Copy the force constants */
560 *ip = getIparams(*idef, il->iatoms[i * 2]);
561 int a_molb = il->iatoms[i * 2 + 1] - a_offset;
562 if (molb->posres_xA.empty())
564 gmx_incons("Position restraint coordinates are missing");
566 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
567 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
568 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
569 if (!molb->posres_xB.empty())
571 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
572 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
573 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
577 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
578 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
579 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
581 /* Set the parameter index for idef->iparams_posre */
582 il->iatoms[i * 2] = i;
586 template<typename IdefType>
587 static void set_fbposres_params(IdefType* idef, const gmx_molblock_t* molb, int i0, int a_offset)
589 auto* il = &idef->il[F_FBPOSRES];
590 int i1 = il->size() / 2;
591 resizeIParams(&idef->iparams_fbposres, i1);
592 for (int i = i0; i < i1; i++)
594 t_iparams* ip = &idef->iparams_fbposres[i];
595 /* Copy the force constants */
596 *ip = getIparams(*idef, il->iatoms[i * 2]);
597 int a_molb = il->iatoms[i * 2 + 1] - a_offset;
598 if (molb->posres_xA.empty())
600 gmx_incons("Position restraint coordinates are missing");
602 /* Take flat-bottom posres reference from normal position restraints */
603 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
604 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
605 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
606 /* Note: no B-type for flat-bottom posres */
608 /* Set the parameter index for idef->iparams_posre */
609 il->iatoms[i * 2] = i;
613 /*! \brief Copy parameters to idef structure from mtop.
615 * Makes a deep copy of the force field parameters data structure from a gmx_mtop_t.
616 * Used to initialize legacy topology types.
618 * \param[in] mtop Reference to input mtop.
619 * \param[in] idef Pointer to idef to populate.
621 static void copyFFParametersFromMtop(const gmx_mtop_t& mtop, t_idef* idef)
623 const gmx_ffparams_t* ffp = &mtop.ffparams;
625 idef->ntypes = ffp->numTypes();
626 idef->atnr = ffp->atnr;
627 /* we can no longer copy the pointers to the mtop members,
628 * because they will become invalid as soon as mtop gets free'd.
629 * We also need to make sure to only operate on valid data!
632 if (!ffp->functype.empty())
634 snew(idef->functype, ffp->functype.size());
635 std::copy(ffp->functype.data(), ffp->functype.data() + ffp->functype.size(), idef->functype);
639 idef->functype = nullptr;
641 if (!ffp->iparams.empty())
643 snew(idef->iparams, ffp->iparams.size());
644 std::copy(ffp->iparams.data(), ffp->iparams.data() + ffp->iparams.size(), idef->iparams);
648 idef->iparams = nullptr;
650 idef->iparams_posres = nullptr;
651 idef->iparams_fbposres = nullptr;
652 idef->fudgeQQ = ffp->fudgeQQ;
653 idef->ilsort = ilsortUNKNOWN;
656 /*! \brief Copy idef structure from mtop.
658 * Makes a deep copy of an idef data structure from a gmx_mtop_t.
659 * Used to initialize legacy topology types.
661 * \param[in] mtop Reference to input mtop.
662 * \param[in] idef Pointer to idef to populate.
663 * \param[in] mergeConstr Decide if constraints will be merged.
665 template<typename IdefType>
666 static void copyIListsFromMtop(const gmx_mtop_t& mtop, IdefType* idef, bool mergeConstr)
669 for (const gmx_molblock_t& molb : mtop.molblock)
671 const gmx_moltype_t& molt = mtop.moltype[molb.type];
673 int srcnr = molt.atoms.nr;
676 int nposre_old = idef->il[F_POSRES].size();
677 int nfbposre_old = idef->il[F_FBPOSRES].size();
678 for (int ftype = 0; ftype < F_NRE; ftype++)
680 if (mergeConstr && ftype == F_CONSTR && !molt.ilist[F_CONSTRNC].empty())
682 /* Merge all constrains into one ilist.
683 * This simplifies the constraint code.
685 for (int mol = 0; mol < molb.nmol; mol++)
687 ilistcat(ftype, &idef->il[F_CONSTR], molt.ilist[F_CONSTR], 1, destnr + mol * srcnr, srcnr);
688 ilistcat(ftype, &idef->il[F_CONSTR], molt.ilist[F_CONSTRNC], 1, destnr + mol * srcnr, srcnr);
691 else if (!(mergeConstr && ftype == F_CONSTRNC))
693 ilistcat(ftype, &idef->il[ftype], molt.ilist[ftype], molb.nmol, destnr, srcnr);
696 if (idef->il[F_POSRES].size() > nposre_old)
698 /* Executing this line line stops gmxdump -sys working
699 * correctly. I'm not aware there's an elegant fix. */
700 set_posres_params(idef, &molb, nposre_old / 2, natoms);
702 if (idef->il[F_FBPOSRES].size() > nfbposre_old)
704 set_fbposres_params(idef, &molb, nfbposre_old / 2, natoms);
707 natoms += molb.nmol * srcnr;
710 if (mtop.bIntermolecularInteractions)
712 for (int ftype = 0; ftype < F_NRE; ftype++)
714 ilistcat(ftype, &idef->il[ftype], (*mtop.intermolecular_ilist)[ftype], 1, 0, mtop.natoms);
718 // We have not (yet) sorted free-energy interactions to the end of the ilists
719 idef->ilsort = ilsortNO_FE;
722 /*! \brief Copy atomtypes from mtop
724 * Makes a deep copy of t_atomtypes from gmx_mtop_t.
725 * Used to initialize legacy topology types.
727 * \param[in] mtop Reference to input mtop.
728 * \param[in] atomtypes Pointer to atomtypes to populate.
730 static void copyAtomtypesFromMtop(const gmx_mtop_t& mtop, t_atomtypes* atomtypes)
732 atomtypes->nr = mtop.atomtypes.nr;
733 if (mtop.atomtypes.atomnumber)
735 snew(atomtypes->atomnumber, mtop.atomtypes.nr);
736 std::copy(mtop.atomtypes.atomnumber,
737 mtop.atomtypes.atomnumber + mtop.atomtypes.nr,
738 atomtypes->atomnumber);
742 atomtypes->atomnumber = nullptr;
746 /*! \brief Generate a single list of lists of exclusions for the whole system
748 * \param[in] mtop Reference to input mtop.
750 static gmx::ListOfLists<int> globalExclusionLists(const gmx_mtop_t& mtop)
752 gmx::ListOfLists<int> excls;
755 for (const gmx_molblock_t& molb : mtop.molblock)
757 const gmx_moltype_t& molt = mtop.moltype[molb.type];
759 for (int mol = 0; mol < molb.nmol; mol++)
761 excls.appendListOfLists(molt.excls, atomIndex);
763 atomIndex += molt.atoms.nr;
770 /*! \brief Updates inter-molecular exclusion lists
772 * This function updates inter-molecular exclusions to exclude all
773 * non-bonded interactions between a given list of atoms
775 * \param[inout] excls existing exclusions in local topology
776 * \param[in] ids list of global IDs of atoms
778 static void addMimicExclusions(gmx::ListOfLists<int>* excls, const gmx::ArrayRef<const int> ids)
780 t_blocka inter_excl{};
781 init_blocka(&inter_excl);
782 size_t n_q = ids.size();
784 inter_excl.nr = excls->ssize();
785 inter_excl.nra = n_q * n_q;
787 size_t total_nra = n_q * n_q;
789 snew(inter_excl.index, excls->ssize() + 1);
790 snew(inter_excl.a, total_nra);
792 for (int i = 0; i < inter_excl.nr; ++i)
794 inter_excl.index[i] = 0;
797 /* Here we loop over the list of QM atom ids
798 * and create exclusions between all of them resulting in
799 * n_q * n_q sized exclusion list
802 for (int k = 0; k < inter_excl.nr; ++k)
804 inter_excl.index[k] = prev_index;
805 for (long i = 0; i < ids.ssize(); i++)
811 size_t index = n_q * i;
812 inter_excl.index[ids[i]] = index;
813 prev_index = index + n_q;
814 for (size_t j = 0; j < n_q; ++j)
816 inter_excl.a[n_q * i + j] = ids[j];
820 inter_excl.index[ids[n_q - 1] + 1] = n_q * n_q;
822 inter_excl.index[inter_excl.nr] = n_q * n_q;
824 std::vector<gmx::ExclusionBlock> qmexcl2(excls->size());
825 gmx::blockaToExclusionBlocks(&inter_excl, qmexcl2);
827 // Merge the created exclusion list with the existing one
828 gmx::mergeExclusions(excls, qmexcl2);
831 static void sortFreeEnergyInteractionsAtEnd(const gmx_mtop_t& mtop, InteractionDefinitions* idef)
833 std::vector<real> qA(mtop.natoms);
834 std::vector<real> qB(mtop.natoms);
835 for (const AtomProxy atomP : AtomRange(mtop))
837 const t_atom& local = atomP.atom();
838 int index = atomP.globalAtomNumber();
840 qB[index] = local.qB;
842 gmx_sort_ilist_fe(idef, qA.data(), qB.data());
845 static void gen_local_top(const gmx_mtop_t& mtop,
846 bool freeEnergyInteractionsAtEnd,
850 copyIListsFromMtop(mtop, &top->idef, bMergeConstr);
851 if (freeEnergyInteractionsAtEnd)
853 sortFreeEnergyInteractionsAtEnd(mtop, &top->idef);
855 top->excls = globalExclusionLists(mtop);
856 if (!mtop.intermolecularExclusionGroup.empty())
858 addMimicExclusions(&top->excls, mtop.intermolecularExclusionGroup);
862 void gmx_mtop_generate_local_top(const gmx_mtop_t& mtop, gmx_localtop_t* top, bool freeEnergyInteractionsAtEnd)
864 gen_local_top(mtop, freeEnergyInteractionsAtEnd, true, top);
867 /*! \brief Fills an array with molecule begin/end atom indices
869 * \param[in] mtop The global topology
870 * \param[out] index Array of size nr. of molecules + 1 to be filled with molecule begin/end indices
872 static void fillMoleculeIndices(const gmx_mtop_t& mtop, gmx::ArrayRef<int> index)
874 int globalAtomIndex = 0;
875 int globalMolIndex = 0;
876 index[globalMolIndex] = globalAtomIndex;
877 for (const gmx_molblock_t& molb : mtop.molblock)
879 int numAtomsPerMolecule = mtop.moltype[molb.type].atoms.nr;
880 for (int mol = 0; mol < molb.nmol; mol++)
882 globalAtomIndex += numAtomsPerMolecule;
884 index[globalMolIndex] = globalAtomIndex;
889 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t& mtop)
891 gmx::RangePartitioning mols;
893 for (const gmx_molblock_t& molb : mtop.molblock)
895 int numAtomsPerMolecule = mtop.moltype[molb.type].atoms.nr;
896 for (int mol = 0; mol < molb.nmol; mol++)
898 mols.appendBlock(numAtomsPerMolecule);
905 std::vector<gmx::Range<int>> atomRangeOfEachResidue(const gmx_moltype_t& moltype)
907 std::vector<gmx::Range<int>> atomRanges;
908 int currentResidueNumber = moltype.atoms.atom[0].resind;
910 // Go through all atoms in a molecule to store first and last atoms in each residue.
911 for (int i = 0; i < moltype.atoms.nr; i++)
913 int residueOfThisAtom = moltype.atoms.atom[i].resind;
914 if (residueOfThisAtom != currentResidueNumber)
916 // This atom belongs to the next residue, so record the range for the previous residue,
917 // remembering that end points to one place past the last atom.
919 atomRanges.emplace_back(startAtom, endAtom);
920 // Prepare for the current residue
922 currentResidueNumber = residueOfThisAtom;
925 // special treatment for last residue in this molecule.
926 atomRanges.emplace_back(startAtom, moltype.atoms.nr);
931 /*! \brief Creates and returns a deprecated t_block struct with molecule indices
933 * \param[in] mtop The global topology
935 static t_block gmx_mtop_molecules_t_block(const gmx_mtop_t& mtop)
939 mols.nr = gmx_mtop_num_molecules(mtop);
940 mols.nalloc_index = mols.nr + 1;
941 snew(mols.index, mols.nalloc_index);
943 fillMoleculeIndices(mtop, gmx::arrayRefFromArray(mols.index, mols.nr + 1));
948 static void gen_t_topology(const gmx_mtop_t& mtop, bool bMergeConstr, t_topology* top)
950 copyAtomtypesFromMtop(mtop, &top->atomtypes);
951 for (int ftype = 0; ftype < F_NRE; ftype++)
953 top->idef.il[ftype].nr = 0;
954 top->idef.il[ftype].nalloc = 0;
955 top->idef.il[ftype].iatoms = nullptr;
957 copyFFParametersFromMtop(mtop, &top->idef);
958 copyIListsFromMtop(mtop, &top->idef, bMergeConstr);
960 top->name = mtop.name;
961 top->atoms = gmx_mtop_global_atoms(mtop);
962 top->mols = gmx_mtop_molecules_t_block(mtop);
963 top->bIntermolecularInteractions = mtop.bIntermolecularInteractions;
964 top->symtab = mtop.symtab;
967 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t* mtop, bool freeMTop)
971 gen_t_topology(*mtop, false, &top);
975 // Clear pointers and counts, such that the pointers copied to top
976 // keep pointing to valid data after destroying mtop.
977 mtop->symtab.symbuf = nullptr;
983 std::vector<int> get_atom_index(const gmx_mtop_t& mtop)
986 std::vector<int> atom_index;
987 for (const AtomProxy atomP : AtomRange(mtop))
989 const t_atom& local = atomP.atom();
990 int index = atomP.globalAtomNumber();
991 if (local.ptype == ParticleType::Atom)
993 atom_index.push_back(index);
999 void convertAtomsToMtop(t_symtab* symtab, char** name, t_atoms* atoms, gmx_mtop_t* mtop)
1001 mtop->symtab = *symtab;
1005 mtop->moltype.clear();
1006 mtop->moltype.resize(1);
1007 mtop->moltype.back().atoms = *atoms;
1009 mtop->molblock.resize(1);
1010 mtop->molblock[0].type = 0;
1011 mtop->molblock[0].nmol = 1;
1013 mtop->bIntermolecularInteractions = FALSE;
1015 mtop->natoms = atoms->nr;
1017 mtop->haveMoleculeIndices = false;
1022 bool haveFepPerturbedNBInteractions(const gmx_mtop_t& mtop)
1024 for (const gmx_moltype_t& molt : mtop.moltype)
1026 for (int a = 0; a < molt.atoms.nr; a++)
1028 if (PERTURBED(molt.atoms.atom[a]))
1038 bool haveFepPerturbedMasses(const gmx_mtop_t& mtop)
1040 for (const gmx_moltype_t& molt : mtop.moltype)
1042 for (int a = 0; a < molt.atoms.nr; a++)
1044 const t_atom& atom = molt.atoms.atom[a];
1045 if (atom.m != atom.mB)
1055 bool havePerturbedConstraints(const gmx_mtop_t& mtop)
1057 // This code assumes that all perturbed constraints parameters are actually used
1058 const auto& ffparams = mtop.ffparams;
1060 for (gmx::index i = 0; i < gmx::ssize(ffparams.functype); i++)
1062 if (ffparams.functype[i] == F_CONSTR || ffparams.functype[i] == F_CONSTRNC)
1064 const auto& iparams = ffparams.iparams[i];
1065 if (iparams.constr.dA != iparams.constr.dB)