Use more modern clang for linters and sanitizers
[alexxy/gromacs.git] / src / gromacs / topology / mtop_util.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "mtop_util.h"
40
41 #include <climits>
42 #include <cstddef>
43 #include <cstdio>
44 #include <cstdlib>
45 #include <cstring>
46
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/enumerationhelpers.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/real.h"
59 #include "gromacs/utility/smalloc.h"
60
61 void gmx_mtop_count_atomtypes(const gmx_mtop_t& mtop, int state, int typecount[])
62 {
63     for (int i = 0; i < mtop.ffparams.atnr; ++i)
64     {
65         typecount[i] = 0;
66     }
67     for (const gmx_molblock_t& molb : mtop.molblock)
68     {
69         const t_atoms& atoms = mtop.moltype[molb.type].atoms;
70         for (int i = 0; i < atoms.nr; ++i)
71         {
72             const int tpi = (state == 0) ? atoms.atom[i].type : atoms.atom[i].typeB;
73             typecount[tpi] += molb.nmol;
74         }
75     }
76 }
77
78 int gmx_mtop_num_molecules(const gmx_mtop_t& mtop)
79 {
80     int numMolecules = 0;
81     for (const gmx_molblock_t& molb : mtop.molblock)
82     {
83         numMolecules += molb.nmol;
84     }
85     return numMolecules;
86 }
87
88 int gmx_mtop_nres(const gmx_mtop_t& mtop)
89 {
90     int nres = 0;
91     for (const gmx_molblock_t& molb : mtop.molblock)
92     {
93         nres += molb.nmol * mtop.moltype[molb.type].atoms.nres;
94     }
95     return nres;
96 }
97
98 AtomIterator::AtomIterator(const gmx_mtop_t& mtop, int globalAtomNumber) :
99     mtop_(&mtop),
100     mblock_(0),
101     atoms_(&mtop.moltype[mtop.molblock[0].type].atoms),
102     currentMolecule_(0),
103     highestResidueNumber_(mtop.maxResNumberNotRenumbered()),
104     localAtomNumber_(0),
105     globalAtomNumber_(globalAtomNumber)
106 {
107     GMX_ASSERT(globalAtomNumber == 0 || globalAtomNumber == mtop.natoms,
108                "Starting at other atoms not implemented yet");
109 }
110
111 AtomIterator& AtomIterator::operator++()
112 {
113     localAtomNumber_++;
114     globalAtomNumber_++;
115
116     if (localAtomNumber_ >= atoms_->nr)
117     {
118         if (atoms_->nres <= mtop_->maxResiduesPerMoleculeToTriggerRenumber())
119         {
120             /* Single residue molecule, increase the count with one */
121             highestResidueNumber_ += atoms_->nres;
122         }
123         currentMolecule_++;
124         localAtomNumber_ = 0;
125         if (currentMolecule_ >= mtop_->molblock[mblock_].nmol)
126         {
127             mblock_++;
128             if (mblock_ >= mtop_->molblock.size())
129             {
130                 return *this;
131             }
132             atoms_           = &mtop_->moltype[mtop_->molblock[mblock_].type].atoms;
133             currentMolecule_ = 0;
134         }
135     }
136     return *this;
137 }
138
139 bool AtomIterator::operator==(const AtomIterator& o) const
140 {
141     return mtop_ == o.mtop_ && globalAtomNumber_ == o.globalAtomNumber_;
142 }
143
144 const t_atom& AtomProxy::atom() const
145 {
146     return it_->atoms_->atom[it_->localAtomNumber_];
147 }
148
149 int AtomProxy::globalAtomNumber() const
150 {
151     return it_->globalAtomNumber_;
152 }
153
154 const char* AtomProxy::atomName() const
155 {
156     return *(it_->atoms_->atomname[it_->localAtomNumber_]);
157 }
158
159 const char* AtomProxy::residueName() const
160 {
161     int residueIndexInMolecule = it_->atoms_->atom[it_->localAtomNumber_].resind;
162     return *(it_->atoms_->resinfo[residueIndexInMolecule].name);
163 }
164
165 int AtomProxy::residueNumber() const
166 {
167     int residueIndexInMolecule = it_->atoms_->atom[it_->localAtomNumber_].resind;
168     if (it_->atoms_->nres <= it_->mtop_->maxResiduesPerMoleculeToTriggerRenumber())
169     {
170         return it_->highestResidueNumber_ + 1 + residueIndexInMolecule;
171     }
172     else
173     {
174         return it_->atoms_->resinfo[residueIndexInMolecule].nr;
175     }
176 }
177
178 const gmx_moltype_t& AtomProxy::moleculeType() const
179 {
180     return it_->mtop_->moltype[it_->mtop_->molblock[it_->mblock_].type];
181 }
182
183 int AtomProxy::atomNumberInMol() const
184 {
185     return it_->localAtomNumber_;
186 }
187
188 struct gmx_mtop_atomloop_block
189 {
190     const gmx_mtop_t* mtop;
191     size_t            mblock;
192     const t_atoms*    atoms;
193     int               at_local;
194 };
195
196 gmx_mtop_atomloop_block_t gmx_mtop_atomloop_block_init(const gmx_mtop_t& mtop)
197 {
198     struct gmx_mtop_atomloop_block* aloop = nullptr;
199
200     snew(aloop, 1);
201
202     aloop->mtop     = &mtop;
203     aloop->mblock   = 0;
204     aloop->atoms    = &mtop.moltype[mtop.molblock[aloop->mblock].type].atoms;
205     aloop->at_local = -1;
206
207     return aloop;
208 }
209
210 static void gmx_mtop_atomloop_block_destroy(gmx_mtop_atomloop_block_t aloop)
211 {
212     sfree(aloop);
213 }
214
215 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop, const t_atom** atom, int* nmol)
216 {
217     if (aloop == nullptr)
218     {
219         gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
220     }
221
222     aloop->at_local++;
223
224     if (aloop->at_local >= aloop->atoms->nr)
225     {
226         aloop->mblock++;
227         if (aloop->mblock >= aloop->mtop->molblock.size())
228         {
229             gmx_mtop_atomloop_block_destroy(aloop);
230             return FALSE;
231         }
232         aloop->atoms    = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
233         aloop->at_local = 0;
234     }
235
236     *atom = &aloop->atoms->atom[aloop->at_local];
237     *nmol = aloop->mtop->molblock[aloop->mblock].nmol;
238
239     return TRUE;
240 }
241
242 IListIterator::IListIterator(const gmx_mtop_t& mtop, size_t molblock) :
243     mtop_(&mtop),
244     mblock_(molblock)
245 {
246 }
247
248 IListIterator& IListIterator::operator++()
249 {
250     mblock_++;
251     return *this;
252 }
253
254 bool IListIterator::operator==(const IListIterator& o) const
255 {
256     return mtop_ == o.mtop_ && mblock_ == o.mblock_;
257 }
258
259 const InteractionLists& IListProxy::list() const
260 {
261     // one past the end means we want to take the
262     // intermolecular list instead.
263     if (it_->mblock_ == it_->mtop_->molblock.size())
264     {
265         return *it_->mtop_->intermolecular_ilist;
266     }
267     else
268     {
269         return it_->mtop_->moltype[it_->mtop_->molblock[it_->mblock_].type].ilist;
270     }
271 }
272
273 int IListProxy::nmol() const
274 {
275     // one past the end means we want to take the
276     // intermolecular list instead.
277     if (it_->mblock_ == it_->mtop_->molblock.size())
278     {
279         return 1;
280     }
281     else
282     {
283         return it_->mtop_->molblock[it_->mblock_].nmol;
284     }
285 }
286
287 IListRange::IListRange(const gmx_mtop_t& mtop) : begin_(mtop), end_(mtop, mtop.molblock.size())
288 {
289     if (mtop.bIntermolecularInteractions)
290     {
291         end_ = IListIterator(mtop, mtop.molblock.size() + 1);
292     }
293 }
294
295 int gmx_mtop_ftype_count(const gmx_mtop_t& mtop, int ftype)
296 {
297     int n = 0;
298
299     for (const IListProxy il : IListRange(mtop))
300     {
301         n += il.nmol() * il.list()[ftype].size() / (1 + NRAL(ftype));
302     }
303
304     if (mtop.bIntermolecularInteractions)
305     {
306         n += (*mtop.intermolecular_ilist)[ftype].size() / (1 + NRAL(ftype));
307     }
308
309     return n;
310 }
311
312 int gmx_mtop_interaction_count(const gmx_mtop_t& mtop, const int unsigned if_flags)
313 {
314     int n = 0;
315
316     for (const IListProxy il : IListRange(mtop))
317     {
318         for (int ftype = 0; ftype < F_NRE; ftype++)
319         {
320             if ((interaction_function[ftype].flags & if_flags) == if_flags)
321             {
322                 n += il.nmol() * il.list()[ftype].size() / (1 + NRAL(ftype));
323             }
324         }
325     }
326
327     if (mtop.bIntermolecularInteractions)
328     {
329         for (int ftype = 0; ftype < F_NRE; ftype++)
330         {
331             if ((interaction_function[ftype].flags & if_flags) == if_flags)
332             {
333                 n += (*mtop.intermolecular_ilist)[ftype].size() / (1 + NRAL(ftype));
334             }
335         }
336     }
337
338     return n;
339 }
340
341 gmx::EnumerationArray<ParticleType, int> gmx_mtop_particletype_count(const gmx_mtop_t& mtop)
342 {
343     gmx::EnumerationArray<ParticleType, int> count = { { 0 } };
344
345     for (const auto& molblock : mtop.molblock)
346     {
347         const t_atoms& atoms = mtop.moltype[molblock.type].atoms;
348         for (int a = 0; a < atoms.nr; a++)
349         {
350             count[atoms.atom[a].ptype] += molblock.nmol;
351         }
352     }
353
354     return count;
355 }
356
357 static void atomcat(t_atoms* dest, const t_atoms* src, int copies, int maxres_renum, int* maxresnr)
358 {
359     int i = 0, j = 0, l = 0, size = 0;
360     int srcnr  = src->nr;
361     int destnr = dest->nr;
362
363     if (dest->nr == 0)
364     {
365         dest->haveMass    = src->haveMass;
366         dest->haveType    = src->haveType;
367         dest->haveCharge  = src->haveCharge;
368         dest->haveBState  = src->haveBState;
369         dest->havePdbInfo = src->havePdbInfo;
370     }
371     else
372     {
373         dest->haveMass    = dest->haveMass && src->haveMass;
374         dest->haveType    = dest->haveType && src->haveType;
375         dest->haveCharge  = dest->haveCharge && src->haveCharge;
376         dest->haveBState  = dest->haveBState && src->haveBState;
377         dest->havePdbInfo = dest->havePdbInfo && src->havePdbInfo;
378     }
379
380     if (srcnr)
381     {
382         size = destnr + copies * srcnr;
383         srenew(dest->atom, size);
384         srenew(dest->atomname, size);
385         if (dest->haveType)
386         {
387             srenew(dest->atomtype, size);
388             if (dest->haveBState)
389             {
390                 srenew(dest->atomtypeB, size);
391             }
392         }
393         if (dest->havePdbInfo)
394         {
395             srenew(dest->pdbinfo, size);
396         }
397     }
398     if (src->nres)
399     {
400         size = dest->nres + copies * src->nres;
401         srenew(dest->resinfo, size);
402     }
403
404     /* residue information */
405     for (l = dest->nres, j = 0; (j < copies); j++, l += src->nres)
406     {
407         memcpy(reinterpret_cast<char*>(&(dest->resinfo[l])),
408                reinterpret_cast<char*>(&(src->resinfo[0])),
409                static_cast<size_t>(src->nres * sizeof(src->resinfo[0])));
410     }
411
412     for (l = destnr, j = 0; (j < copies); j++, l += srcnr)
413     {
414         memcpy(reinterpret_cast<char*>(&(dest->atom[l])),
415                reinterpret_cast<char*>(&(src->atom[0])),
416                static_cast<size_t>(srcnr * sizeof(src->atom[0])));
417         memcpy(reinterpret_cast<char*>(&(dest->atomname[l])),
418                reinterpret_cast<char*>(&(src->atomname[0])),
419                static_cast<size_t>(srcnr * sizeof(src->atomname[0])));
420         if (dest->haveType)
421         {
422             memcpy(reinterpret_cast<char*>(&(dest->atomtype[l])),
423                    reinterpret_cast<char*>(&(src->atomtype[0])),
424                    static_cast<size_t>(srcnr * sizeof(src->atomtype[0])));
425             if (dest->haveBState)
426             {
427                 memcpy(reinterpret_cast<char*>(&(dest->atomtypeB[l])),
428                        reinterpret_cast<char*>(&(src->atomtypeB[0])),
429                        static_cast<size_t>(srcnr * sizeof(src->atomtypeB[0])));
430             }
431         }
432         if (dest->havePdbInfo)
433         {
434             memcpy(reinterpret_cast<char*>(&(dest->pdbinfo[l])),
435                    reinterpret_cast<char*>(&(src->pdbinfo[0])),
436                    static_cast<size_t>(srcnr * sizeof(src->pdbinfo[0])));
437         }
438     }
439
440     /* Increment residue indices */
441     for (l = destnr, j = 0; (j < copies); j++)
442     {
443         for (i = 0; (i < srcnr); i++, l++)
444         {
445             dest->atom[l].resind = dest->nres + j * src->nres + src->atom[i].resind;
446         }
447     }
448
449     if (src->nres <= maxres_renum)
450     {
451         /* Single residue molecule, continue counting residues */
452         for (j = 0; (j < copies); j++)
453         {
454             for (l = 0; l < src->nres; l++)
455             {
456                 (*maxresnr)++;
457                 dest->resinfo[dest->nres + j * src->nres + l].nr = *maxresnr;
458             }
459         }
460     }
461
462     dest->nres += copies * src->nres;
463     dest->nr += copies * src->nr;
464 }
465
466 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t& mtop)
467 {
468     t_atoms atoms;
469
470     init_t_atoms(&atoms, 0, FALSE);
471
472     int maxresnr = mtop.maxResNumberNotRenumbered();
473     for (const gmx_molblock_t& molb : mtop.molblock)
474     {
475         atomcat(&atoms,
476                 &mtop.moltype[molb.type].atoms,
477                 molb.nmol,
478                 mtop.maxResiduesPerMoleculeToTriggerRenumber(),
479                 &maxresnr);
480     }
481
482     return atoms;
483 }
484
485 /*
486  * The cat routines below are old code from src/kernel/topcat.c
487  */
488
489 static void ilistcat(int ftype, InteractionList* dest, const InteractionList& src, int copies, int dnum, int snum)
490 {
491     const int nral = NRAL(ftype);
492
493     size_t destIndex = dest->iatoms.size();
494     dest->iatoms.resize(dest->iatoms.size() + copies * src.size());
495
496     for (int c = 0; c < copies; c++)
497     {
498         for (int i = 0; i < src.size();)
499         {
500             dest->iatoms[destIndex++] = src.iatoms[i++];
501             for (int a = 0; a < nral; a++)
502             {
503                 dest->iatoms[destIndex++] = dnum + src.iatoms[i++];
504             }
505         }
506         dnum += snum;
507     }
508 }
509
510 static void ilistcat(int ftype, t_ilist* dest, const InteractionList& src, int copies, int dnum, int snum)
511 {
512     const int nral = NRAL(ftype);
513
514     dest->nalloc = dest->nr + copies * src.size();
515     srenew(dest->iatoms, dest->nalloc);
516
517     for (int c = 0; c < copies; c++)
518     {
519         for (int i = 0; i < src.size();)
520         {
521             dest->iatoms[dest->nr++] = src.iatoms[i++];
522             for (int a = 0; a < nral; a++)
523             {
524                 dest->iatoms[dest->nr++] = dnum + src.iatoms[i++];
525             }
526         }
527         dnum += snum;
528     }
529 }
530
531 static const t_iparams& getIparams(const InteractionDefinitions& idef, const int index)
532 {
533     return idef.iparams[index];
534 }
535
536 static const t_iparams& getIparams(const t_idef& idef, const int index)
537 {
538     return idef.iparams[index];
539 }
540
541 static void resizeIParams(std::vector<t_iparams>* iparams, const int newSize)
542 {
543     iparams->resize(newSize);
544 }
545
546 static void resizeIParams(t_iparams** iparams, const int newSize)
547 {
548     srenew(*iparams, newSize);
549 }
550
551 template<typename IdefType>
552 static void set_posres_params(IdefType* idef, const gmx_molblock_t* molb, int i0, int a_offset)
553 {
554     auto* il = &idef->il[F_POSRES];
555     int   i1 = il->size() / 2;
556     resizeIParams(&idef->iparams_posres, i1);
557     for (int i = i0; i < i1; i++)
558     {
559         t_iparams* ip = &idef->iparams_posres[i];
560         /* Copy the force constants */
561         *ip        = getIparams(*idef, il->iatoms[i * 2]);
562         int a_molb = il->iatoms[i * 2 + 1] - a_offset;
563         if (molb->posres_xA.empty())
564         {
565             gmx_incons("Position restraint coordinates are missing");
566         }
567         ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
568         ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
569         ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
570         if (!molb->posres_xB.empty())
571         {
572             ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
573             ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
574             ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
575         }
576         else
577         {
578             ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
579             ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
580             ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
581         }
582         /* Set the parameter index for idef->iparams_posre */
583         il->iatoms[i * 2] = i;
584     }
585 }
586
587 template<typename IdefType>
588 static void set_fbposres_params(IdefType* idef, const gmx_molblock_t* molb, int i0, int a_offset)
589 {
590     auto* il = &idef->il[F_FBPOSRES];
591     int   i1 = il->size() / 2;
592     resizeIParams(&idef->iparams_fbposres, i1);
593     for (int i = i0; i < i1; i++)
594     {
595         t_iparams* ip = &idef->iparams_fbposres[i];
596         /* Copy the force constants */
597         *ip        = getIparams(*idef, il->iatoms[i * 2]);
598         int a_molb = il->iatoms[i * 2 + 1] - a_offset;
599         if (molb->posres_xA.empty())
600         {
601             gmx_incons("Position restraint coordinates are missing");
602         }
603         /* Take flat-bottom posres reference from normal position restraints */
604         ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
605         ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
606         ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
607         /* Note: no B-type for flat-bottom posres */
608
609         /* Set the parameter index for idef->iparams_posre */
610         il->iatoms[i * 2] = i;
611     }
612 }
613
614 /*! \brief Copy parameters to idef structure from mtop.
615  *
616  * Makes a deep copy of the force field parameters data structure from a gmx_mtop_t.
617  * Used to initialize legacy topology types.
618  *
619  * \param[in] mtop Reference to input mtop.
620  * \param[in] idef Pointer to idef to populate.
621  */
622 static void copyFFParametersFromMtop(const gmx_mtop_t& mtop, t_idef* idef)
623 {
624     const gmx_ffparams_t* ffp = &mtop.ffparams;
625
626     idef->ntypes = ffp->numTypes();
627     idef->atnr   = ffp->atnr;
628     /* we can no longer copy the pointers to the mtop members,
629      * because they will become invalid as soon as mtop gets free'd.
630      * We also need to make sure to only operate on valid data!
631      */
632
633     if (!ffp->functype.empty())
634     {
635         snew(idef->functype, ffp->functype.size());
636         std::copy(ffp->functype.data(), ffp->functype.data() + ffp->functype.size(), idef->functype);
637     }
638     else
639     {
640         idef->functype = nullptr;
641     }
642     if (!ffp->iparams.empty())
643     {
644         snew(idef->iparams, ffp->iparams.size());
645         std::copy(ffp->iparams.data(), ffp->iparams.data() + ffp->iparams.size(), idef->iparams);
646     }
647     else
648     {
649         idef->iparams = nullptr;
650     }
651     idef->iparams_posres   = nullptr;
652     idef->iparams_fbposres = nullptr;
653     idef->fudgeQQ          = ffp->fudgeQQ;
654     idef->ilsort           = ilsortUNKNOWN;
655 }
656
657 /*! \brief Copy idef structure from mtop.
658  *
659  * Makes a deep copy of an idef data structure from a gmx_mtop_t.
660  * Used to initialize legacy topology types.
661  *
662  * \param[in] mtop Reference to input mtop.
663  * \param[in] idef Pointer to idef to populate.
664  * \param[in] mergeConstr Decide if constraints will be merged.
665  */
666 template<typename IdefType>
667 static void copyIListsFromMtop(const gmx_mtop_t& mtop, IdefType* idef, bool mergeConstr)
668 {
669     int natoms = 0;
670     for (const gmx_molblock_t& molb : mtop.molblock)
671     {
672         const gmx_moltype_t& molt = mtop.moltype[molb.type];
673
674         int srcnr  = molt.atoms.nr;
675         int destnr = natoms;
676
677         int nposre_old   = idef->il[F_POSRES].size();
678         int nfbposre_old = idef->il[F_FBPOSRES].size();
679         for (int ftype = 0; ftype < F_NRE; ftype++)
680         {
681             if (mergeConstr && ftype == F_CONSTR && !molt.ilist[F_CONSTRNC].empty())
682             {
683                 /* Merge all constrains into one ilist.
684                  * This simplifies the constraint code.
685                  */
686                 for (int mol = 0; mol < molb.nmol; mol++)
687                 {
688                     ilistcat(ftype, &idef->il[F_CONSTR], molt.ilist[F_CONSTR], 1, destnr + mol * srcnr, srcnr);
689                     ilistcat(ftype, &idef->il[F_CONSTR], molt.ilist[F_CONSTRNC], 1, destnr + mol * srcnr, srcnr);
690                 }
691             }
692             else if (!(mergeConstr && ftype == F_CONSTRNC))
693             {
694                 ilistcat(ftype, &idef->il[ftype], molt.ilist[ftype], molb.nmol, destnr, srcnr);
695             }
696         }
697         if (idef->il[F_POSRES].size() > nposre_old)
698         {
699             /* Executing this line line stops gmxdump -sys working
700              * correctly. I'm not aware there's an elegant fix. */
701             set_posres_params(idef, &molb, nposre_old / 2, natoms);
702         }
703         if (idef->il[F_FBPOSRES].size() > nfbposre_old)
704         {
705             set_fbposres_params(idef, &molb, nfbposre_old / 2, natoms);
706         }
707
708         natoms += molb.nmol * srcnr;
709     }
710
711     if (mtop.bIntermolecularInteractions)
712     {
713         for (int ftype = 0; ftype < F_NRE; ftype++)
714         {
715             ilistcat(ftype, &idef->il[ftype], (*mtop.intermolecular_ilist)[ftype], 1, 0, mtop.natoms);
716         }
717     }
718
719     // We have not (yet) sorted free-energy interactions to the end of the ilists
720     idef->ilsort = ilsortNO_FE;
721 }
722
723 /*! \brief Copy atomtypes from mtop
724  *
725  * Makes a deep copy of t_atomtypes from gmx_mtop_t.
726  * Used to initialize legacy topology types.
727  *
728  * \param[in] mtop Reference to input mtop.
729  * \param[in] atomtypes Pointer to atomtypes to populate.
730  */
731 static void copyAtomtypesFromMtop(const gmx_mtop_t& mtop, t_atomtypes* atomtypes)
732 {
733     atomtypes->nr = mtop.atomtypes.nr;
734     if (mtop.atomtypes.atomnumber)
735     {
736         snew(atomtypes->atomnumber, mtop.atomtypes.nr);
737         std::copy(mtop.atomtypes.atomnumber,
738                   mtop.atomtypes.atomnumber + mtop.atomtypes.nr,
739                   atomtypes->atomnumber);
740     }
741     else
742     {
743         atomtypes->atomnumber = nullptr;
744     }
745 }
746
747 /*! \brief Generate a single list of lists of exclusions for the whole system
748  *
749  * \param[in] mtop  Reference to input mtop.
750  */
751 static gmx::ListOfLists<int> globalExclusionLists(const gmx_mtop_t& mtop)
752 {
753     gmx::ListOfLists<int> excls;
754
755     int atomIndex = 0;
756     for (const gmx_molblock_t& molb : mtop.molblock)
757     {
758         const gmx_moltype_t& molt = mtop.moltype[molb.type];
759
760         for (int mol = 0; mol < molb.nmol; mol++)
761         {
762             excls.appendListOfLists(molt.excls, atomIndex);
763
764             atomIndex += molt.atoms.nr;
765         }
766     }
767
768     return excls;
769 }
770
771 /*! \brief Updates inter-molecular exclusion lists
772  *
773  * This function updates inter-molecular exclusions to exclude all
774  * non-bonded interactions between a given list of atoms
775  *
776  * \param[inout]    excls   existing exclusions in local topology
777  * \param[in]       ids     list of global IDs of atoms
778  */
779 static void addMimicExclusions(gmx::ListOfLists<int>* excls, const gmx::ArrayRef<const int> ids)
780 {
781     t_blocka inter_excl{};
782     init_blocka(&inter_excl);
783     size_t n_q = ids.size();
784
785     inter_excl.nr  = excls->ssize();
786     inter_excl.nra = n_q * n_q;
787
788     size_t total_nra = n_q * n_q;
789
790     snew(inter_excl.index, excls->ssize() + 1);
791     snew(inter_excl.a, total_nra);
792
793     for (int i = 0; i < inter_excl.nr; ++i)
794     {
795         inter_excl.index[i] = 0;
796     }
797
798     /* Here we loop over the list of QM atom ids
799      *  and create exclusions between all of them resulting in
800      *  n_q * n_q sized exclusion list
801      */
802     int prev_index = 0;
803     for (int k = 0; k < inter_excl.nr; ++k)
804     {
805         inter_excl.index[k] = prev_index;
806         for (long i = 0; i < ids.ssize(); i++)
807         {
808             if (k != ids[i])
809             {
810                 continue;
811             }
812             size_t index             = n_q * i;
813             inter_excl.index[ids[i]] = index;
814             prev_index               = index + n_q;
815             for (size_t j = 0; j < n_q; ++j)
816             {
817                 inter_excl.a[n_q * i + j] = ids[j];
818             }
819         }
820     }
821     inter_excl.index[ids[n_q - 1] + 1] = n_q * n_q;
822
823     inter_excl.index[inter_excl.nr] = n_q * n_q;
824
825     std::vector<gmx::ExclusionBlock> qmexcl2(excls->size());
826     gmx::blockaToExclusionBlocks(&inter_excl, qmexcl2);
827
828     // Merge the created exclusion list with the existing one
829     gmx::mergeExclusions(excls, qmexcl2);
830 }
831
832 static void sortFreeEnergyInteractionsAtEnd(const gmx_mtop_t& mtop, InteractionDefinitions* idef)
833 {
834     std::vector<real> qA(mtop.natoms);
835     std::vector<real> qB(mtop.natoms);
836     for (const AtomProxy atomP : AtomRange(mtop))
837     {
838         const t_atom& local = atomP.atom();
839         int           index = atomP.globalAtomNumber();
840         qA[index]           = local.q;
841         qB[index]           = local.qB;
842     }
843     gmx_sort_ilist_fe(idef, qA.data(), qB.data());
844 }
845
846 static void gen_local_top(const gmx_mtop_t& mtop,
847                           bool              freeEnergyInteractionsAtEnd,
848                           bool              bMergeConstr,
849                           gmx_localtop_t*   top)
850 {
851     copyIListsFromMtop(mtop, &top->idef, bMergeConstr);
852     if (freeEnergyInteractionsAtEnd)
853     {
854         sortFreeEnergyInteractionsAtEnd(mtop, &top->idef);
855     }
856     top->excls = globalExclusionLists(mtop);
857     if (!mtop.intermolecularExclusionGroup.empty())
858     {
859         addMimicExclusions(&top->excls, mtop.intermolecularExclusionGroup);
860     }
861 }
862
863 void gmx_mtop_generate_local_top(const gmx_mtop_t& mtop, gmx_localtop_t* top, bool freeEnergyInteractionsAtEnd)
864 {
865     gen_local_top(mtop, freeEnergyInteractionsAtEnd, true, top);
866 }
867
868 /*! \brief Fills an array with molecule begin/end atom indices
869  *
870  * \param[in]  mtop   The global topology
871  * \param[out] index  Array of size nr. of molecules + 1 to be filled with molecule begin/end indices
872  */
873 static void fillMoleculeIndices(const gmx_mtop_t& mtop, gmx::ArrayRef<int> index)
874 {
875     int globalAtomIndex   = 0;
876     int globalMolIndex    = 0;
877     index[globalMolIndex] = globalAtomIndex;
878     for (const gmx_molblock_t& molb : mtop.molblock)
879     {
880         int numAtomsPerMolecule = mtop.moltype[molb.type].atoms.nr;
881         for (int mol = 0; mol < molb.nmol; mol++)
882         {
883             globalAtomIndex += numAtomsPerMolecule;
884             globalMolIndex += 1;
885             index[globalMolIndex] = globalAtomIndex;
886         }
887     }
888 }
889
890 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t& mtop)
891 {
892     gmx::RangePartitioning mols;
893
894     for (const gmx_molblock_t& molb : mtop.molblock)
895     {
896         int numAtomsPerMolecule = mtop.moltype[molb.type].atoms.nr;
897         for (int mol = 0; mol < molb.nmol; mol++)
898         {
899             mols.appendBlock(numAtomsPerMolecule);
900         }
901     }
902
903     return mols;
904 }
905
906 std::vector<gmx::Range<int>> atomRangeOfEachResidue(const gmx_moltype_t& moltype)
907 {
908     std::vector<gmx::Range<int>> atomRanges;
909     int                          currentResidueNumber = moltype.atoms.atom[0].resind;
910     int                          startAtom            = 0;
911     // Go through all atoms in a molecule to store first and last atoms in each residue.
912     for (int i = 0; i < moltype.atoms.nr; i++)
913     {
914         int residueOfThisAtom = moltype.atoms.atom[i].resind;
915         if (residueOfThisAtom != currentResidueNumber)
916         {
917             // This atom belongs to the next residue, so record the range for the previous residue,
918             // remembering that end points to one place past the last atom.
919             int endAtom = i;
920             atomRanges.emplace_back(startAtom, endAtom);
921             // Prepare for the current residue
922             startAtom            = endAtom;
923             currentResidueNumber = residueOfThisAtom;
924         }
925     }
926     // special treatment for last residue in this molecule.
927     atomRanges.emplace_back(startAtom, moltype.atoms.nr);
928
929     return atomRanges;
930 }
931
932 /*! \brief Creates and returns a deprecated t_block struct with molecule indices
933  *
934  * \param[in] mtop  The global topology
935  */
936 static t_block gmx_mtop_molecules_t_block(const gmx_mtop_t& mtop)
937 {
938     t_block mols;
939
940     mols.nr           = gmx_mtop_num_molecules(mtop);
941     mols.nalloc_index = mols.nr + 1;
942     snew(mols.index, mols.nalloc_index);
943
944     fillMoleculeIndices(mtop, gmx::arrayRefFromArray(mols.index, mols.nr + 1));
945
946     return mols;
947 }
948
949 static void gen_t_topology(const gmx_mtop_t& mtop, bool bMergeConstr, t_topology* top)
950 {
951     copyAtomtypesFromMtop(mtop, &top->atomtypes);
952     for (int ftype = 0; ftype < F_NRE; ftype++)
953     {
954         top->idef.il[ftype].nr     = 0;
955         top->idef.il[ftype].nalloc = 0;
956         top->idef.il[ftype].iatoms = nullptr;
957     }
958     copyFFParametersFromMtop(mtop, &top->idef);
959     copyIListsFromMtop(mtop, &top->idef, bMergeConstr);
960
961     top->name                        = mtop.name;
962     top->atoms                       = gmx_mtop_global_atoms(mtop);
963     top->mols                        = gmx_mtop_molecules_t_block(mtop);
964     top->bIntermolecularInteractions = mtop.bIntermolecularInteractions;
965     top->symtab                      = mtop.symtab;
966 }
967
968 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t* mtop, bool freeMTop)
969 {
970     t_topology top;
971
972     gen_t_topology(*mtop, false, &top);
973
974     if (freeMTop)
975     {
976         // Clear pointers and counts, such that the pointers copied to top
977         // keep pointing to valid data after destroying mtop.
978         mtop->symtab.symbuf = nullptr;
979         mtop->symtab.nr     = 0;
980     }
981     return top;
982 }
983
984 std::vector<int> get_atom_index(const gmx_mtop_t& mtop)
985 {
986
987     std::vector<int> atom_index;
988     for (const AtomProxy atomP : AtomRange(mtop))
989     {
990         const t_atom& local = atomP.atom();
991         int           index = atomP.globalAtomNumber();
992         if (local.ptype == ParticleType::Atom)
993         {
994             atom_index.push_back(index);
995         }
996     }
997     return atom_index;
998 }
999
1000 void convertAtomsToMtop(t_symtab* symtab, char** name, t_atoms* atoms, gmx_mtop_t* mtop)
1001 {
1002     mtop->symtab = *symtab;
1003
1004     mtop->name = name;
1005
1006     mtop->moltype.clear();
1007     mtop->moltype.resize(1);
1008     mtop->moltype.back().atoms = *atoms;
1009
1010     mtop->molblock.resize(1);
1011     mtop->molblock[0].type = 0;
1012     mtop->molblock[0].nmol = 1;
1013
1014     mtop->bIntermolecularInteractions = FALSE;
1015
1016     mtop->natoms = atoms->nr;
1017
1018     mtop->haveMoleculeIndices = false;
1019
1020     mtop->finalize();
1021 }
1022
1023 bool haveFepPerturbedNBInteractions(const gmx_mtop_t& mtop)
1024 {
1025     for (const gmx_moltype_t& molt : mtop.moltype)
1026     {
1027         for (int a = 0; a < molt.atoms.nr; a++)
1028         {
1029             if (PERTURBED(molt.atoms.atom[a]))
1030             {
1031                 return true;
1032             }
1033         }
1034     }
1035
1036     return false;
1037 }
1038
1039 bool haveFepPerturbedMasses(const gmx_mtop_t& mtop)
1040 {
1041     for (const gmx_moltype_t& molt : mtop.moltype)
1042     {
1043         for (int a = 0; a < molt.atoms.nr; a++)
1044         {
1045             const t_atom& atom = molt.atoms.atom[a];
1046             if (atom.m != atom.mB)
1047             {
1048                 return true;
1049             }
1050         }
1051     }
1052
1053     return false;
1054 }
1055
1056 bool havePerturbedConstraints(const gmx_mtop_t& mtop)
1057 {
1058     // This code assumes that all perturbed constraints parameters are actually used
1059     const auto& ffparams = mtop.ffparams;
1060
1061     for (gmx::index i = 0; i < gmx::ssize(ffparams.functype); i++)
1062     {
1063         if (ffparams.functype[i] == F_CONSTR || ffparams.functype[i] == F_CONSTRNC)
1064         {
1065             const auto& iparams = ffparams.iparams[i];
1066             if (iparams.constr.dA != iparams.constr.dB)
1067             {
1068                 return true;
1069             }
1070         }
1071     }
1072
1073     return false;
1074 }