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 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)
74 tpi = atoms.atom[i].type;
78 tpi = atoms.atom[i].typeB;
80 typecount[tpi] += molb.nmol;
85 int gmx_mtop_num_molecules(const gmx_mtop_t& mtop)
88 for (const gmx_molblock_t& molb : mtop.molblock)
90 numMolecules += molb.nmol;
95 int gmx_mtop_nres(const gmx_mtop_t* mtop)
98 for (const gmx_molblock_t& molb : mtop->molblock)
100 nres += molb.nmol * mtop->moltype[molb.type].atoms.nres;
105 AtomIterator::AtomIterator(const gmx_mtop_t& mtop, int globalAtomNumber) :
108 atoms_(&mtop.moltype[mtop.molblock[0].type].atoms),
110 highestResidueNumber_(mtop.maxResNumberNotRenumbered()),
112 globalAtomNumber_(globalAtomNumber)
114 GMX_ASSERT(globalAtomNumber == 0 || globalAtomNumber == mtop.natoms,
115 "Starting at other atoms not implemented yet");
118 AtomIterator& AtomIterator::operator++()
123 if (localAtomNumber_ >= atoms_->nr)
125 if (atoms_->nres <= mtop_->maxResiduesPerMoleculeToTriggerRenumber())
127 /* Single residue molecule, increase the count with one */
128 highestResidueNumber_ += atoms_->nres;
131 localAtomNumber_ = 0;
132 if (currentMolecule_ >= mtop_->molblock[mblock_].nmol)
135 if (mblock_ >= mtop_->molblock.size())
139 atoms_ = &mtop_->moltype[mtop_->molblock[mblock_].type].atoms;
140 currentMolecule_ = 0;
146 bool AtomIterator::operator==(const AtomIterator& o) const
148 return mtop_ == o.mtop_ && globalAtomNumber_ == o.globalAtomNumber_;
151 const t_atom& AtomProxy::atom() const
153 return it_->atoms_->atom[it_->localAtomNumber_];
156 int AtomProxy::globalAtomNumber() const
158 return it_->globalAtomNumber_;
161 const char* AtomProxy::atomName() const
163 return *(it_->atoms_->atomname[it_->localAtomNumber_]);
166 const char* AtomProxy::residueName() const
168 int residueIndexInMolecule = it_->atoms_->atom[it_->localAtomNumber_].resind;
169 return *(it_->atoms_->resinfo[residueIndexInMolecule].name);
172 int AtomProxy::residueNumber() const
174 int residueIndexInMolecule = it_->atoms_->atom[it_->localAtomNumber_].resind;
175 if (it_->atoms_->nres <= it_->mtop_->maxResiduesPerMoleculeToTriggerRenumber())
177 return it_->highestResidueNumber_ + 1 + residueIndexInMolecule;
181 return it_->atoms_->resinfo[residueIndexInMolecule].nr;
185 const gmx_moltype_t& AtomProxy::moleculeType() const
187 return it_->mtop_->moltype[it_->mtop_->molblock[it_->mblock_].type];
190 int AtomProxy::atomNumberInMol() const
192 return it_->localAtomNumber_;
195 typedef struct gmx_mtop_atomloop_block
197 const gmx_mtop_t* mtop;
199 const t_atoms* atoms;
201 } t_gmx_mtop_atomloop_block;
203 gmx_mtop_atomloop_block_t gmx_mtop_atomloop_block_init(const gmx_mtop_t* mtop)
205 struct gmx_mtop_atomloop_block* aloop;
211 aloop->atoms = &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
212 aloop->at_local = -1;
217 static void gmx_mtop_atomloop_block_destroy(gmx_mtop_atomloop_block_t aloop)
222 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop, const t_atom** atom, int* nmol)
224 if (aloop == nullptr)
226 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
231 if (aloop->at_local >= aloop->atoms->nr)
234 if (aloop->mblock >= aloop->mtop->molblock.size())
236 gmx_mtop_atomloop_block_destroy(aloop);
239 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
243 *atom = &aloop->atoms->atom[aloop->at_local];
244 *nmol = aloop->mtop->molblock[aloop->mblock].nmol;
249 typedef struct gmx_mtop_ilistloop
251 const gmx_mtop_t* mtop;
255 gmx_mtop_ilistloop_t gmx_mtop_ilistloop_init(const gmx_mtop_t* mtop)
257 struct gmx_mtop_ilistloop* iloop;
267 gmx_mtop_ilistloop_t gmx_mtop_ilistloop_init(const gmx_mtop_t& mtop)
269 return gmx_mtop_ilistloop_init(&mtop);
272 static void gmx_mtop_ilistloop_destroy(gmx_mtop_ilistloop_t iloop)
277 const InteractionLists* gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop, int* nmol)
279 if (iloop == nullptr)
281 gmx_incons("gmx_mtop_ilistloop_next called without calling gmx_mtop_ilistloop_init");
285 if (iloop->mblock >= gmx::ssize(iloop->mtop->molblock))
287 if (iloop->mblock == gmx::ssize(iloop->mtop->molblock) && iloop->mtop->bIntermolecularInteractions)
290 return iloop->mtop->intermolecular_ilist.get();
293 gmx_mtop_ilistloop_destroy(iloop);
297 *nmol = iloop->mtop->molblock[iloop->mblock].nmol;
299 return &iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
301 typedef struct gmx_mtop_ilistloop_all
303 const gmx_mtop_t* mtop;
307 } t_gmx_mtop_ilist_all;
309 int gmx_mtop_ftype_count(const gmx_mtop_t* mtop, int ftype)
311 gmx_mtop_ilistloop_t iloop;
316 iloop = gmx_mtop_ilistloop_init(mtop);
317 while (const InteractionLists* il = gmx_mtop_ilistloop_next(iloop, &nmol))
319 n += nmol * (*il)[ftype].size() / (1 + NRAL(ftype));
322 if (mtop->bIntermolecularInteractions)
324 n += (*mtop->intermolecular_ilist)[ftype].size() / (1 + NRAL(ftype));
330 int gmx_mtop_ftype_count(const gmx_mtop_t& mtop, int ftype)
332 return gmx_mtop_ftype_count(&mtop, ftype);
335 int gmx_mtop_interaction_count(const gmx_mtop_t& mtop, const int unsigned if_flags)
339 gmx_mtop_ilistloop_t iloop = gmx_mtop_ilistloop_init(mtop);
341 while (const InteractionLists* il = gmx_mtop_ilistloop_next(iloop, &nmol))
343 for (int ftype = 0; ftype < F_NRE; ftype++)
345 if ((interaction_function[ftype].flags & if_flags) == if_flags)
347 n += nmol * (*il)[ftype].size() / (1 + NRAL(ftype));
352 if (mtop.bIntermolecularInteractions)
354 for (int ftype = 0; ftype < F_NRE; ftype++)
356 if ((interaction_function[ftype].flags & if_flags) == if_flags)
358 n += (*mtop.intermolecular_ilist)[ftype].size() / (1 + NRAL(ftype));
366 std::array<int, eptNR> gmx_mtop_particletype_count(const gmx_mtop_t& mtop)
368 std::array<int, eptNR> count = { { 0 } };
370 for (const auto& molblock : mtop.molblock)
372 const t_atoms& atoms = mtop.moltype[molblock.type].atoms;
373 for (int a = 0; a < atoms.nr; a++)
375 count[atoms.atom[a].ptype] += molblock.nmol;
382 static void atomcat(t_atoms* dest, const t_atoms* src, int copies, int maxres_renum, int* maxresnr)
386 int destnr = dest->nr;
390 dest->haveMass = src->haveMass;
391 dest->haveType = src->haveType;
392 dest->haveCharge = src->haveCharge;
393 dest->haveBState = src->haveBState;
394 dest->havePdbInfo = src->havePdbInfo;
398 dest->haveMass = dest->haveMass && src->haveMass;
399 dest->haveType = dest->haveType && src->haveType;
400 dest->haveCharge = dest->haveCharge && src->haveCharge;
401 dest->haveBState = dest->haveBState && src->haveBState;
402 dest->havePdbInfo = dest->havePdbInfo && src->havePdbInfo;
407 size = destnr + copies * srcnr;
408 srenew(dest->atom, size);
409 srenew(dest->atomname, size);
412 srenew(dest->atomtype, size);
413 if (dest->haveBState)
415 srenew(dest->atomtypeB, size);
418 if (dest->havePdbInfo)
420 srenew(dest->pdbinfo, size);
425 size = dest->nres + copies * src->nres;
426 srenew(dest->resinfo, size);
429 /* residue information */
430 for (l = dest->nres, j = 0; (j < copies); j++, l += src->nres)
432 memcpy(reinterpret_cast<char*>(&(dest->resinfo[l])), reinterpret_cast<char*>(&(src->resinfo[0])),
433 static_cast<size_t>(src->nres * sizeof(src->resinfo[0])));
436 for (l = destnr, j = 0; (j < copies); j++, l += srcnr)
438 memcpy(reinterpret_cast<char*>(&(dest->atom[l])), reinterpret_cast<char*>(&(src->atom[0])),
439 static_cast<size_t>(srcnr * sizeof(src->atom[0])));
440 memcpy(reinterpret_cast<char*>(&(dest->atomname[l])),
441 reinterpret_cast<char*>(&(src->atomname[0])),
442 static_cast<size_t>(srcnr * sizeof(src->atomname[0])));
445 memcpy(reinterpret_cast<char*>(&(dest->atomtype[l])),
446 reinterpret_cast<char*>(&(src->atomtype[0])),
447 static_cast<size_t>(srcnr * sizeof(src->atomtype[0])));
448 if (dest->haveBState)
450 memcpy(reinterpret_cast<char*>(&(dest->atomtypeB[l])),
451 reinterpret_cast<char*>(&(src->atomtypeB[0])),
452 static_cast<size_t>(srcnr * sizeof(src->atomtypeB[0])));
455 if (dest->havePdbInfo)
457 memcpy(reinterpret_cast<char*>(&(dest->pdbinfo[l])),
458 reinterpret_cast<char*>(&(src->pdbinfo[0])),
459 static_cast<size_t>(srcnr * sizeof(src->pdbinfo[0])));
463 /* Increment residue indices */
464 for (l = destnr, j = 0; (j < copies); j++)
466 for (i = 0; (i < srcnr); i++, l++)
468 dest->atom[l].resind = dest->nres + j * src->nres + src->atom[i].resind;
472 if (src->nres <= maxres_renum)
474 /* Single residue molecule, continue counting residues */
475 for (j = 0; (j < copies); j++)
477 for (l = 0; l < src->nres; l++)
480 dest->resinfo[dest->nres + j * src->nres + l].nr = *maxresnr;
485 dest->nres += copies * src->nres;
486 dest->nr += copies * src->nr;
489 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t* mtop)
493 init_t_atoms(&atoms, 0, FALSE);
495 int maxresnr = mtop->maxResNumberNotRenumbered();
496 for (const gmx_molblock_t& molb : mtop->molblock)
498 atomcat(&atoms, &mtop->moltype[molb.type].atoms, molb.nmol,
499 mtop->maxResiduesPerMoleculeToTriggerRenumber(), &maxresnr);
506 * The cat routines below are old code from src/kernel/topcat.c
509 static void ilistcat(int ftype, InteractionList* dest, const InteractionList& src, int copies, int dnum, int snum)
515 size_t destIndex = dest->iatoms.size();
516 dest->iatoms.resize(dest->iatoms.size() + copies * src.size());
518 for (c = 0; c < copies; c++)
520 for (i = 0; i < src.size();)
522 dest->iatoms[destIndex++] = src.iatoms[i++];
523 for (a = 0; a < nral; a++)
525 dest->iatoms[destIndex++] = dnum + src.iatoms[i++];
532 static void ilistcat(int ftype, t_ilist* dest, const InteractionList& src, int copies, int dnum, int snum)
538 dest->nalloc = dest->nr + copies * src.size();
539 srenew(dest->iatoms, dest->nalloc);
541 for (c = 0; c < copies; c++)
543 for (i = 0; i < src.size();)
545 dest->iatoms[dest->nr++] = src.iatoms[i++];
546 for (a = 0; a < nral; a++)
548 dest->iatoms[dest->nr++] = dnum + src.iatoms[i++];
555 static const t_iparams& getIparams(const InteractionDefinitions& idef, const int index)
557 return idef.iparams[index];
560 static const t_iparams& getIparams(const t_idef& idef, const int index)
562 return idef.iparams[index];
565 static void resizeIParams(std::vector<t_iparams>* iparams, const int newSize)
567 iparams->resize(newSize);
570 static void resizeIParams(t_iparams** iparams, const int newSize)
572 srenew(*iparams, newSize);
575 template<typename IdefType>
576 static void set_posres_params(IdefType* idef, const gmx_molblock_t* molb, int i0, int a_offset)
581 auto* il = &idef->il[F_POSRES];
583 resizeIParams(&idef->iparams_posres, i1);
584 for (i = i0; i < i1; i++)
586 ip = &idef->iparams_posres[i];
587 /* Copy the force constants */
588 *ip = getIparams(*idef, il->iatoms[i * 2]);
589 a_molb = il->iatoms[i * 2 + 1] - a_offset;
590 if (molb->posres_xA.empty())
592 gmx_incons("Position restraint coordinates are missing");
594 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
595 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
596 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
597 if (!molb->posres_xB.empty())
599 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
600 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
601 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
605 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
606 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
607 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
609 /* Set the parameter index for idef->iparams_posre */
610 il->iatoms[i * 2] = i;
614 template<typename IdefType>
615 static void set_fbposres_params(IdefType* idef, const gmx_molblock_t* molb, int i0, int a_offset)
620 auto* il = &idef->il[F_FBPOSRES];
622 resizeIParams(&idef->iparams_fbposres, i1);
623 for (i = i0; i < i1; i++)
625 ip = &idef->iparams_fbposres[i];
626 /* Copy the force constants */
627 *ip = getIparams(*idef, il->iatoms[i * 2]);
628 a_molb = il->iatoms[i * 2 + 1] - a_offset;
629 if (molb->posres_xA.empty())
631 gmx_incons("Position restraint coordinates are missing");
633 /* Take flat-bottom posres reference from normal position restraints */
634 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
635 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
636 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
637 /* Note: no B-type for flat-bottom posres */
639 /* Set the parameter index for idef->iparams_posre */
640 il->iatoms[i * 2] = i;
644 /*! \brief Copy parameters to idef structure from mtop.
646 * Makes a deep copy of the force field parameters data structure from a gmx_mtop_t.
647 * Used to initialize legacy topology types.
649 * \param[in] mtop Reference to input mtop.
650 * \param[in] idef Pointer to idef to populate.
652 static void copyFFParametersFromMtop(const gmx_mtop_t& mtop, t_idef* idef)
654 const gmx_ffparams_t* ffp = &mtop.ffparams;
656 idef->ntypes = ffp->numTypes();
657 idef->atnr = ffp->atnr;
658 /* we can no longer copy the pointers to the mtop members,
659 * because they will become invalid as soon as mtop gets free'd.
660 * We also need to make sure to only operate on valid data!
663 if (!ffp->functype.empty())
665 snew(idef->functype, ffp->functype.size());
666 std::copy(ffp->functype.data(), ffp->functype.data() + ffp->functype.size(), idef->functype);
670 idef->functype = nullptr;
672 if (!ffp->iparams.empty())
674 snew(idef->iparams, ffp->iparams.size());
675 std::copy(ffp->iparams.data(), ffp->iparams.data() + ffp->iparams.size(), idef->iparams);
679 idef->iparams = nullptr;
681 idef->iparams_posres = nullptr;
682 idef->iparams_fbposres = nullptr;
683 idef->fudgeQQ = ffp->fudgeQQ;
684 idef->ilsort = ilsortUNKNOWN;
687 /*! \brief Copy idef structure from mtop.
689 * Makes a deep copy of an idef data structure from a gmx_mtop_t.
690 * Used to initialize legacy topology types.
692 * \param[in] mtop Reference to input mtop.
693 * \param[in] idef Pointer to idef to populate.
694 * \param[in] mergeConstr Decide if constraints will be merged.
696 template<typename IdefType>
697 static void copyIListsFromMtop(const gmx_mtop_t& mtop, IdefType* idef, bool mergeConstr)
700 for (const gmx_molblock_t& molb : mtop.molblock)
702 const gmx_moltype_t& molt = mtop.moltype[molb.type];
704 int srcnr = molt.atoms.nr;
707 int nposre_old = idef->il[F_POSRES].size();
708 int nfbposre_old = idef->il[F_FBPOSRES].size();
709 for (int ftype = 0; ftype < F_NRE; ftype++)
711 if (mergeConstr && ftype == F_CONSTR && !molt.ilist[F_CONSTRNC].empty())
713 /* Merge all constrains into one ilist.
714 * This simplifies the constraint code.
716 for (int mol = 0; mol < molb.nmol; mol++)
718 ilistcat(ftype, &idef->il[F_CONSTR], molt.ilist[F_CONSTR], 1,
719 destnr + mol * srcnr, srcnr);
720 ilistcat(ftype, &idef->il[F_CONSTR], molt.ilist[F_CONSTRNC], 1,
721 destnr + mol * srcnr, srcnr);
724 else if (!(mergeConstr && ftype == F_CONSTRNC))
726 ilistcat(ftype, &idef->il[ftype], molt.ilist[ftype], molb.nmol, destnr, srcnr);
729 if (idef->il[F_POSRES].size() > nposre_old)
731 /* Executing this line line stops gmxdump -sys working
732 * correctly. I'm not aware there's an elegant fix. */
733 set_posres_params(idef, &molb, nposre_old / 2, natoms);
735 if (idef->il[F_FBPOSRES].size() > nfbposre_old)
737 set_fbposres_params(idef, &molb, nfbposre_old / 2, natoms);
740 natoms += molb.nmol * srcnr;
743 if (mtop.bIntermolecularInteractions)
745 for (int ftype = 0; ftype < F_NRE; ftype++)
747 ilistcat(ftype, &idef->il[ftype], (*mtop.intermolecular_ilist)[ftype], 1, 0, mtop.natoms);
751 // We have not (yet) sorted free-energy interactions to the end of the ilists
752 idef->ilsort = ilsortNO_FE;
755 /*! \brief Copy atomtypes from mtop
757 * Makes a deep copy of t_atomtypes from gmx_mtop_t.
758 * Used to initialize legacy topology types.
760 * \param[in] mtop Reference to input mtop.
761 * \param[in] atomtypes Pointer to atomtypes to populate.
763 static void copyAtomtypesFromMtop(const gmx_mtop_t& mtop, t_atomtypes* atomtypes)
765 atomtypes->nr = mtop.atomtypes.nr;
766 if (mtop.atomtypes.atomnumber)
768 snew(atomtypes->atomnumber, mtop.atomtypes.nr);
769 std::copy(mtop.atomtypes.atomnumber, mtop.atomtypes.atomnumber + mtop.atomtypes.nr,
770 atomtypes->atomnumber);
774 atomtypes->atomnumber = nullptr;
778 /*! \brief Generate a single list of lists of exclusions for the whole system
780 * \param[in] mtop Reference to input mtop.
782 static gmx::ListOfLists<int> globalExclusionLists(const gmx_mtop_t& mtop)
784 gmx::ListOfLists<int> excls;
787 for (const gmx_molblock_t& molb : mtop.molblock)
789 const gmx_moltype_t& molt = mtop.moltype[molb.type];
791 for (int mol = 0; mol < molb.nmol; mol++)
793 excls.appendListOfLists(molt.excls, atomIndex);
795 atomIndex += molt.atoms.nr;
802 /*! \brief Updates inter-molecular exclusion lists
804 * This function updates inter-molecular exclusions to exclude all
805 * non-bonded interactions between a given list of atoms
807 * \param[inout] excls existing exclusions in local topology
808 * \param[in] ids list of global IDs of atoms
810 static void addMimicExclusions(gmx::ListOfLists<int>* excls, const gmx::ArrayRef<const int> ids)
812 t_blocka inter_excl{};
813 init_blocka(&inter_excl);
814 size_t n_q = ids.size();
816 inter_excl.nr = excls->ssize();
817 inter_excl.nra = n_q * n_q;
819 size_t total_nra = n_q * n_q;
821 snew(inter_excl.index, excls->ssize() + 1);
822 snew(inter_excl.a, total_nra);
824 for (int i = 0; i < inter_excl.nr; ++i)
826 inter_excl.index[i] = 0;
829 /* Here we loop over the list of QM atom ids
830 * and create exclusions between all of them resulting in
831 * n_q * n_q sized exclusion list
834 for (int k = 0; k < inter_excl.nr; ++k)
836 inter_excl.index[k] = prev_index;
837 for (long i = 0; i < ids.ssize(); i++)
843 size_t index = n_q * i;
844 inter_excl.index[ids[i]] = index;
845 prev_index = index + n_q;
846 for (size_t j = 0; j < n_q; ++j)
848 inter_excl.a[n_q * i + j] = ids[j];
852 inter_excl.index[ids[n_q - 1] + 1] = n_q * n_q;
854 inter_excl.index[inter_excl.nr] = n_q * n_q;
856 std::vector<gmx::ExclusionBlock> qmexcl2(excls->size());
857 gmx::blockaToExclusionBlocks(&inter_excl, qmexcl2);
859 // Merge the created exclusion list with the existing one
860 gmx::mergeExclusions(excls, qmexcl2);
863 static void sortFreeEnergyInteractionsAtEnd(const gmx_mtop_t& mtop, InteractionDefinitions* idef)
865 std::vector<real> qA(mtop.natoms);
866 std::vector<real> qB(mtop.natoms);
867 for (const AtomProxy atomP : AtomRange(mtop))
869 const t_atom& local = atomP.atom();
870 int index = atomP.globalAtomNumber();
872 qB[index] = local.qB;
874 gmx_sort_ilist_fe(idef, qA.data(), qB.data());
877 static void gen_local_top(const gmx_mtop_t& mtop,
878 bool freeEnergyInteractionsAtEnd,
882 copyIListsFromMtop(mtop, &top->idef, bMergeConstr);
883 if (freeEnergyInteractionsAtEnd)
885 sortFreeEnergyInteractionsAtEnd(mtop, &top->idef);
887 top->excls = globalExclusionLists(mtop);
888 if (!mtop.intermolecularExclusionGroup.empty())
890 addMimicExclusions(&top->excls, mtop.intermolecularExclusionGroup);
894 void gmx_mtop_generate_local_top(const gmx_mtop_t& mtop, gmx_localtop_t* top, bool freeEnergyInteractionsAtEnd)
896 gen_local_top(mtop, freeEnergyInteractionsAtEnd, true, top);
899 /*! \brief Fills an array with molecule begin/end atom indices
901 * \param[in] mtop The global topology
902 * \param[out] index Array of size nr. of molecules + 1 to be filled with molecule begin/end indices
904 static void fillMoleculeIndices(const gmx_mtop_t& mtop, gmx::ArrayRef<int> index)
906 int globalAtomIndex = 0;
907 int globalMolIndex = 0;
908 index[globalMolIndex] = globalAtomIndex;
909 for (const gmx_molblock_t& molb : mtop.molblock)
911 int numAtomsPerMolecule = mtop.moltype[molb.type].atoms.nr;
912 for (int mol = 0; mol < molb.nmol; mol++)
914 globalAtomIndex += numAtomsPerMolecule;
916 index[globalMolIndex] = globalAtomIndex;
921 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t& mtop)
923 gmx::RangePartitioning mols;
925 for (const gmx_molblock_t& molb : mtop.molblock)
927 int numAtomsPerMolecule = mtop.moltype[molb.type].atoms.nr;
928 for (int mol = 0; mol < molb.nmol; mol++)
930 mols.appendBlock(numAtomsPerMolecule);
937 std::vector<gmx::Range<int>> atomRangeOfEachResidue(const gmx_moltype_t& moltype)
939 std::vector<gmx::Range<int>> atomRanges;
940 int currentResidueNumber = moltype.atoms.atom[0].resind;
942 // Go through all atoms in a molecule to store first and last atoms in each residue.
943 for (int i = 0; i < moltype.atoms.nr; i++)
945 int residueOfThisAtom = moltype.atoms.atom[i].resind;
946 if (residueOfThisAtom != currentResidueNumber)
948 // This atom belongs to the next residue, so record the range for the previous residue,
949 // remembering that end points to one place past the last atom.
951 atomRanges.emplace_back(startAtom, endAtom);
952 // Prepare for the current residue
954 currentResidueNumber = residueOfThisAtom;
957 // special treatment for last residue in this molecule.
958 atomRanges.emplace_back(startAtom, moltype.atoms.nr);
963 /*! \brief Creates and returns a deprecated t_block struct with molecule indices
965 * \param[in] mtop The global topology
967 static t_block gmx_mtop_molecules_t_block(const gmx_mtop_t& mtop)
971 mols.nr = gmx_mtop_num_molecules(mtop);
972 mols.nalloc_index = mols.nr + 1;
973 snew(mols.index, mols.nalloc_index);
975 fillMoleculeIndices(mtop, gmx::arrayRefFromArray(mols.index, mols.nr + 1));
980 static void gen_t_topology(const gmx_mtop_t& mtop, bool bMergeConstr, t_topology* top)
982 copyAtomtypesFromMtop(mtop, &top->atomtypes);
983 for (int ftype = 0; ftype < F_NRE; ftype++)
985 top->idef.il[ftype].nr = 0;
986 top->idef.il[ftype].nalloc = 0;
987 top->idef.il[ftype].iatoms = nullptr;
989 copyFFParametersFromMtop(mtop, &top->idef);
990 copyIListsFromMtop(mtop, &top->idef, bMergeConstr);
992 top->name = mtop.name;
993 top->atoms = gmx_mtop_global_atoms(&mtop);
994 top->mols = gmx_mtop_molecules_t_block(mtop);
995 top->bIntermolecularInteractions = mtop.bIntermolecularInteractions;
996 top->symtab = mtop.symtab;
999 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t* mtop, bool freeMTop)
1003 gen_t_topology(*mtop, false, &top);
1007 // Clear pointers and counts, such that the pointers copied to top
1008 // keep pointing to valid data after destroying mtop.
1009 mtop->symtab.symbuf = nullptr;
1010 mtop->symtab.nr = 0;
1015 std::vector<int> get_atom_index(const gmx_mtop_t* mtop)
1018 std::vector<int> atom_index;
1019 for (const AtomProxy atomP : AtomRange(*mtop))
1021 const t_atom& local = atomP.atom();
1022 int index = atomP.globalAtomNumber();
1023 if (local.ptype == eptAtom)
1025 atom_index.push_back(index);
1031 void convertAtomsToMtop(t_symtab* symtab, char** name, t_atoms* atoms, gmx_mtop_t* mtop)
1033 mtop->symtab = *symtab;
1037 mtop->moltype.clear();
1038 mtop->moltype.resize(1);
1039 mtop->moltype.back().atoms = *atoms;
1041 mtop->molblock.resize(1);
1042 mtop->molblock[0].type = 0;
1043 mtop->molblock[0].nmol = 1;
1045 mtop->bIntermolecularInteractions = FALSE;
1047 mtop->natoms = atoms->nr;
1049 mtop->haveMoleculeIndices = false;
1054 bool haveFepPerturbedNBInteractions(const gmx_mtop_t& mtop)
1056 for (const gmx_moltype_t& molt : mtop.moltype)
1058 for (int a = 0; a < molt.atoms.nr; a++)
1060 if (PERTURBED(molt.atoms.atom[a]))
1070 bool haveFepPerturbedMasses(const gmx_mtop_t& mtop)
1072 for (const gmx_moltype_t& molt : mtop.moltype)
1074 for (int a = 0; a < molt.atoms.nr; a++)
1076 const t_atom& atom = molt.atoms.atom[a];
1077 if (atom.m != atom.mB)
1087 bool havePerturbedConstraints(const gmx_mtop_t& mtop)
1089 // This code assumes that all perturbed constraints parameters are actually used
1090 const auto& ffparams = mtop.ffparams;
1092 for (gmx::index i = 0; i < gmx::ssize(ffparams.functype); i++)
1094 if (ffparams.functype[i] == F_CONSTR || ffparams.functype[i] == F_CONSTRNC)
1096 const auto& iparams = ffparams.iparams[i];
1097 if (iparams.constr.dA != iparams.constr.dB)