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"
42 #include "gromacs/legacyheaders/types/ifunc.h"
43 #include "gromacs/legacyheaders/types/inputrec.h"
45 #include "gromacs/topology/block.h"
46 #include "gromacs/topology/symtab.h"
47 #include "gromacs/topology/topology.h"
48 #include "gromacs/topology/topsort.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/smalloc.h"
52 static int gmx_mtop_maxresnr(const gmx_mtop_t *mtop, int maxres_renum)
59 for (mt = 0; mt < mtop->nmoltype; mt++)
61 atoms = &mtop->moltype[mt].atoms;
62 if (atoms->nres > maxres_renum)
64 for (r = 0; r < atoms->nres; r++)
66 if (atoms->resinfo[r].nr > maxresnr)
68 maxresnr = atoms->resinfo[r].nr;
77 void gmx_mtop_finalize(gmx_mtop_t *mtop)
81 mtop->maxres_renum = 1;
83 env = getenv("GMX_MAXRESRENUM");
86 sscanf(env, "%d", &mtop->maxres_renum);
88 if (mtop->maxres_renum == -1)
90 /* -1 signals renumber residues in all molecules */
91 mtop->maxres_renum = INT_MAX;
94 mtop->maxresnr = gmx_mtop_maxresnr(mtop, mtop->maxres_renum);
97 void gmx_mtop_count_atomtypes(const gmx_mtop_t *mtop, int state, int typecount[])
102 for (i = 0; i < mtop->ffparams.atnr; ++i)
106 for (mb = 0; mb < mtop->nmolblock; ++mb)
108 nmol = mtop->molblock[mb].nmol;
109 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
110 for (i = 0; i < atoms->nr; ++i)
114 tpi = atoms->atom[i].type;
118 tpi = atoms->atom[i].typeB;
120 typecount[tpi] += nmol;
125 int ncg_mtop(const gmx_mtop_t *mtop)
131 for (mb = 0; mb < mtop->nmolblock; mb++)
134 mtop->molblock[mb].nmol*
135 mtop->moltype[mtop->molblock[mb].type].cgs.nr;
141 void gmx_mtop_remove_chargegroups(gmx_mtop_t *mtop)
147 for (mt = 0; mt < mtop->nmoltype; mt++)
149 cgs = &mtop->moltype[mt].cgs;
150 if (cgs->nr < mtop->moltype[mt].atoms.nr)
152 cgs->nr = mtop->moltype[mt].atoms.nr;
153 srenew(cgs->index, cgs->nr+1);
154 for (i = 0; i < cgs->nr+1; i++)
170 typedef struct gmx_mtop_atomlookup
172 const gmx_mtop_t *mtop;
176 } t_gmx_mtop_atomlookup;
179 gmx_mtop_atomlookup_t
180 gmx_mtop_atomlookup_init(const gmx_mtop_t *mtop)
182 t_gmx_mtop_atomlookup *alook;
184 int a_start, a_end, na, na_start = -1;
189 alook->nmb = mtop->nmolblock;
191 snew(alook->mba, alook->nmb);
194 for (mb = 0; mb < mtop->nmolblock; mb++)
196 na = mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
197 a_end = a_start + na;
199 alook->mba[mb].a_start = a_start;
200 alook->mba[mb].a_end = a_end;
201 alook->mba[mb].na_mol = mtop->molblock[mb].natoms_mol;
203 /* We start the binary search with the largest block */
204 if (mb == 0 || na > na_start)
206 alook->mb_start = mb;
216 gmx_mtop_atomlookup_t
217 gmx_mtop_atomlookup_settle_init(const gmx_mtop_t *mtop)
219 t_gmx_mtop_atomlookup *alook;
221 int na, na_start = -1;
223 alook = gmx_mtop_atomlookup_init(mtop);
225 /* Check if the starting molblock has settle */
226 if (mtop->moltype[mtop->molblock[alook->mb_start].type].ilist[F_SETTLE].nr == 0)
228 /* Search the largest molblock with settle */
229 alook->mb_start = -1;
230 for (mb = 0; mb < mtop->nmolblock; mb++)
232 if (mtop->moltype[mtop->molblock[mb].type].ilist[F_SETTLE].nr > 0)
234 na = alook->mba[mb].a_end - alook->mba[mb].a_start;
235 if (alook->mb_start == -1 || na > na_start)
237 alook->mb_start = mb;
243 if (alook->mb_start == -1)
245 gmx_incons("gmx_mtop_atomlookup_settle_init called without settles");
253 gmx_mtop_atomlookup_destroy(gmx_mtop_atomlookup_t alook)
259 void gmx_mtop_atomnr_to_atom(const gmx_mtop_atomlookup_t alook,
264 int a_start, atnr_mol;
267 if (atnr_global < 0 || atnr_global >= mtop->natoms)
269 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)",
270 atnr_global, 0, mtop->natoms-1);
276 mb = alook->mb_start;
280 a_start = alook->mba[mb].a_start;
281 if (atnr_global < a_start)
285 else if (atnr_global >= alook->mba[mb].a_end)
293 mb = ((mb0 + mb1 + 1)>>1);
296 atnr_mol = (atnr_global - a_start) % alook->mba[mb].na_mol;
298 *atom = &alook->mtop->moltype[alook->mtop->molblock[mb].type].atoms.atom[atnr_mol];
301 void gmx_mtop_atomnr_to_ilist(const gmx_mtop_atomlookup_t alook,
303 t_ilist **ilist_mol, int *atnr_offset)
306 int a_start, atnr_local;
309 if (atnr_global < 0 || atnr_global >= mtop->natoms)
311 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)",
312 atnr_global, 0, mtop->natoms-1);
318 mb = alook->mb_start;
322 a_start = alook->mba[mb].a_start;
323 if (atnr_global < a_start)
327 else if (atnr_global >= alook->mba[mb].a_end)
335 mb = ((mb0 + mb1 + 1)>>1);
338 *ilist_mol = alook->mtop->moltype[alook->mtop->molblock[mb].type].ilist;
340 atnr_local = (atnr_global - a_start) % alook->mba[mb].na_mol;
342 *atnr_offset = atnr_global - atnr_local;
345 void gmx_mtop_atomnr_to_molblock_ind(const gmx_mtop_atomlookup_t alook,
347 int *molb, int *molnr, int *atnr_mol)
353 if (atnr_global < 0 || atnr_global >= mtop->natoms)
355 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)",
356 atnr_global, 0, mtop->natoms-1);
362 mb = alook->mb_start;
366 a_start = alook->mba[mb].a_start;
367 if (atnr_global < a_start)
371 else if (atnr_global >= alook->mba[mb].a_end)
379 mb = ((mb0 + mb1 + 1)>>1);
383 *molnr = (atnr_global - a_start) / alook->mba[mb].na_mol;
384 *atnr_mol = atnr_global - a_start - (*molnr)*alook->mba[mb].na_mol;
387 void gmx_mtop_atominfo_global(const gmx_mtop_t *mtop, int atnr_global,
388 char **atomname, int *resnr, char **resname)
390 int mb, a_start, a_end, maxresnr, at_loc;
391 gmx_molblock_t *molb;
392 t_atoms *atoms = NULL;
394 if (atnr_global < 0 || atnr_global >= mtop->natoms)
396 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)",
397 atnr_global, 0, mtop->natoms-1);
402 maxresnr = mtop->maxresnr;
407 if (atoms->nres <= mtop->maxres_renum)
409 /* Single residue molecule, keep counting */
410 maxresnr += mtop->molblock[mb].nmol*atoms->nres;
414 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
416 a_end = a_start + mtop->molblock[mb].nmol*atoms->nr;
418 while (atnr_global >= a_end);
420 at_loc = (atnr_global - a_start) % atoms->nr;
421 *atomname = *(atoms->atomname[at_loc]);
422 if (atoms->nres > mtop->maxres_renum)
424 *resnr = atoms->resinfo[atoms->atom[at_loc].resind].nr;
428 /* Single residue molecule, keep counting */
429 *resnr = maxresnr + 1 + (atnr_global - a_start)/atoms->nr*atoms->nres + atoms->atom[at_loc].resind;
431 *resname = *(atoms->resinfo[atoms->atom[at_loc].resind].name);
434 typedef struct gmx_mtop_atomloop_all
436 const gmx_mtop_t *mtop;
443 } t_gmx_mtop_atomloop_all;
445 gmx_mtop_atomloop_all_t
446 gmx_mtop_atomloop_all_init(const gmx_mtop_t *mtop)
448 struct gmx_mtop_atomloop_all *aloop;
455 &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
457 aloop->maxresnr = mtop->maxresnr;
458 aloop->at_local = -1;
459 aloop->at_global = -1;
464 static void gmx_mtop_atomloop_all_destroy(gmx_mtop_atomloop_all_t aloop)
469 gmx_bool gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop,
470 int *at_global, t_atom **atom)
474 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
480 if (aloop->at_local >= aloop->atoms->nr)
482 if (aloop->atoms->nres <= aloop->mtop->maxres_renum)
484 /* Single residue molecule, increase the count with one */
485 aloop->maxresnr += aloop->atoms->nres;
489 if (aloop->mol >= aloop->mtop->molblock[aloop->mblock].nmol)
492 if (aloop->mblock >= aloop->mtop->nmolblock)
494 gmx_mtop_atomloop_all_destroy(aloop);
497 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
502 *at_global = aloop->at_global;
503 *atom = &aloop->atoms->atom[aloop->at_local];
508 void gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop,
509 char **atomname, int *resnr, char **resname)
513 *atomname = *(aloop->atoms->atomname[aloop->at_local]);
514 resind_mol = aloop->atoms->atom[aloop->at_local].resind;
515 *resnr = aloop->atoms->resinfo[resind_mol].nr;
516 if (aloop->atoms->nres <= aloop->mtop->maxres_renum)
518 *resnr = aloop->maxresnr + 1 + resind_mol;
520 *resname = *(aloop->atoms->resinfo[resind_mol].name);
523 void gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop,
524 gmx_moltype_t **moltype, int *at_mol)
526 *moltype = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type];
527 *at_mol = aloop->at_local;
530 typedef struct gmx_mtop_atomloop_block
532 const gmx_mtop_t *mtop;
536 } t_gmx_mtop_atomloop_block;
538 gmx_mtop_atomloop_block_t
539 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop)
541 struct gmx_mtop_atomloop_block *aloop;
547 aloop->atoms = &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
548 aloop->at_local = -1;
553 static void gmx_mtop_atomloop_block_destroy(gmx_mtop_atomloop_block_t aloop)
558 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
559 t_atom **atom, int *nmol)
563 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
568 if (aloop->at_local >= aloop->atoms->nr)
571 if (aloop->mblock >= aloop->mtop->nmolblock)
573 gmx_mtop_atomloop_block_destroy(aloop);
576 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
580 *atom = &aloop->atoms->atom[aloop->at_local];
581 *nmol = aloop->mtop->molblock[aloop->mblock].nmol;
586 typedef struct gmx_mtop_ilistloop
588 const gmx_mtop_t *mtop;
593 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop)
595 struct gmx_mtop_ilistloop *iloop;
605 static void gmx_mtop_ilistloop_destroy(gmx_mtop_ilistloop_t iloop)
610 gmx_bool gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
611 t_ilist **ilist_mol, int *nmol)
615 gmx_incons("gmx_mtop_ilistloop_next called without calling gmx_mtop_ilistloop_init");
619 if (iloop->mblock == iloop->mtop->nmolblock)
621 gmx_mtop_ilistloop_destroy(iloop);
626 iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
628 *nmol = iloop->mtop->molblock[iloop->mblock].nmol;
632 typedef struct gmx_mtop_ilistloop_all
634 const gmx_mtop_t *mtop;
638 } t_gmx_mtop_ilist_all;
640 gmx_mtop_ilistloop_all_t
641 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop)
643 struct gmx_mtop_ilistloop_all *iloop;
655 static void gmx_mtop_ilistloop_all_destroy(gmx_mtop_ilistloop_all_t iloop)
660 gmx_bool gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
661 t_ilist **ilist_mol, int *atnr_offset)
663 gmx_molblock_t *molb;
667 gmx_incons("gmx_mtop_ilistloop_all_next called without calling gmx_mtop_ilistloop_all_init");
672 iloop->a_offset += iloop->mtop->molblock[iloop->mblock].natoms_mol;
677 if (iloop->mol >= iloop->mtop->molblock[iloop->mblock].nmol)
681 if (iloop->mblock == iloop->mtop->nmolblock)
683 gmx_mtop_ilistloop_all_destroy(iloop);
689 iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
691 *atnr_offset = iloop->a_offset;
696 int gmx_mtop_ftype_count(const gmx_mtop_t *mtop, int ftype)
698 gmx_mtop_ilistloop_t iloop;
704 iloop = gmx_mtop_ilistloop_init(mtop);
705 while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
707 n += nmol*il[ftype].nr/(1+NRAL(ftype));
713 t_block gmx_mtop_global_cgs(const gmx_mtop_t *mtop)
715 t_block cgs_gl, *cgs_mol;
717 gmx_molblock_t *molb;
720 /* In most cases this is too much, but we realloc at the end */
721 snew(cgs_gl.index, mtop->natoms+1);
725 for (mb = 0; mb < mtop->nmolblock; mb++)
727 molb = &mtop->molblock[mb];
728 cgs_mol = &mtop->moltype[molb->type].cgs;
729 for (mol = 0; mol < molb->nmol; mol++)
731 for (cg = 0; cg < cgs_mol->nr; cg++)
733 cgs_gl.index[cgs_gl.nr+1] =
734 cgs_gl.index[cgs_gl.nr] +
735 cgs_mol->index[cg+1] - cgs_mol->index[cg];
740 cgs_gl.nalloc_index = cgs_gl.nr + 1;
741 srenew(cgs_gl.index, cgs_gl.nalloc_index);
746 static void atomcat(t_atoms *dest, t_atoms *src, int copies,
747 int maxres_renum, int *maxresnr)
751 int destnr = dest->nr;
755 size = destnr+copies*srcnr;
756 srenew(dest->atom, size);
757 srenew(dest->atomname, size);
758 srenew(dest->atomtype, size);
759 srenew(dest->atomtypeB, size);
763 size = dest->nres+copies*src->nres;
764 srenew(dest->resinfo, size);
767 /* residue information */
768 for (l = dest->nres, j = 0; (j < copies); j++, l += src->nres)
770 memcpy((char *) &(dest->resinfo[l]), (char *) &(src->resinfo[0]),
771 (size_t)(src->nres*sizeof(src->resinfo[0])));
774 for (l = destnr, j = 0; (j < copies); j++, l += srcnr)
776 memcpy((char *) &(dest->atomname[l]), (char *) &(src->atomname[0]),
777 (size_t)(srcnr*sizeof(src->atomname[0])));
778 memcpy((char *) &(dest->atomtype[l]), (char *) &(src->atomtype[0]),
779 (size_t)(srcnr*sizeof(src->atomtype[0])));
780 memcpy((char *) &(dest->atomtypeB[l]), (char *) &(src->atomtypeB[0]),
781 (size_t)(srcnr*sizeof(src->atomtypeB[0])));
782 memcpy((char *) &(dest->atom[l]), (char *) &(src->atom[0]),
783 (size_t)(srcnr*sizeof(src->atom[0])));
786 /* Increment residue indices */
787 for (l = destnr, j = 0; (j < copies); j++)
789 for (i = 0; (i < srcnr); i++, l++)
791 dest->atom[l].resind = dest->nres+j*src->nres+src->atom[i].resind;
795 if (src->nres <= maxres_renum)
797 /* Single residue molecule, continue counting residues */
798 for (j = 0; (j < copies); j++)
800 for (l = 0; l < src->nres; l++)
803 dest->resinfo[dest->nres+j*src->nres+l].nr = *maxresnr;
808 dest->nres += copies*src->nres;
809 dest->nr += copies*src->nr;
812 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t *mtop)
816 gmx_molblock_t *molb;
818 init_t_atoms(&atoms, 0, FALSE);
820 maxresnr = mtop->maxresnr;
821 for (mb = 0; mb < mtop->nmolblock; mb++)
823 molb = &mtop->molblock[mb];
824 atomcat(&atoms, &mtop->moltype[molb->type].atoms, molb->nmol,
825 mtop->maxres_renum, &maxresnr);
831 void gmx_mtop_make_atomic_charge_groups(gmx_mtop_t *mtop,
832 gmx_bool bKeepSingleMolCG)
837 for (mb = 0; mb < mtop->nmolblock; mb++)
839 cgs_mol = &mtop->moltype[mtop->molblock[mb].type].cgs;
840 if (!(bKeepSingleMolCG && cgs_mol->nr == 1))
842 cgs_mol->nr = mtop->molblock[mb].natoms_mol;
843 cgs_mol->nalloc_index = cgs_mol->nr + 1;
844 srenew(cgs_mol->index, cgs_mol->nalloc_index);
845 for (cg = 0; cg < cgs_mol->nr+1; cg++)
847 cgs_mol->index[cg] = cg;
854 * The cat routines below are old code from src/kernel/topcat.c
857 static void blockcat(t_block *dest, t_block *src, int copies)
859 int i, j, l, nra, size;
863 size = (dest->nr+copies*src->nr+1);
864 srenew(dest->index, size);
867 nra = dest->index[dest->nr];
868 for (l = dest->nr, j = 0; (j < copies); j++)
870 for (i = 0; (i < src->nr); i++)
872 dest->index[l++] = nra + src->index[i];
874 nra += src->index[src->nr];
876 dest->nr += copies*src->nr;
877 dest->index[dest->nr] = nra;
880 static void blockacat(t_blocka *dest, t_blocka *src, int copies,
884 int destnr = dest->nr;
885 int destnra = dest->nra;
889 size = (dest->nr+copies*src->nr+1);
890 srenew(dest->index, size);
894 size = (dest->nra+copies*src->nra);
895 srenew(dest->a, size);
898 for (l = destnr, j = 0; (j < copies); j++)
900 for (i = 0; (i < src->nr); i++)
902 dest->index[l++] = dest->nra+src->index[i];
904 dest->nra += src->nra;
906 for (l = destnra, j = 0; (j < copies); j++)
908 for (i = 0; (i < src->nra); i++)
910 dest->a[l++] = dnum+src->a[i];
915 dest->index[dest->nr] = dest->nra;
918 static void ilistcat(int ftype, t_ilist *dest, t_ilist *src, int copies,
925 dest->nalloc = dest->nr + copies*src->nr;
926 srenew(dest->iatoms, dest->nalloc);
928 for (c = 0; c < copies; c++)
930 for (i = 0; i < src->nr; )
932 dest->iatoms[dest->nr++] = src->iatoms[i++];
933 for (a = 0; a < nral; a++)
935 dest->iatoms[dest->nr++] = dnum + src->iatoms[i++];
942 static void set_posres_params(t_idef *idef, gmx_molblock_t *molb,
943 int i0, int a_offset)
949 il = &idef->il[F_POSRES];
951 idef->iparams_posres_nalloc = i1;
952 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
953 for (i = i0; i < i1; i++)
955 ip = &idef->iparams_posres[i];
956 /* Copy the force constants */
957 *ip = idef->iparams[il->iatoms[i*2]];
958 a_molb = il->iatoms[i*2+1] - a_offset;
959 if (molb->nposres_xA == 0)
961 gmx_incons("Position restraint coordinates are missing");
963 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
964 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
965 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
966 if (molb->nposres_xB > 0)
968 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
969 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
970 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
974 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
975 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
976 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
978 /* Set the parameter index for idef->iparams_posre */
983 static void set_fbposres_params(t_idef *idef, gmx_molblock_t *molb,
984 int i0, int a_offset)
990 il = &idef->il[F_FBPOSRES];
992 idef->iparams_fbposres_nalloc = i1;
993 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
994 for (i = i0; i < i1; i++)
996 ip = &idef->iparams_fbposres[i];
997 /* Copy the force constants */
998 *ip = idef->iparams[il->iatoms[i*2]];
999 a_molb = il->iatoms[i*2+1] - a_offset;
1000 if (molb->nposres_xA == 0)
1002 gmx_incons("Position restraint coordinates are missing");
1004 /* Take flat-bottom posres reference from normal position restraints */
1005 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
1006 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
1007 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
1008 /* Note: no B-type for flat-bottom posres */
1010 /* Set the parameter index for idef->iparams_posre */
1011 il->iatoms[i*2] = i;
1015 static void gen_local_top(const gmx_mtop_t *mtop, const t_inputrec *ir,
1016 gmx_bool bMergeConstr,
1017 gmx_localtop_t *top)
1019 int mb, srcnr, destnr, ftype, ftype_dest, mt, natoms, mol, nposre_old, nfbposre_old;
1020 gmx_molblock_t *molb;
1021 gmx_moltype_t *molt;
1022 const gmx_ffparams_t *ffp;
1025 gmx_mtop_atomloop_all_t aloop;
1029 top->atomtypes = mtop->atomtypes;
1031 ffp = &mtop->ffparams;
1034 idef->ntypes = ffp->ntypes;
1035 idef->atnr = ffp->atnr;
1036 idef->functype = ffp->functype;
1037 idef->iparams = ffp->iparams;
1038 idef->iparams_posres = NULL;
1039 idef->iparams_posres_nalloc = 0;
1040 idef->iparams_fbposres = NULL;
1041 idef->iparams_fbposres_nalloc = 0;
1042 idef->fudgeQQ = ffp->fudgeQQ;
1043 idef->cmap_grid = ffp->cmap_grid;
1044 idef->ilsort = ilsortUNKNOWN;
1046 init_block(&top->cgs);
1047 init_blocka(&top->excls);
1048 for (ftype = 0; ftype < F_NRE; ftype++)
1050 idef->il[ftype].nr = 0;
1051 idef->il[ftype].nalloc = 0;
1052 idef->il[ftype].iatoms = NULL;
1056 for (mb = 0; mb < mtop->nmolblock; mb++)
1058 molb = &mtop->molblock[mb];
1059 molt = &mtop->moltype[molb->type];
1061 srcnr = molt->atoms.nr;
1064 blockcat(&top->cgs, &molt->cgs, molb->nmol);
1066 blockacat(&top->excls, &molt->excls, molb->nmol, destnr, srcnr);
1068 nposre_old = idef->il[F_POSRES].nr;
1069 nfbposre_old = idef->il[F_FBPOSRES].nr;
1070 for (ftype = 0; ftype < F_NRE; ftype++)
1073 ftype == F_CONSTR && molt->ilist[F_CONSTRNC].nr > 0)
1075 /* Merge all constrains into one ilist.
1076 * This simplifies the constraint code.
1078 for (mol = 0; mol < molb->nmol; mol++)
1080 ilistcat(ftype, &idef->il[F_CONSTR], &molt->ilist[F_CONSTR],
1081 1, destnr+mol*srcnr, srcnr);
1082 ilistcat(ftype, &idef->il[F_CONSTR], &molt->ilist[F_CONSTRNC],
1083 1, destnr+mol*srcnr, srcnr);
1086 else if (!(bMergeConstr && ftype == F_CONSTRNC))
1088 ilistcat(ftype, &idef->il[ftype], &molt->ilist[ftype],
1089 molb->nmol, destnr, srcnr);
1092 if (idef->il[F_POSRES].nr > nposre_old)
1094 /* Executing this line line stops gmxdump -sys working
1095 * correctly. I'm not aware there's an elegant fix. */
1096 set_posres_params(idef, molb, nposre_old/2, natoms);
1098 if (idef->il[F_FBPOSRES].nr > nfbposre_old)
1100 set_fbposres_params(idef, molb, nfbposre_old/2, natoms);
1103 natoms += molb->nmol*srcnr;
1108 top->idef.ilsort = ilsortUNKNOWN;
1112 if (ir->efep != efepNO && gmx_mtop_bondeds_free_energy(mtop))
1114 snew(qA, mtop->natoms);
1115 snew(qB, mtop->natoms);
1116 aloop = gmx_mtop_atomloop_all_init(mtop);
1117 while (gmx_mtop_atomloop_all_next(aloop, &ag, &atom))
1122 gmx_sort_ilist_fe(&top->idef, qA, qB);
1128 top->idef.ilsort = ilsortNO_FE;
1133 gmx_localtop_t *gmx_mtop_generate_local_top(const gmx_mtop_t *mtop,
1134 const t_inputrec *ir)
1136 gmx_localtop_t *top;
1140 gen_local_top(mtop, ir, TRUE, top);
1145 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop)
1148 gmx_localtop_t ltop;
1151 gen_local_top(mtop, NULL, FALSE, <op);
1153 top.name = mtop->name;
1154 top.idef = ltop.idef;
1155 top.atomtypes = ltop.atomtypes;
1157 top.excls = ltop.excls;
1158 top.atoms = gmx_mtop_global_atoms(mtop);
1159 top.mols = mtop->mols;
1160 top.symtab = mtop->symtab;
1162 /* We only need to free the moltype and molblock data,
1163 * all other pointers have been copied to top.
1165 * Well, except for the group data, but we can't free those, because they
1166 * are used somewhere even after a call to this function.
1168 for (mt = 0; mt < mtop->nmoltype; mt++)
1170 done_moltype(&mtop->moltype[mt]);
1172 sfree(mtop->moltype);
1174 for (mb = 0; mb < mtop->nmolblock; mb++)
1176 done_molblock(&mtop->molblock[mb]);
1178 sfree(mtop->molblock);