2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2008,2009,2010.
5 * Copyright (c) 2012,2013,2014,2015,2016 by 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 gmx_mtop_ilistloop_all_t gmx_mtop_ilistloop_all_init(const gmx_mtop_t* mtop)
411 struct gmx_mtop_ilistloop_all* iloop;
423 static void gmx_mtop_ilistloop_all_destroy(gmx_mtop_ilistloop_all_t iloop)
428 const InteractionLists* gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop, int* atnr_offset)
431 if (iloop == nullptr)
434 "gmx_mtop_ilistloop_all_next called without calling gmx_mtop_ilistloop_all_init");
439 iloop->a_offset += iloop->mtop->moleculeBlockIndices[iloop->mblock].numAtomsPerMolecule;
444 /* Inter-molecular interactions, if present, are indexed with
445 * iloop->mblock == iloop->mtop->nmolblock, thus we should separately
446 * check for this value in this conditional.
448 if (iloop->mblock == iloop->mtop->molblock.size()
449 || iloop->mol >= iloop->mtop->molblock[iloop->mblock].nmol)
453 if (iloop->mblock >= iloop->mtop->molblock.size())
455 if (iloop->mblock == iloop->mtop->molblock.size() && iloop->mtop->bIntermolecularInteractions)
458 return iloop->mtop->intermolecular_ilist.get();
461 gmx_mtop_ilistloop_all_destroy(iloop);
466 *atnr_offset = iloop->a_offset;
468 return &iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
471 int gmx_mtop_ftype_count(const gmx_mtop_t* mtop, int ftype)
473 gmx_mtop_ilistloop_t iloop;
478 iloop = gmx_mtop_ilistloop_init(mtop);
479 while (const InteractionLists* il = gmx_mtop_ilistloop_next(iloop, &nmol))
481 n += nmol * (*il)[ftype].size() / (1 + NRAL(ftype));
484 if (mtop->bIntermolecularInteractions)
486 n += (*mtop->intermolecular_ilist)[ftype].size() / (1 + NRAL(ftype));
492 int gmx_mtop_ftype_count(const gmx_mtop_t& mtop, int ftype)
494 return gmx_mtop_ftype_count(&mtop, ftype);
497 int gmx_mtop_interaction_count(const gmx_mtop_t& mtop, const int unsigned if_flags)
501 gmx_mtop_ilistloop_t iloop = gmx_mtop_ilistloop_init(mtop);
503 while (const InteractionLists* il = gmx_mtop_ilistloop_next(iloop, &nmol))
505 for (int ftype = 0; ftype < F_NRE; ftype++)
507 if ((interaction_function[ftype].flags & if_flags) == if_flags)
509 n += nmol * (*il)[ftype].size() / (1 + NRAL(ftype));
514 if (mtop.bIntermolecularInteractions)
516 for (int ftype = 0; ftype < F_NRE; ftype++)
518 if ((interaction_function[ftype].flags & if_flags) == if_flags)
520 n += (*mtop.intermolecular_ilist)[ftype].size() / (1 + NRAL(ftype));
528 static void atomcat(t_atoms* dest, const t_atoms* src, int copies, int maxres_renum, int* maxresnr)
532 int destnr = dest->nr;
536 dest->haveMass = src->haveMass;
537 dest->haveType = src->haveType;
538 dest->haveCharge = src->haveCharge;
539 dest->haveBState = src->haveBState;
540 dest->havePdbInfo = src->havePdbInfo;
544 dest->haveMass = dest->haveMass && src->haveMass;
545 dest->haveType = dest->haveType && src->haveType;
546 dest->haveCharge = dest->haveCharge && src->haveCharge;
547 dest->haveBState = dest->haveBState && src->haveBState;
548 dest->havePdbInfo = dest->havePdbInfo && src->havePdbInfo;
553 size = destnr + copies * srcnr;
554 srenew(dest->atom, size);
555 srenew(dest->atomname, size);
558 srenew(dest->atomtype, size);
559 if (dest->haveBState)
561 srenew(dest->atomtypeB, size);
564 if (dest->havePdbInfo)
566 srenew(dest->pdbinfo, size);
571 size = dest->nres + copies * src->nres;
572 srenew(dest->resinfo, size);
575 /* residue information */
576 for (l = dest->nres, j = 0; (j < copies); j++, l += src->nres)
578 memcpy(reinterpret_cast<char*>(&(dest->resinfo[l])), reinterpret_cast<char*>(&(src->resinfo[0])),
579 static_cast<size_t>(src->nres * sizeof(src->resinfo[0])));
582 for (l = destnr, j = 0; (j < copies); j++, l += srcnr)
584 memcpy(reinterpret_cast<char*>(&(dest->atom[l])), reinterpret_cast<char*>(&(src->atom[0])),
585 static_cast<size_t>(srcnr * sizeof(src->atom[0])));
586 memcpy(reinterpret_cast<char*>(&(dest->atomname[l])),
587 reinterpret_cast<char*>(&(src->atomname[0])),
588 static_cast<size_t>(srcnr * sizeof(src->atomname[0])));
591 memcpy(reinterpret_cast<char*>(&(dest->atomtype[l])),
592 reinterpret_cast<char*>(&(src->atomtype[0])),
593 static_cast<size_t>(srcnr * sizeof(src->atomtype[0])));
594 if (dest->haveBState)
596 memcpy(reinterpret_cast<char*>(&(dest->atomtypeB[l])),
597 reinterpret_cast<char*>(&(src->atomtypeB[0])),
598 static_cast<size_t>(srcnr * sizeof(src->atomtypeB[0])));
601 if (dest->havePdbInfo)
603 memcpy(reinterpret_cast<char*>(&(dest->pdbinfo[l])),
604 reinterpret_cast<char*>(&(src->pdbinfo[0])),
605 static_cast<size_t>(srcnr * sizeof(src->pdbinfo[0])));
609 /* Increment residue indices */
610 for (l = destnr, j = 0; (j < copies); j++)
612 for (i = 0; (i < srcnr); i++, l++)
614 dest->atom[l].resind = dest->nres + j * src->nres + src->atom[i].resind;
618 if (src->nres <= maxres_renum)
620 /* Single residue molecule, continue counting residues */
621 for (j = 0; (j < copies); j++)
623 for (l = 0; l < src->nres; l++)
626 dest->resinfo[dest->nres + j * src->nres + l].nr = *maxresnr;
631 dest->nres += copies * src->nres;
632 dest->nr += copies * src->nr;
635 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t* mtop)
639 init_t_atoms(&atoms, 0, FALSE);
641 int maxresnr = mtop->maxresnr;
642 for (const gmx_molblock_t& molb : mtop->molblock)
644 atomcat(&atoms, &mtop->moltype[molb.type].atoms, molb.nmol, mtop->maxres_renum, &maxresnr);
651 * The cat routines below are old code from src/kernel/topcat.c
654 static void ilistcat(int ftype, t_ilist* dest, const InteractionList& src, int copies, int dnum, int snum)
660 dest->nalloc = dest->nr + copies * src.size();
661 srenew(dest->iatoms, dest->nalloc);
663 for (c = 0; c < copies; c++)
665 for (i = 0; i < src.size();)
667 dest->iatoms[dest->nr++] = src.iatoms[i++];
668 for (a = 0; a < nral; a++)
670 dest->iatoms[dest->nr++] = dnum + src.iatoms[i++];
677 static void set_posres_params(t_idef* idef, const gmx_molblock_t* molb, int i0, int a_offset)
683 il = &idef->il[F_POSRES];
685 idef->iparams_posres_nalloc = i1;
686 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
687 for (i = i0; i < i1; i++)
689 ip = &idef->iparams_posres[i];
690 /* Copy the force constants */
691 *ip = idef->iparams[il->iatoms[i * 2]];
692 a_molb = il->iatoms[i * 2 + 1] - a_offset;
693 if (molb->posres_xA.empty())
695 gmx_incons("Position restraint coordinates are missing");
697 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
698 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
699 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
700 if (!molb->posres_xB.empty())
702 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
703 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
704 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
708 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
709 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
710 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
712 /* Set the parameter index for idef->iparams_posre */
713 il->iatoms[i * 2] = i;
717 static void set_fbposres_params(t_idef* idef, const gmx_molblock_t* molb, int i0, int a_offset)
723 il = &idef->il[F_FBPOSRES];
725 idef->iparams_fbposres_nalloc = i1;
726 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
727 for (i = i0; i < i1; i++)
729 ip = &idef->iparams_fbposres[i];
730 /* Copy the force constants */
731 *ip = idef->iparams[il->iatoms[i * 2]];
732 a_molb = il->iatoms[i * 2 + 1] - a_offset;
733 if (molb->posres_xA.empty())
735 gmx_incons("Position restraint coordinates are missing");
737 /* Take flat-bottom posres reference from normal position restraints */
738 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
739 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
740 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
741 /* Note: no B-type for flat-bottom posres */
743 /* Set the parameter index for idef->iparams_posre */
744 il->iatoms[i * 2] = i;
748 /*! \brief Copy idef structure from mtop.
750 * Makes a deep copy of an idef data structure from a gmx_mtop_t.
751 * Used to initialize legacy topology types.
753 * \param[in] mtop Reference to input mtop.
754 * \param[in] idef Pointer to idef to populate.
755 * \param[in] mergeConstr Decide if constraints will be merged.
756 * \param[in] freeEnergyInteractionsAtEnd Decide if free energy stuff should
757 * be added at the end.
759 static void copyIdefFromMtop(const gmx_mtop_t& mtop, t_idef* idef, bool freeEnergyInteractionsAtEnd, bool mergeConstr)
761 const gmx_ffparams_t* ffp = &mtop.ffparams;
763 idef->ntypes = ffp->numTypes();
764 idef->atnr = ffp->atnr;
765 /* we can no longer copy the pointers to the mtop members,
766 * because they will become invalid as soon as mtop gets free'd.
767 * We also need to make sure to only operate on valid data!
770 if (!ffp->functype.empty())
772 snew(idef->functype, ffp->functype.size());
773 std::copy(ffp->functype.data(), ffp->functype.data() + ffp->functype.size(), idef->functype);
777 idef->functype = nullptr;
779 if (!ffp->iparams.empty())
781 snew(idef->iparams, ffp->iparams.size());
782 std::copy(ffp->iparams.data(), ffp->iparams.data() + ffp->iparams.size(), idef->iparams);
786 idef->iparams = nullptr;
788 idef->iparams_posres = nullptr;
789 idef->iparams_posres_nalloc = 0;
790 idef->iparams_fbposres = nullptr;
791 idef->iparams_fbposres_nalloc = 0;
792 idef->fudgeQQ = ffp->fudgeQQ;
793 idef->cmap_grid = new gmx_cmap_t;
794 *idef->cmap_grid = ffp->cmap_grid;
795 idef->ilsort = ilsortUNKNOWN;
797 for (int ftype = 0; ftype < F_NRE; ftype++)
799 idef->il[ftype].nr = 0;
800 idef->il[ftype].nalloc = 0;
801 idef->il[ftype].iatoms = nullptr;
805 for (const gmx_molblock_t& molb : mtop.molblock)
807 const gmx_moltype_t& molt = mtop.moltype[molb.type];
809 int srcnr = molt.atoms.nr;
812 int nposre_old = idef->il[F_POSRES].nr;
813 int nfbposre_old = idef->il[F_FBPOSRES].nr;
814 for (int ftype = 0; ftype < F_NRE; ftype++)
816 if (mergeConstr && ftype == F_CONSTR && molt.ilist[F_CONSTRNC].size() > 0)
818 /* Merge all constrains into one ilist.
819 * This simplifies the constraint code.
821 for (int mol = 0; mol < molb.nmol; mol++)
823 ilistcat(ftype, &idef->il[F_CONSTR], molt.ilist[F_CONSTR], 1,
824 destnr + mol * srcnr, srcnr);
825 ilistcat(ftype, &idef->il[F_CONSTR], molt.ilist[F_CONSTRNC], 1,
826 destnr + mol * srcnr, srcnr);
829 else if (!(mergeConstr && ftype == F_CONSTRNC))
831 ilistcat(ftype, &idef->il[ftype], molt.ilist[ftype], molb.nmol, destnr, srcnr);
834 if (idef->il[F_POSRES].nr > nposre_old)
836 /* Executing this line line stops gmxdump -sys working
837 * correctly. I'm not aware there's an elegant fix. */
838 set_posres_params(idef, &molb, nposre_old / 2, natoms);
840 if (idef->il[F_FBPOSRES].nr > nfbposre_old)
842 set_fbposres_params(idef, &molb, nfbposre_old / 2, natoms);
845 natoms += molb.nmol * srcnr;
848 if (mtop.bIntermolecularInteractions)
850 for (int ftype = 0; ftype < F_NRE; ftype++)
852 ilistcat(ftype, &idef->il[ftype], (*mtop.intermolecular_ilist)[ftype], 1, 0, mtop.natoms);
856 if (freeEnergyInteractionsAtEnd && gmx_mtop_bondeds_free_energy(&mtop))
858 std::vector<real> qA(mtop.natoms);
859 std::vector<real> qB(mtop.natoms);
860 for (const AtomProxy atomP : AtomRange(mtop))
862 const t_atom& local = atomP.atom();
863 int index = atomP.globalAtomNumber();
865 qB[index] = local.qB;
867 gmx_sort_ilist_fe(idef, qA.data(), qB.data());
871 idef->ilsort = ilsortNO_FE;
875 /*! \brief Copy atomtypes from mtop
877 * Makes a deep copy of t_atomtypes from gmx_mtop_t.
878 * Used to initialize legacy topology types.
880 * \param[in] mtop Reference to input mtop.
881 * \param[in] atomtypes Pointer to atomtypes to populate.
883 static void copyAtomtypesFromMtop(const gmx_mtop_t& mtop, t_atomtypes* atomtypes)
885 atomtypes->nr = mtop.atomtypes.nr;
886 if (mtop.atomtypes.atomnumber)
888 snew(atomtypes->atomnumber, mtop.atomtypes.nr);
889 std::copy(mtop.atomtypes.atomnumber, mtop.atomtypes.atomnumber + mtop.atomtypes.nr,
890 atomtypes->atomnumber);
894 atomtypes->atomnumber = nullptr;
898 /*! \brief Generate a single list of lists of exclusions for the whole system
900 * \param[in] mtop Reference to input mtop.
902 static gmx::ListOfLists<int> globalExclusionLists(const gmx_mtop_t& mtop)
904 gmx::ListOfLists<int> excls;
907 for (const gmx_molblock_t& molb : mtop.molblock)
909 const gmx_moltype_t& molt = mtop.moltype[molb.type];
911 for (int mol = 0; mol < molb.nmol; mol++)
913 excls.appendListOfLists(molt.excls, atomIndex);
915 atomIndex += molt.atoms.nr;
922 /*! \brief Updates inter-molecular exclusion lists
924 * This function updates inter-molecular exclusions to exclude all
925 * non-bonded interactions between a given list of atoms
927 * \param[inout] excls existing exclusions in local topology
928 * \param[in] ids list of global IDs of atoms
930 static void addMimicExclusions(gmx::ListOfLists<int>* excls, const gmx::ArrayRef<const int> ids)
932 t_blocka inter_excl{};
933 init_blocka(&inter_excl);
934 size_t n_q = ids.size();
936 inter_excl.nr = excls->ssize();
937 inter_excl.nra = n_q * n_q;
939 size_t total_nra = n_q * n_q;
941 snew(inter_excl.index, excls->ssize() + 1);
942 snew(inter_excl.a, total_nra);
944 for (int i = 0; i < inter_excl.nr; ++i)
946 inter_excl.index[i] = 0;
949 /* Here we loop over the list of QM atom ids
950 * and create exclusions between all of them resulting in
951 * n_q * n_q sized exclusion list
954 for (int k = 0; k < inter_excl.nr; ++k)
956 inter_excl.index[k] = prev_index;
957 for (long i = 0; i < ids.ssize(); i++)
963 size_t index = n_q * i;
964 inter_excl.index[ids[i]] = index;
965 prev_index = index + n_q;
966 for (size_t j = 0; j < n_q; ++j)
968 inter_excl.a[n_q * i + j] = ids[j];
972 inter_excl.index[ids[n_q - 1] + 1] = n_q * n_q;
974 inter_excl.index[inter_excl.nr] = n_q * n_q;
976 std::vector<gmx::ExclusionBlock> qmexcl2(excls->size());
977 gmx::blockaToExclusionBlocks(&inter_excl, qmexcl2);
979 // Merge the created exclusion list with the existing one
980 gmx::mergeExclusions(excls, qmexcl2);
983 static void gen_local_top(const gmx_mtop_t& mtop,
984 bool freeEnergyInteractionsAtEnd,
988 copyAtomtypesFromMtop(mtop, &top->atomtypes);
989 copyIdefFromMtop(mtop, &top->idef, freeEnergyInteractionsAtEnd, bMergeConstr);
990 top->excls = globalExclusionLists(mtop);
991 if (!mtop.intermolecularExclusionGroup.empty())
993 addMimicExclusions(&top->excls, mtop.intermolecularExclusionGroup);
997 void gmx_mtop_generate_local_top(const gmx_mtop_t& mtop, gmx_localtop_t* top, bool freeEnergyInteractionsAtEnd)
999 gen_local_top(mtop, freeEnergyInteractionsAtEnd, true, top);
1002 /*! \brief Fills an array with molecule begin/end atom indices
1004 * \param[in] mtop The global topology
1005 * \param[out] index Array of size nr. of molecules + 1 to be filled with molecule begin/end indices
1007 static void fillMoleculeIndices(const gmx_mtop_t& mtop, gmx::ArrayRef<int> index)
1009 int globalAtomIndex = 0;
1010 int globalMolIndex = 0;
1011 index[globalMolIndex] = globalAtomIndex;
1012 for (const gmx_molblock_t& molb : mtop.molblock)
1014 int numAtomsPerMolecule = mtop.moltype[molb.type].atoms.nr;
1015 for (int mol = 0; mol < molb.nmol; mol++)
1017 globalAtomIndex += numAtomsPerMolecule;
1018 globalMolIndex += 1;
1019 index[globalMolIndex] = globalAtomIndex;
1024 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t& mtop)
1026 gmx::RangePartitioning mols;
1028 for (const gmx_molblock_t& molb : mtop.molblock)
1030 int numAtomsPerMolecule = mtop.moltype[molb.type].atoms.nr;
1031 for (int mol = 0; mol < molb.nmol; mol++)
1033 mols.appendBlock(numAtomsPerMolecule);
1040 /*! \brief Creates and returns a deprecated t_block struct with molecule indices
1042 * \param[in] mtop The global topology
1044 static t_block gmx_mtop_molecules_t_block(const gmx_mtop_t& mtop)
1048 mols.nr = gmx_mtop_num_molecules(mtop);
1049 mols.nalloc_index = mols.nr + 1;
1050 snew(mols.index, mols.nalloc_index);
1052 fillMoleculeIndices(mtop, gmx::arrayRefFromArray(mols.index, mols.nr + 1));
1057 static void gen_t_topology(const gmx_mtop_t& mtop,
1058 bool freeEnergyInteractionsAtEnd,
1062 copyAtomtypesFromMtop(mtop, &top->atomtypes);
1063 copyIdefFromMtop(mtop, &top->idef, freeEnergyInteractionsAtEnd, 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, 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);