2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2008,2009,2010,2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "mtop_util.h"
45 #include "gromacs/legacyheaders/types/enums.h"
46 #include "gromacs/legacyheaders/types/ifunc.h"
47 #include "gromacs/legacyheaders/types/inputrec.h"
48 #include "gromacs/legacyheaders/types/simple.h"
49 #include "gromacs/math/vectypes.h"
50 #include "gromacs/topology/atoms.h"
51 #include "gromacs/topology/block.h"
52 #include "gromacs/topology/idef.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/topology/topsort.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/real.h"
57 #include "gromacs/utility/smalloc.h"
59 static int gmx_mtop_maxresnr(const gmx_mtop_t *mtop, int maxres_renum)
66 for (mt = 0; mt < mtop->nmoltype; mt++)
68 atoms = &mtop->moltype[mt].atoms;
69 if (atoms->nres > maxres_renum)
71 for (r = 0; r < atoms->nres; r++)
73 if (atoms->resinfo[r].nr > maxresnr)
75 maxresnr = atoms->resinfo[r].nr;
84 void gmx_mtop_finalize(gmx_mtop_t *mtop)
88 mtop->maxres_renum = 1;
90 env = getenv("GMX_MAXRESRENUM");
93 sscanf(env, "%d", &mtop->maxres_renum);
95 if (mtop->maxres_renum == -1)
97 /* -1 signals renumber residues in all molecules */
98 mtop->maxres_renum = INT_MAX;
101 mtop->maxresnr = gmx_mtop_maxresnr(mtop, mtop->maxres_renum);
104 void gmx_mtop_count_atomtypes(const gmx_mtop_t *mtop, int state, int typecount[])
106 int i, mb, nmol, tpi;
109 for (i = 0; i < mtop->ffparams.atnr; ++i)
113 for (mb = 0; mb < mtop->nmolblock; ++mb)
115 nmol = mtop->molblock[mb].nmol;
116 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
117 for (i = 0; i < atoms->nr; ++i)
121 tpi = atoms->atom[i].type;
125 tpi = atoms->atom[i].typeB;
127 typecount[tpi] += nmol;
132 int ncg_mtop(const gmx_mtop_t *mtop)
138 for (mb = 0; mb < mtop->nmolblock; mb++)
141 mtop->molblock[mb].nmol*
142 mtop->moltype[mtop->molblock[mb].type].cgs.nr;
148 void gmx_mtop_remove_chargegroups(gmx_mtop_t *mtop)
154 for (mt = 0; mt < mtop->nmoltype; mt++)
156 cgs = &mtop->moltype[mt].cgs;
157 if (cgs->nr < mtop->moltype[mt].atoms.nr)
159 cgs->nr = mtop->moltype[mt].atoms.nr;
160 srenew(cgs->index, cgs->nr+1);
161 for (i = 0; i < cgs->nr+1; i++)
177 typedef struct gmx_mtop_atomlookup
179 const gmx_mtop_t *mtop;
183 } t_gmx_mtop_atomlookup;
186 gmx_mtop_atomlookup_t
187 gmx_mtop_atomlookup_init(const gmx_mtop_t *mtop)
189 t_gmx_mtop_atomlookup *alook;
191 int a_start, a_end, na, na_start = -1;
196 alook->nmb = mtop->nmolblock;
198 snew(alook->mba, alook->nmb);
201 for (mb = 0; mb < mtop->nmolblock; mb++)
203 na = mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
204 a_end = a_start + na;
206 alook->mba[mb].a_start = a_start;
207 alook->mba[mb].a_end = a_end;
208 alook->mba[mb].na_mol = mtop->molblock[mb].natoms_mol;
210 /* We start the binary search with the largest block */
211 if (mb == 0 || na > na_start)
213 alook->mb_start = mb;
223 gmx_mtop_atomlookup_t
224 gmx_mtop_atomlookup_settle_init(const gmx_mtop_t *mtop)
226 t_gmx_mtop_atomlookup *alook;
228 int na, na_start = -1;
230 alook = gmx_mtop_atomlookup_init(mtop);
232 /* Check if the starting molblock has settle */
233 if (mtop->moltype[mtop->molblock[alook->mb_start].type].ilist[F_SETTLE].nr == 0)
235 /* Search the largest molblock with settle */
236 alook->mb_start = -1;
237 for (mb = 0; mb < mtop->nmolblock; mb++)
239 if (mtop->moltype[mtop->molblock[mb].type].ilist[F_SETTLE].nr > 0)
241 na = alook->mba[mb].a_end - alook->mba[mb].a_start;
242 if (alook->mb_start == -1 || na > na_start)
244 alook->mb_start = mb;
250 if (alook->mb_start == -1)
252 gmx_incons("gmx_mtop_atomlookup_settle_init called without settles");
260 gmx_mtop_atomlookup_destroy(gmx_mtop_atomlookup_t alook)
266 void gmx_mtop_atomnr_to_atom(const gmx_mtop_atomlookup_t alook,
271 int a_start, atnr_mol;
274 if (atnr_global < 0 || atnr_global >= mtop->natoms)
276 gmx_fatal(FARGS, "gmx_mtop_atomnr_to_moltype was called with atnr_global=%d which is not in the atom range of this system (%d-%d)",
277 atnr_global, 0, mtop->natoms-1);
283 mb = alook->mb_start;
287 a_start = alook->mba[mb].a_start;
288 if (atnr_global < a_start)
292 else if (atnr_global >= alook->mba[mb].a_end)
300 mb = ((mb0 + mb1 + 1)>>1);
303 atnr_mol = (atnr_global - a_start) % alook->mba[mb].na_mol;
305 *atom = &alook->mtop->moltype[alook->mtop->molblock[mb].type].atoms.atom[atnr_mol];
308 void gmx_mtop_atomnr_to_ilist(const gmx_mtop_atomlookup_t alook,
310 t_ilist **ilist_mol, int *atnr_offset)
313 int a_start, atnr_local;
316 if (atnr_global < 0 || atnr_global >= mtop->natoms)
318 gmx_fatal(FARGS, "gmx_mtop_atomnr_to_moltype was called with atnr_global=%d which is not in the atom range of this system (%d-%d)",
319 atnr_global, 0, mtop->natoms-1);
325 mb = alook->mb_start;
329 a_start = alook->mba[mb].a_start;
330 if (atnr_global < a_start)
334 else if (atnr_global >= alook->mba[mb].a_end)
342 mb = ((mb0 + mb1 + 1)>>1);
345 *ilist_mol = alook->mtop->moltype[alook->mtop->molblock[mb].type].ilist;
347 atnr_local = (atnr_global - a_start) % alook->mba[mb].na_mol;
349 *atnr_offset = atnr_global - atnr_local;
352 void gmx_mtop_atomnr_to_molblock_ind(const gmx_mtop_atomlookup_t alook,
354 int *molb, int *molnr, int *atnr_mol)
360 if (atnr_global < 0 || atnr_global >= mtop->natoms)
362 gmx_fatal(FARGS, "gmx_mtop_atomnr_to_moltype was called with atnr_global=%d which is not in the atom range of this system (%d-%d)",
363 atnr_global, 0, mtop->natoms-1);
369 mb = alook->mb_start;
373 a_start = alook->mba[mb].a_start;
374 if (atnr_global < a_start)
378 else if (atnr_global >= alook->mba[mb].a_end)
386 mb = ((mb0 + mb1 + 1)>>1);
390 *molnr = (atnr_global - a_start) / alook->mba[mb].na_mol;
391 *atnr_mol = atnr_global - a_start - (*molnr)*alook->mba[mb].na_mol;
394 void gmx_mtop_atominfo_global(const gmx_mtop_t *mtop, int atnr_global,
395 char **atomname, int *resnr, char **resname)
397 int mb, a_start, a_end, maxresnr, at_loc;
398 gmx_molblock_t *molb;
399 t_atoms *atoms = NULL;
401 if (atnr_global < 0 || atnr_global >= mtop->natoms)
403 gmx_fatal(FARGS, "gmx_mtop_atominfo_global was called with atnr_global=%d which is not in the atom range of this system (%d-%d)",
404 atnr_global, 0, mtop->natoms-1);
409 maxresnr = mtop->maxresnr;
414 if (atoms->nres <= mtop->maxres_renum)
416 /* Single residue molecule, keep counting */
417 maxresnr += mtop->molblock[mb].nmol*atoms->nres;
421 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
423 a_end = a_start + mtop->molblock[mb].nmol*atoms->nr;
425 while (atnr_global >= a_end);
427 at_loc = (atnr_global - a_start) % atoms->nr;
428 *atomname = *(atoms->atomname[at_loc]);
429 if (atoms->nres > mtop->maxres_renum)
431 *resnr = atoms->resinfo[atoms->atom[at_loc].resind].nr;
435 /* Single residue molecule, keep counting */
436 *resnr = maxresnr + 1 + (atnr_global - a_start)/atoms->nr*atoms->nres + atoms->atom[at_loc].resind;
438 *resname = *(atoms->resinfo[atoms->atom[at_loc].resind].name);
441 typedef struct gmx_mtop_atomloop_all
443 const gmx_mtop_t *mtop;
450 } t_gmx_mtop_atomloop_all;
452 gmx_mtop_atomloop_all_t
453 gmx_mtop_atomloop_all_init(const gmx_mtop_t *mtop)
455 struct gmx_mtop_atomloop_all *aloop;
462 &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
464 aloop->maxresnr = mtop->maxresnr;
465 aloop->at_local = -1;
466 aloop->at_global = -1;
471 static void gmx_mtop_atomloop_all_destroy(gmx_mtop_atomloop_all_t aloop)
476 gmx_bool gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop,
477 int *at_global, t_atom **atom)
481 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
487 if (aloop->at_local >= aloop->atoms->nr)
489 if (aloop->atoms->nres <= aloop->mtop->maxres_renum)
491 /* Single residue molecule, increase the count with one */
492 aloop->maxresnr += aloop->atoms->nres;
496 if (aloop->mol >= aloop->mtop->molblock[aloop->mblock].nmol)
499 if (aloop->mblock >= aloop->mtop->nmolblock)
501 gmx_mtop_atomloop_all_destroy(aloop);
504 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
509 *at_global = aloop->at_global;
510 *atom = &aloop->atoms->atom[aloop->at_local];
515 void gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop,
516 char **atomname, int *resnr, char **resname)
520 *atomname = *(aloop->atoms->atomname[aloop->at_local]);
521 resind_mol = aloop->atoms->atom[aloop->at_local].resind;
522 *resnr = aloop->atoms->resinfo[resind_mol].nr;
523 if (aloop->atoms->nres <= aloop->mtop->maxres_renum)
525 *resnr = aloop->maxresnr + 1 + resind_mol;
527 *resname = *(aloop->atoms->resinfo[resind_mol].name);
530 void gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop,
531 gmx_moltype_t **moltype, int *at_mol)
533 *moltype = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type];
534 *at_mol = aloop->at_local;
537 typedef struct gmx_mtop_atomloop_block
539 const gmx_mtop_t *mtop;
543 } t_gmx_mtop_atomloop_block;
545 gmx_mtop_atomloop_block_t
546 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop)
548 struct gmx_mtop_atomloop_block *aloop;
554 aloop->atoms = &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
555 aloop->at_local = -1;
560 static void gmx_mtop_atomloop_block_destroy(gmx_mtop_atomloop_block_t aloop)
565 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
566 t_atom **atom, int *nmol)
570 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
575 if (aloop->at_local >= aloop->atoms->nr)
578 if (aloop->mblock >= aloop->mtop->nmolblock)
580 gmx_mtop_atomloop_block_destroy(aloop);
583 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
587 *atom = &aloop->atoms->atom[aloop->at_local];
588 *nmol = aloop->mtop->molblock[aloop->mblock].nmol;
593 typedef struct gmx_mtop_ilistloop
595 const gmx_mtop_t *mtop;
600 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop)
602 struct gmx_mtop_ilistloop *iloop;
612 static void gmx_mtop_ilistloop_destroy(gmx_mtop_ilistloop_t iloop)
617 gmx_bool gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
618 t_ilist **ilist_mol, int *nmol)
622 gmx_incons("gmx_mtop_ilistloop_next called without calling gmx_mtop_ilistloop_init");
626 if (iloop->mblock == iloop->mtop->nmolblock)
628 gmx_mtop_ilistloop_destroy(iloop);
633 iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
635 *nmol = iloop->mtop->molblock[iloop->mblock].nmol;
639 typedef struct gmx_mtop_ilistloop_all
641 const gmx_mtop_t *mtop;
645 } t_gmx_mtop_ilist_all;
647 gmx_mtop_ilistloop_all_t
648 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop)
650 struct gmx_mtop_ilistloop_all *iloop;
662 static void gmx_mtop_ilistloop_all_destroy(gmx_mtop_ilistloop_all_t iloop)
667 gmx_bool gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
668 t_ilist **ilist_mol, int *atnr_offset)
670 gmx_molblock_t *molb;
674 gmx_incons("gmx_mtop_ilistloop_all_next called without calling gmx_mtop_ilistloop_all_init");
679 iloop->a_offset += iloop->mtop->molblock[iloop->mblock].natoms_mol;
684 if (iloop->mol >= iloop->mtop->molblock[iloop->mblock].nmol)
688 if (iloop->mblock == iloop->mtop->nmolblock)
690 gmx_mtop_ilistloop_all_destroy(iloop);
696 iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
698 *atnr_offset = iloop->a_offset;
703 int gmx_mtop_ftype_count(const gmx_mtop_t *mtop, int ftype)
705 gmx_mtop_ilistloop_t iloop;
711 iloop = gmx_mtop_ilistloop_init(mtop);
712 while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
714 n += nmol*il[ftype].nr/(1+NRAL(ftype));
720 t_block gmx_mtop_global_cgs(const gmx_mtop_t *mtop)
722 t_block cgs_gl, *cgs_mol;
724 gmx_molblock_t *molb;
727 /* In most cases this is too much, but we realloc at the end */
728 snew(cgs_gl.index, mtop->natoms+1);
732 for (mb = 0; mb < mtop->nmolblock; mb++)
734 molb = &mtop->molblock[mb];
735 cgs_mol = &mtop->moltype[molb->type].cgs;
736 for (mol = 0; mol < molb->nmol; mol++)
738 for (cg = 0; cg < cgs_mol->nr; cg++)
740 cgs_gl.index[cgs_gl.nr+1] =
741 cgs_gl.index[cgs_gl.nr] +
742 cgs_mol->index[cg+1] - cgs_mol->index[cg];
747 cgs_gl.nalloc_index = cgs_gl.nr + 1;
748 srenew(cgs_gl.index, cgs_gl.nalloc_index);
753 static void atomcat(t_atoms *dest, t_atoms *src, int copies,
754 int maxres_renum, int *maxresnr)
758 int destnr = dest->nr;
762 size = destnr+copies*srcnr;
763 srenew(dest->atom, size);
764 srenew(dest->atomname, size);
765 srenew(dest->atomtype, size);
766 srenew(dest->atomtypeB, size);
770 size = dest->nres+copies*src->nres;
771 srenew(dest->resinfo, size);
774 /* residue information */
775 for (l = dest->nres, j = 0; (j < copies); j++, l += src->nres)
777 memcpy((char *) &(dest->resinfo[l]), (char *) &(src->resinfo[0]),
778 (size_t)(src->nres*sizeof(src->resinfo[0])));
781 for (l = destnr, j = 0; (j < copies); j++, l += srcnr)
783 memcpy((char *) &(dest->atomname[l]), (char *) &(src->atomname[0]),
784 (size_t)(srcnr*sizeof(src->atomname[0])));
785 memcpy((char *) &(dest->atomtype[l]), (char *) &(src->atomtype[0]),
786 (size_t)(srcnr*sizeof(src->atomtype[0])));
787 memcpy((char *) &(dest->atomtypeB[l]), (char *) &(src->atomtypeB[0]),
788 (size_t)(srcnr*sizeof(src->atomtypeB[0])));
789 memcpy((char *) &(dest->atom[l]), (char *) &(src->atom[0]),
790 (size_t)(srcnr*sizeof(src->atom[0])));
793 /* Increment residue indices */
794 for (l = destnr, j = 0; (j < copies); j++)
796 for (i = 0; (i < srcnr); i++, l++)
798 dest->atom[l].resind = dest->nres+j*src->nres+src->atom[i].resind;
802 if (src->nres <= maxres_renum)
804 /* Single residue molecule, continue counting residues */
805 for (j = 0; (j < copies); j++)
807 for (l = 0; l < src->nres; l++)
810 dest->resinfo[dest->nres+j*src->nres+l].nr = *maxresnr;
815 dest->nres += copies*src->nres;
816 dest->nr += copies*src->nr;
819 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t *mtop)
823 gmx_molblock_t *molb;
825 init_t_atoms(&atoms, 0, FALSE);
827 maxresnr = mtop->maxresnr;
828 for (mb = 0; mb < mtop->nmolblock; mb++)
830 molb = &mtop->molblock[mb];
831 atomcat(&atoms, &mtop->moltype[molb->type].atoms, molb->nmol,
832 mtop->maxres_renum, &maxresnr);
838 void gmx_mtop_make_atomic_charge_groups(gmx_mtop_t *mtop,
839 gmx_bool bKeepSingleMolCG)
844 for (mb = 0; mb < mtop->nmolblock; mb++)
846 cgs_mol = &mtop->moltype[mtop->molblock[mb].type].cgs;
847 if (!(bKeepSingleMolCG && cgs_mol->nr == 1))
849 cgs_mol->nr = mtop->molblock[mb].natoms_mol;
850 cgs_mol->nalloc_index = cgs_mol->nr + 1;
851 srenew(cgs_mol->index, cgs_mol->nalloc_index);
852 for (cg = 0; cg < cgs_mol->nr+1; cg++)
854 cgs_mol->index[cg] = cg;
861 * The cat routines below are old code from src/kernel/topcat.c
864 static void blockcat(t_block *dest, t_block *src, int copies)
866 int i, j, l, nra, size;
870 size = (dest->nr+copies*src->nr+1);
871 srenew(dest->index, size);
874 nra = dest->index[dest->nr];
875 for (l = dest->nr, j = 0; (j < copies); j++)
877 for (i = 0; (i < src->nr); i++)
879 dest->index[l++] = nra + src->index[i];
881 nra += src->index[src->nr];
883 dest->nr += copies*src->nr;
884 dest->index[dest->nr] = nra;
887 static void blockacat(t_blocka *dest, t_blocka *src, int copies,
891 int destnr = dest->nr;
892 int destnra = dest->nra;
896 size = (dest->nr+copies*src->nr+1);
897 srenew(dest->index, size);
901 size = (dest->nra+copies*src->nra);
902 srenew(dest->a, size);
905 for (l = destnr, j = 0; (j < copies); j++)
907 for (i = 0; (i < src->nr); i++)
909 dest->index[l++] = dest->nra+src->index[i];
911 dest->nra += src->nra;
913 for (l = destnra, j = 0; (j < copies); j++)
915 for (i = 0; (i < src->nra); i++)
917 dest->a[l++] = dnum+src->a[i];
922 dest->index[dest->nr] = dest->nra;
925 static void ilistcat(int ftype, t_ilist *dest, t_ilist *src, int copies,
932 dest->nalloc = dest->nr + copies*src->nr;
933 srenew(dest->iatoms, dest->nalloc);
935 for (c = 0; c < copies; c++)
937 for (i = 0; i < src->nr; )
939 dest->iatoms[dest->nr++] = src->iatoms[i++];
940 for (a = 0; a < nral; a++)
942 dest->iatoms[dest->nr++] = dnum + src->iatoms[i++];
949 static void set_posres_params(t_idef *idef, gmx_molblock_t *molb,
950 int i0, int a_offset)
956 il = &idef->il[F_POSRES];
958 idef->iparams_posres_nalloc = i1;
959 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
960 for (i = i0; i < i1; i++)
962 ip = &idef->iparams_posres[i];
963 /* Copy the force constants */
964 *ip = idef->iparams[il->iatoms[i*2]];
965 a_molb = il->iatoms[i*2+1] - a_offset;
966 if (molb->nposres_xA == 0)
968 gmx_incons("Position restraint coordinates are missing");
970 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
971 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
972 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
973 if (molb->nposres_xB > 0)
975 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
976 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
977 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
981 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
982 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
983 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
985 /* Set the parameter index for idef->iparams_posre */
990 static void set_fbposres_params(t_idef *idef, gmx_molblock_t *molb,
991 int i0, int a_offset)
997 il = &idef->il[F_FBPOSRES];
999 idef->iparams_fbposres_nalloc = i1;
1000 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
1001 for (i = i0; i < i1; i++)
1003 ip = &idef->iparams_fbposres[i];
1004 /* Copy the force constants */
1005 *ip = idef->iparams[il->iatoms[i*2]];
1006 a_molb = il->iatoms[i*2+1] - a_offset;
1007 if (molb->nposres_xA == 0)
1009 gmx_incons("Position restraint coordinates are missing");
1011 /* Take flat-bottom posres reference from normal position restraints */
1012 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
1013 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
1014 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
1015 /* Note: no B-type for flat-bottom posres */
1017 /* Set the parameter index for idef->iparams_posre */
1018 il->iatoms[i*2] = i;
1022 static void gen_local_top(const gmx_mtop_t *mtop, const t_inputrec *ir,
1023 gmx_bool bMergeConstr,
1024 gmx_localtop_t *top)
1026 int mb, srcnr, destnr, ftype, ftype_dest, mt, natoms, mol, nposre_old, nfbposre_old;
1027 gmx_molblock_t *molb;
1028 gmx_moltype_t *molt;
1029 const gmx_ffparams_t *ffp;
1032 gmx_mtop_atomloop_all_t aloop;
1036 top->atomtypes = mtop->atomtypes;
1038 ffp = &mtop->ffparams;
1041 idef->ntypes = ffp->ntypes;
1042 idef->atnr = ffp->atnr;
1043 idef->functype = ffp->functype;
1044 idef->iparams = ffp->iparams;
1045 idef->iparams_posres = NULL;
1046 idef->iparams_posres_nalloc = 0;
1047 idef->iparams_fbposres = NULL;
1048 idef->iparams_fbposres_nalloc = 0;
1049 idef->fudgeQQ = ffp->fudgeQQ;
1050 idef->cmap_grid = ffp->cmap_grid;
1051 idef->ilsort = ilsortUNKNOWN;
1053 init_block(&top->cgs);
1054 init_blocka(&top->excls);
1055 for (ftype = 0; ftype < F_NRE; ftype++)
1057 idef->il[ftype].nr = 0;
1058 idef->il[ftype].nalloc = 0;
1059 idef->il[ftype].iatoms = NULL;
1063 for (mb = 0; mb < mtop->nmolblock; mb++)
1065 molb = &mtop->molblock[mb];
1066 molt = &mtop->moltype[molb->type];
1068 srcnr = molt->atoms.nr;
1071 blockcat(&top->cgs, &molt->cgs, molb->nmol);
1073 blockacat(&top->excls, &molt->excls, molb->nmol, destnr, srcnr);
1075 nposre_old = idef->il[F_POSRES].nr;
1076 nfbposre_old = idef->il[F_FBPOSRES].nr;
1077 for (ftype = 0; ftype < F_NRE; ftype++)
1080 ftype == F_CONSTR && molt->ilist[F_CONSTRNC].nr > 0)
1082 /* Merge all constrains into one ilist.
1083 * This simplifies the constraint code.
1085 for (mol = 0; mol < molb->nmol; mol++)
1087 ilistcat(ftype, &idef->il[F_CONSTR], &molt->ilist[F_CONSTR],
1088 1, destnr+mol*srcnr, srcnr);
1089 ilistcat(ftype, &idef->il[F_CONSTR], &molt->ilist[F_CONSTRNC],
1090 1, destnr+mol*srcnr, srcnr);
1093 else if (!(bMergeConstr && ftype == F_CONSTRNC))
1095 ilistcat(ftype, &idef->il[ftype], &molt->ilist[ftype],
1096 molb->nmol, destnr, srcnr);
1099 if (idef->il[F_POSRES].nr > nposre_old)
1101 /* Executing this line line stops gmxdump -sys working
1102 * correctly. I'm not aware there's an elegant fix. */
1103 set_posres_params(idef, molb, nposre_old/2, natoms);
1105 if (idef->il[F_FBPOSRES].nr > nfbposre_old)
1107 set_fbposres_params(idef, molb, nfbposre_old/2, natoms);
1110 natoms += molb->nmol*srcnr;
1115 top->idef.ilsort = ilsortUNKNOWN;
1119 if (ir->efep != efepNO && gmx_mtop_bondeds_free_energy(mtop))
1121 snew(qA, mtop->natoms);
1122 snew(qB, mtop->natoms);
1123 aloop = gmx_mtop_atomloop_all_init(mtop);
1124 while (gmx_mtop_atomloop_all_next(aloop, &ag, &atom))
1129 gmx_sort_ilist_fe(&top->idef, qA, qB);
1135 top->idef.ilsort = ilsortNO_FE;
1140 gmx_localtop_t *gmx_mtop_generate_local_top(const gmx_mtop_t *mtop,
1141 const t_inputrec *ir)
1143 gmx_localtop_t *top;
1147 gen_local_top(mtop, ir, TRUE, top);
1152 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop)
1155 gmx_localtop_t ltop;
1158 gen_local_top(mtop, NULL, FALSE, <op);
1160 top.name = mtop->name;
1161 top.idef = ltop.idef;
1162 top.atomtypes = ltop.atomtypes;
1164 top.excls = ltop.excls;
1165 top.atoms = gmx_mtop_global_atoms(mtop);
1166 top.mols = mtop->mols;
1167 top.symtab = mtop->symtab;
1169 /* We only need to free the moltype and molblock data,
1170 * all other pointers have been copied to top.
1172 * Well, except for the group data, but we can't free those, because they
1173 * are used somewhere even after a call to this function.
1175 for (mt = 0; mt < mtop->nmoltype; mt++)
1177 done_moltype(&mtop->moltype[mt]);
1179 sfree(mtop->moltype);
1181 for (mb = 0; mb < mtop->nmolblock; mb++)
1183 done_molblock(&mtop->molblock[mb]);
1185 sfree(mtop->molblock);