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 "gromacs/topology/mtop_util.h"
44 #include "gromacs/legacyheaders/types/ifunc.h"
45 #include "gromacs/legacyheaders/types/inputrec.h"
47 #include "gromacs/topology/block.h"
48 #include "gromacs/topology/symtab.h"
49 #include "gromacs/topology/topology.h"
50 #include "gromacs/topology/topsort.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/smalloc.h"
54 static int gmx_mtop_maxresnr(const gmx_mtop_t *mtop, int maxres_renum)
61 for (mt = 0; mt < mtop->nmoltype; mt++)
63 atoms = &mtop->moltype[mt].atoms;
64 if (atoms->nres > maxres_renum)
66 for (r = 0; r < atoms->nres; r++)
68 if (atoms->resinfo[r].nr > maxresnr)
70 maxresnr = atoms->resinfo[r].nr;
79 void gmx_mtop_finalize(gmx_mtop_t *mtop)
83 mtop->maxres_renum = 1;
85 env = getenv("GMX_MAXRESRENUM");
88 sscanf(env, "%d", &mtop->maxres_renum);
90 if (mtop->maxres_renum == -1)
92 /* -1 signals renumber residues in all molecules */
93 mtop->maxres_renum = INT_MAX;
96 mtop->maxresnr = gmx_mtop_maxresnr(mtop, mtop->maxres_renum);
99 void gmx_mtop_count_atomtypes(const gmx_mtop_t *mtop, int state, int typecount[])
101 int i, mb, nmol, tpi;
104 for (i = 0; i < mtop->ffparams.atnr; ++i)
108 for (mb = 0; mb < mtop->nmolblock; ++mb)
110 nmol = mtop->molblock[mb].nmol;
111 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
112 for (i = 0; i < atoms->nr; ++i)
116 tpi = atoms->atom[i].type;
120 tpi = atoms->atom[i].typeB;
122 typecount[tpi] += nmol;
127 int ncg_mtop(const gmx_mtop_t *mtop)
133 for (mb = 0; mb < mtop->nmolblock; mb++)
136 mtop->molblock[mb].nmol*
137 mtop->moltype[mtop->molblock[mb].type].cgs.nr;
143 void gmx_mtop_remove_chargegroups(gmx_mtop_t *mtop)
149 for (mt = 0; mt < mtop->nmoltype; mt++)
151 cgs = &mtop->moltype[mt].cgs;
152 if (cgs->nr < mtop->moltype[mt].atoms.nr)
154 cgs->nr = mtop->moltype[mt].atoms.nr;
155 srenew(cgs->index, cgs->nr+1);
156 for (i = 0; i < cgs->nr+1; i++)
172 typedef struct gmx_mtop_atomlookup
174 const gmx_mtop_t *mtop;
178 } t_gmx_mtop_atomlookup;
181 gmx_mtop_atomlookup_t
182 gmx_mtop_atomlookup_init(const gmx_mtop_t *mtop)
184 t_gmx_mtop_atomlookup *alook;
186 int a_start, a_end, na, na_start = -1;
191 alook->nmb = mtop->nmolblock;
193 snew(alook->mba, alook->nmb);
196 for (mb = 0; mb < mtop->nmolblock; mb++)
198 na = mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
199 a_end = a_start + na;
201 alook->mba[mb].a_start = a_start;
202 alook->mba[mb].a_end = a_end;
203 alook->mba[mb].na_mol = mtop->molblock[mb].natoms_mol;
205 /* We start the binary search with the largest block */
206 if (mb == 0 || na > na_start)
208 alook->mb_start = mb;
218 gmx_mtop_atomlookup_t
219 gmx_mtop_atomlookup_settle_init(const gmx_mtop_t *mtop)
221 t_gmx_mtop_atomlookup *alook;
223 int na, na_start = -1;
225 alook = gmx_mtop_atomlookup_init(mtop);
227 /* Check if the starting molblock has settle */
228 if (mtop->moltype[mtop->molblock[alook->mb_start].type].ilist[F_SETTLE].nr == 0)
230 /* Search the largest molblock with settle */
231 alook->mb_start = -1;
232 for (mb = 0; mb < mtop->nmolblock; mb++)
234 if (mtop->moltype[mtop->molblock[mb].type].ilist[F_SETTLE].nr > 0)
236 na = alook->mba[mb].a_end - alook->mba[mb].a_start;
237 if (alook->mb_start == -1 || na > na_start)
239 alook->mb_start = mb;
245 if (alook->mb_start == -1)
247 gmx_incons("gmx_mtop_atomlookup_settle_init called without settles");
255 gmx_mtop_atomlookup_destroy(gmx_mtop_atomlookup_t alook)
261 void gmx_mtop_atomnr_to_atom(const gmx_mtop_atomlookup_t alook,
266 int a_start, atnr_mol;
269 if (atnr_global < 0 || atnr_global >= mtop->natoms)
271 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)",
272 atnr_global, 0, mtop->natoms-1);
278 mb = alook->mb_start;
282 a_start = alook->mba[mb].a_start;
283 if (atnr_global < a_start)
287 else if (atnr_global >= alook->mba[mb].a_end)
295 mb = ((mb0 + mb1 + 1)>>1);
298 atnr_mol = (atnr_global - a_start) % alook->mba[mb].na_mol;
300 *atom = &alook->mtop->moltype[alook->mtop->molblock[mb].type].atoms.atom[atnr_mol];
303 void gmx_mtop_atomnr_to_ilist(const gmx_mtop_atomlookup_t alook,
305 t_ilist **ilist_mol, int *atnr_offset)
308 int a_start, atnr_local;
311 if (atnr_global < 0 || atnr_global >= mtop->natoms)
313 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)",
314 atnr_global, 0, mtop->natoms-1);
320 mb = alook->mb_start;
324 a_start = alook->mba[mb].a_start;
325 if (atnr_global < a_start)
329 else if (atnr_global >= alook->mba[mb].a_end)
337 mb = ((mb0 + mb1 + 1)>>1);
340 *ilist_mol = alook->mtop->moltype[alook->mtop->molblock[mb].type].ilist;
342 atnr_local = (atnr_global - a_start) % alook->mba[mb].na_mol;
344 *atnr_offset = atnr_global - atnr_local;
347 void gmx_mtop_atomnr_to_molblock_ind(const gmx_mtop_atomlookup_t alook,
349 int *molb, int *molnr, int *atnr_mol)
355 if (atnr_global < 0 || atnr_global >= mtop->natoms)
357 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)",
358 atnr_global, 0, mtop->natoms-1);
364 mb = alook->mb_start;
368 a_start = alook->mba[mb].a_start;
369 if (atnr_global < a_start)
373 else if (atnr_global >= alook->mba[mb].a_end)
381 mb = ((mb0 + mb1 + 1)>>1);
385 *molnr = (atnr_global - a_start) / alook->mba[mb].na_mol;
386 *atnr_mol = atnr_global - a_start - (*molnr)*alook->mba[mb].na_mol;
389 void gmx_mtop_atominfo_global(const gmx_mtop_t *mtop, int atnr_global,
390 char **atomname, int *resnr, char **resname)
392 int mb, a_start, a_end, maxresnr, at_loc;
393 gmx_molblock_t *molb;
394 t_atoms *atoms = NULL;
396 if (atnr_global < 0 || atnr_global >= mtop->natoms)
398 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)",
399 atnr_global, 0, mtop->natoms-1);
404 maxresnr = mtop->maxresnr;
409 if (atoms->nres <= mtop->maxres_renum)
411 /* Single residue molecule, keep counting */
412 maxresnr += mtop->molblock[mb].nmol*atoms->nres;
416 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
418 a_end = a_start + mtop->molblock[mb].nmol*atoms->nr;
420 while (atnr_global >= a_end);
422 at_loc = (atnr_global - a_start) % atoms->nr;
423 *atomname = *(atoms->atomname[at_loc]);
424 if (atoms->nres > mtop->maxres_renum)
426 *resnr = atoms->resinfo[atoms->atom[at_loc].resind].nr;
430 /* Single residue molecule, keep counting */
431 *resnr = maxresnr + 1 + (atnr_global - a_start)/atoms->nr*atoms->nres + atoms->atom[at_loc].resind;
433 *resname = *(atoms->resinfo[atoms->atom[at_loc].resind].name);
436 typedef struct gmx_mtop_atomloop_all
438 const gmx_mtop_t *mtop;
445 } t_gmx_mtop_atomloop_all;
447 gmx_mtop_atomloop_all_t
448 gmx_mtop_atomloop_all_init(const gmx_mtop_t *mtop)
450 struct gmx_mtop_atomloop_all *aloop;
457 &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
459 aloop->maxresnr = mtop->maxresnr;
460 aloop->at_local = -1;
461 aloop->at_global = -1;
466 static void gmx_mtop_atomloop_all_destroy(gmx_mtop_atomloop_all_t aloop)
471 gmx_bool gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop,
472 int *at_global, t_atom **atom)
476 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
482 if (aloop->at_local >= aloop->atoms->nr)
484 if (aloop->atoms->nres <= aloop->mtop->maxres_renum)
486 /* Single residue molecule, increase the count with one */
487 aloop->maxresnr += aloop->atoms->nres;
491 if (aloop->mol >= aloop->mtop->molblock[aloop->mblock].nmol)
494 if (aloop->mblock >= aloop->mtop->nmolblock)
496 gmx_mtop_atomloop_all_destroy(aloop);
499 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
504 *at_global = aloop->at_global;
505 *atom = &aloop->atoms->atom[aloop->at_local];
510 void gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop,
511 char **atomname, int *resnr, char **resname)
515 *atomname = *(aloop->atoms->atomname[aloop->at_local]);
516 resind_mol = aloop->atoms->atom[aloop->at_local].resind;
517 *resnr = aloop->atoms->resinfo[resind_mol].nr;
518 if (aloop->atoms->nres <= aloop->mtop->maxres_renum)
520 *resnr = aloop->maxresnr + 1 + resind_mol;
522 *resname = *(aloop->atoms->resinfo[resind_mol].name);
525 void gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop,
526 gmx_moltype_t **moltype, int *at_mol)
528 *moltype = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type];
529 *at_mol = aloop->at_local;
532 typedef struct gmx_mtop_atomloop_block
534 const gmx_mtop_t *mtop;
538 } t_gmx_mtop_atomloop_block;
540 gmx_mtop_atomloop_block_t
541 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop)
543 struct gmx_mtop_atomloop_block *aloop;
549 aloop->atoms = &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
550 aloop->at_local = -1;
555 static void gmx_mtop_atomloop_block_destroy(gmx_mtop_atomloop_block_t aloop)
560 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
561 t_atom **atom, int *nmol)
565 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
570 if (aloop->at_local >= aloop->atoms->nr)
573 if (aloop->mblock >= aloop->mtop->nmolblock)
575 gmx_mtop_atomloop_block_destroy(aloop);
578 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
582 *atom = &aloop->atoms->atom[aloop->at_local];
583 *nmol = aloop->mtop->molblock[aloop->mblock].nmol;
588 typedef struct gmx_mtop_ilistloop
590 const gmx_mtop_t *mtop;
595 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop)
597 struct gmx_mtop_ilistloop *iloop;
607 static void gmx_mtop_ilistloop_destroy(gmx_mtop_ilistloop_t iloop)
612 gmx_bool gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
613 t_ilist **ilist_mol, int *nmol)
617 gmx_incons("gmx_mtop_ilistloop_next called without calling gmx_mtop_ilistloop_init");
621 if (iloop->mblock == iloop->mtop->nmolblock)
623 gmx_mtop_ilistloop_destroy(iloop);
628 iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
630 *nmol = iloop->mtop->molblock[iloop->mblock].nmol;
634 typedef struct gmx_mtop_ilistloop_all
636 const gmx_mtop_t *mtop;
640 } t_gmx_mtop_ilist_all;
642 gmx_mtop_ilistloop_all_t
643 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop)
645 struct gmx_mtop_ilistloop_all *iloop;
657 static void gmx_mtop_ilistloop_all_destroy(gmx_mtop_ilistloop_all_t iloop)
662 gmx_bool gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
663 t_ilist **ilist_mol, int *atnr_offset)
665 gmx_molblock_t *molb;
669 gmx_incons("gmx_mtop_ilistloop_all_next called without calling gmx_mtop_ilistloop_all_init");
674 iloop->a_offset += iloop->mtop->molblock[iloop->mblock].natoms_mol;
679 if (iloop->mol >= iloop->mtop->molblock[iloop->mblock].nmol)
683 if (iloop->mblock == iloop->mtop->nmolblock)
685 gmx_mtop_ilistloop_all_destroy(iloop);
691 iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
693 *atnr_offset = iloop->a_offset;
698 int gmx_mtop_ftype_count(const gmx_mtop_t *mtop, int ftype)
700 gmx_mtop_ilistloop_t iloop;
706 iloop = gmx_mtop_ilistloop_init(mtop);
707 while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
709 n += nmol*il[ftype].nr/(1+NRAL(ftype));
715 t_block gmx_mtop_global_cgs(const gmx_mtop_t *mtop)
717 t_block cgs_gl, *cgs_mol;
719 gmx_molblock_t *molb;
722 /* In most cases this is too much, but we realloc at the end */
723 snew(cgs_gl.index, mtop->natoms+1);
727 for (mb = 0; mb < mtop->nmolblock; mb++)
729 molb = &mtop->molblock[mb];
730 cgs_mol = &mtop->moltype[molb->type].cgs;
731 for (mol = 0; mol < molb->nmol; mol++)
733 for (cg = 0; cg < cgs_mol->nr; cg++)
735 cgs_gl.index[cgs_gl.nr+1] =
736 cgs_gl.index[cgs_gl.nr] +
737 cgs_mol->index[cg+1] - cgs_mol->index[cg];
742 cgs_gl.nalloc_index = cgs_gl.nr + 1;
743 srenew(cgs_gl.index, cgs_gl.nalloc_index);
748 static void atomcat(t_atoms *dest, t_atoms *src, int copies,
749 int maxres_renum, int *maxresnr)
753 int destnr = dest->nr;
757 size = destnr+copies*srcnr;
758 srenew(dest->atom, size);
759 srenew(dest->atomname, size);
760 srenew(dest->atomtype, size);
761 srenew(dest->atomtypeB, size);
765 size = dest->nres+copies*src->nres;
766 srenew(dest->resinfo, size);
769 /* residue information */
770 for (l = dest->nres, j = 0; (j < copies); j++, l += src->nres)
772 memcpy((char *) &(dest->resinfo[l]), (char *) &(src->resinfo[0]),
773 (size_t)(src->nres*sizeof(src->resinfo[0])));
776 for (l = destnr, j = 0; (j < copies); j++, l += srcnr)
778 memcpy((char *) &(dest->atomname[l]), (char *) &(src->atomname[0]),
779 (size_t)(srcnr*sizeof(src->atomname[0])));
780 memcpy((char *) &(dest->atomtype[l]), (char *) &(src->atomtype[0]),
781 (size_t)(srcnr*sizeof(src->atomtype[0])));
782 memcpy((char *) &(dest->atomtypeB[l]), (char *) &(src->atomtypeB[0]),
783 (size_t)(srcnr*sizeof(src->atomtypeB[0])));
784 memcpy((char *) &(dest->atom[l]), (char *) &(src->atom[0]),
785 (size_t)(srcnr*sizeof(src->atom[0])));
788 /* Increment residue indices */
789 for (l = destnr, j = 0; (j < copies); j++)
791 for (i = 0; (i < srcnr); i++, l++)
793 dest->atom[l].resind = dest->nres+j*src->nres+src->atom[i].resind;
797 if (src->nres <= maxres_renum)
799 /* Single residue molecule, continue counting residues */
800 for (j = 0; (j < copies); j++)
802 for (l = 0; l < src->nres; l++)
805 dest->resinfo[dest->nres+j*src->nres+l].nr = *maxresnr;
810 dest->nres += copies*src->nres;
811 dest->nr += copies*src->nr;
814 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t *mtop)
818 gmx_molblock_t *molb;
820 init_t_atoms(&atoms, 0, FALSE);
822 maxresnr = mtop->maxresnr;
823 for (mb = 0; mb < mtop->nmolblock; mb++)
825 molb = &mtop->molblock[mb];
826 atomcat(&atoms, &mtop->moltype[molb->type].atoms, molb->nmol,
827 mtop->maxres_renum, &maxresnr);
833 void gmx_mtop_make_atomic_charge_groups(gmx_mtop_t *mtop,
834 gmx_bool bKeepSingleMolCG)
839 for (mb = 0; mb < mtop->nmolblock; mb++)
841 cgs_mol = &mtop->moltype[mtop->molblock[mb].type].cgs;
842 if (!(bKeepSingleMolCG && cgs_mol->nr == 1))
844 cgs_mol->nr = mtop->molblock[mb].natoms_mol;
845 cgs_mol->nalloc_index = cgs_mol->nr + 1;
846 srenew(cgs_mol->index, cgs_mol->nalloc_index);
847 for (cg = 0; cg < cgs_mol->nr+1; cg++)
849 cgs_mol->index[cg] = cg;
856 * The cat routines below are old code from src/kernel/topcat.c
859 static void blockcat(t_block *dest, t_block *src, int copies)
861 int i, j, l, nra, size;
865 size = (dest->nr+copies*src->nr+1);
866 srenew(dest->index, size);
869 nra = dest->index[dest->nr];
870 for (l = dest->nr, j = 0; (j < copies); j++)
872 for (i = 0; (i < src->nr); i++)
874 dest->index[l++] = nra + src->index[i];
876 nra += src->index[src->nr];
878 dest->nr += copies*src->nr;
879 dest->index[dest->nr] = nra;
882 static void blockacat(t_blocka *dest, t_blocka *src, int copies,
886 int destnr = dest->nr;
887 int destnra = dest->nra;
891 size = (dest->nr+copies*src->nr+1);
892 srenew(dest->index, size);
896 size = (dest->nra+copies*src->nra);
897 srenew(dest->a, size);
900 for (l = destnr, j = 0; (j < copies); j++)
902 for (i = 0; (i < src->nr); i++)
904 dest->index[l++] = dest->nra+src->index[i];
906 dest->nra += src->nra;
908 for (l = destnra, j = 0; (j < copies); j++)
910 for (i = 0; (i < src->nra); i++)
912 dest->a[l++] = dnum+src->a[i];
917 dest->index[dest->nr] = dest->nra;
920 static void ilistcat(int ftype, t_ilist *dest, t_ilist *src, int copies,
927 dest->nalloc = dest->nr + copies*src->nr;
928 srenew(dest->iatoms, dest->nalloc);
930 for (c = 0; c < copies; c++)
932 for (i = 0; i < src->nr; )
934 dest->iatoms[dest->nr++] = src->iatoms[i++];
935 for (a = 0; a < nral; a++)
937 dest->iatoms[dest->nr++] = dnum + src->iatoms[i++];
944 static void set_posres_params(t_idef *idef, gmx_molblock_t *molb,
945 int i0, int a_offset)
951 il = &idef->il[F_POSRES];
953 idef->iparams_posres_nalloc = i1;
954 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
955 for (i = i0; i < i1; i++)
957 ip = &idef->iparams_posres[i];
958 /* Copy the force constants */
959 *ip = idef->iparams[il->iatoms[i*2]];
960 a_molb = il->iatoms[i*2+1] - a_offset;
961 if (molb->nposres_xA == 0)
963 gmx_incons("Position restraint coordinates are missing");
965 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
966 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
967 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
968 if (molb->nposres_xB > 0)
970 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
971 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
972 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
976 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
977 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
978 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
980 /* Set the parameter index for idef->iparams_posre */
985 static void set_fbposres_params(t_idef *idef, gmx_molblock_t *molb,
986 int i0, int a_offset)
992 il = &idef->il[F_FBPOSRES];
994 idef->iparams_fbposres_nalloc = i1;
995 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
996 for (i = i0; i < i1; i++)
998 ip = &idef->iparams_fbposres[i];
999 /* Copy the force constants */
1000 *ip = idef->iparams[il->iatoms[i*2]];
1001 a_molb = il->iatoms[i*2+1] - a_offset;
1002 if (molb->nposres_xA == 0)
1004 gmx_incons("Position restraint coordinates are missing");
1006 /* Take flat-bottom posres reference from normal position restraints */
1007 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
1008 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
1009 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
1010 /* Note: no B-type for flat-bottom posres */
1012 /* Set the parameter index for idef->iparams_posre */
1013 il->iatoms[i*2] = i;
1017 static void gen_local_top(const gmx_mtop_t *mtop, const t_inputrec *ir,
1018 gmx_bool bMergeConstr,
1019 gmx_localtop_t *top)
1021 int mb, srcnr, destnr, ftype, ftype_dest, mt, natoms, mol, nposre_old, nfbposre_old;
1022 gmx_molblock_t *molb;
1023 gmx_moltype_t *molt;
1024 const gmx_ffparams_t *ffp;
1027 gmx_mtop_atomloop_all_t aloop;
1031 top->atomtypes = mtop->atomtypes;
1033 ffp = &mtop->ffparams;
1036 idef->ntypes = ffp->ntypes;
1037 idef->atnr = ffp->atnr;
1038 idef->functype = ffp->functype;
1039 idef->iparams = ffp->iparams;
1040 idef->iparams_posres = NULL;
1041 idef->iparams_posres_nalloc = 0;
1042 idef->iparams_fbposres = NULL;
1043 idef->iparams_fbposres_nalloc = 0;
1044 idef->fudgeQQ = ffp->fudgeQQ;
1045 idef->cmap_grid = ffp->cmap_grid;
1046 idef->ilsort = ilsortUNKNOWN;
1048 init_block(&top->cgs);
1049 init_blocka(&top->excls);
1050 for (ftype = 0; ftype < F_NRE; ftype++)
1052 idef->il[ftype].nr = 0;
1053 idef->il[ftype].nalloc = 0;
1054 idef->il[ftype].iatoms = NULL;
1058 for (mb = 0; mb < mtop->nmolblock; mb++)
1060 molb = &mtop->molblock[mb];
1061 molt = &mtop->moltype[molb->type];
1063 srcnr = molt->atoms.nr;
1066 blockcat(&top->cgs, &molt->cgs, molb->nmol);
1068 blockacat(&top->excls, &molt->excls, molb->nmol, destnr, srcnr);
1070 nposre_old = idef->il[F_POSRES].nr;
1071 nfbposre_old = idef->il[F_FBPOSRES].nr;
1072 for (ftype = 0; ftype < F_NRE; ftype++)
1075 ftype == F_CONSTR && molt->ilist[F_CONSTRNC].nr > 0)
1077 /* Merge all constrains into one ilist.
1078 * This simplifies the constraint code.
1080 for (mol = 0; mol < molb->nmol; mol++)
1082 ilistcat(ftype, &idef->il[F_CONSTR], &molt->ilist[F_CONSTR],
1083 1, destnr+mol*srcnr, srcnr);
1084 ilistcat(ftype, &idef->il[F_CONSTR], &molt->ilist[F_CONSTRNC],
1085 1, destnr+mol*srcnr, srcnr);
1088 else if (!(bMergeConstr && ftype == F_CONSTRNC))
1090 ilistcat(ftype, &idef->il[ftype], &molt->ilist[ftype],
1091 molb->nmol, destnr, srcnr);
1094 if (idef->il[F_POSRES].nr > nposre_old)
1096 /* Executing this line line stops gmxdump -sys working
1097 * correctly. I'm not aware there's an elegant fix. */
1098 set_posres_params(idef, molb, nposre_old/2, natoms);
1100 if (idef->il[F_FBPOSRES].nr > nfbposre_old)
1102 set_fbposres_params(idef, molb, nfbposre_old/2, natoms);
1105 natoms += molb->nmol*srcnr;
1110 top->idef.ilsort = ilsortUNKNOWN;
1114 if (ir->efep != efepNO && gmx_mtop_bondeds_free_energy(mtop))
1116 snew(qA, mtop->natoms);
1117 snew(qB, mtop->natoms);
1118 aloop = gmx_mtop_atomloop_all_init(mtop);
1119 while (gmx_mtop_atomloop_all_next(aloop, &ag, &atom))
1124 gmx_sort_ilist_fe(&top->idef, qA, qB);
1130 top->idef.ilsort = ilsortNO_FE;
1135 gmx_localtop_t *gmx_mtop_generate_local_top(const gmx_mtop_t *mtop,
1136 const t_inputrec *ir)
1138 gmx_localtop_t *top;
1142 gen_local_top(mtop, ir, TRUE, top);
1147 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop)
1150 gmx_localtop_t ltop;
1153 gen_local_top(mtop, NULL, FALSE, <op);
1155 top.name = mtop->name;
1156 top.idef = ltop.idef;
1157 top.atomtypes = ltop.atomtypes;
1159 top.excls = ltop.excls;
1160 top.atoms = gmx_mtop_global_atoms(mtop);
1161 top.mols = mtop->mols;
1162 top.symtab = mtop->symtab;
1164 /* We only need to free the moltype and molblock data,
1165 * all other pointers have been copied to top.
1167 * Well, except for the group data, but we can't free those, because they
1168 * are used somewhere even after a call to this function.
1170 for (mt = 0; mt < mtop->nmoltype; mt++)
1172 done_moltype(&mtop->moltype[mt]);
1174 sfree(mtop->moltype);
1176 for (mb = 0; mb < mtop->nmolblock; mb++)
1178 done_molblock(&mtop->molblock[mb]);
1180 sfree(mtop->molblock);