Fixed typos in GB inner dielectric commit for some kernels.
[gromacs.git] / src / gmxlib / mtop_util.c
blobee4e5705eb413eee8874a9060e0c4e5536385cb3
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This file is part of Gromacs Copyright (c) 1991-2008
5 * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
7 * This program is free software; you can redistribute it and/or
8 * modify it under the terms of the GNU General Public License
9 * as published by the Free Software Foundation; either version 2
10 * of the License, or (at your option) any later version.
12 * To help us fund GROMACS development, we humbly ask that you cite
13 * the research papers on the package. Check out http://www.gromacs.org
15 * And Hey:
16 * Gnomes, ROck Monsters And Chili Sauce
19 #ifdef HAVE_CONFIG_H
20 #include <config.h>
21 #endif
23 #include <string.h>
24 #include "smalloc.h"
25 #include "typedefs.h"
26 #include "mtop_util.h"
27 #include "topsort.h"
28 #include "symtab.h"
29 #include "gmx_fatal.h"
31 static int gmx_mtop_maxresnr(const gmx_mtop_t *mtop,int maxres_renum)
33 int maxresnr,mt,r;
34 const t_atoms *atoms;
36 maxresnr = 0;
38 for(mt=0; mt<mtop->nmoltype; mt++)
40 atoms = &mtop->moltype[mt].atoms;
41 if (atoms->nres > maxres_renum)
43 for(r=0; r<atoms->nres; r++)
45 if (atoms->resinfo[r].nr > maxresnr)
47 maxresnr = atoms->resinfo[r].nr;
53 return maxresnr;
56 void gmx_mtop_finalize(gmx_mtop_t *mtop)
58 char *env;
60 mtop->maxres_renum = 1;
62 env = getenv("GMX_MAXRESRENUM");
63 if (env != NULL)
65 sscanf(env,"%d",&mtop->maxres_renum);
67 if (mtop->maxres_renum == -1)
69 /* -1 signals renumber residues in all molecules */
70 mtop->maxres_renum = INT_MAX;
73 mtop->maxresnr = gmx_mtop_maxresnr(mtop,mtop->maxres_renum);
76 int ncg_mtop(const gmx_mtop_t *mtop)
78 int ncg;
79 int mb;
81 ncg = 0;
82 for(mb=0; mb<mtop->nmolblock; mb++)
84 ncg +=
85 mtop->molblock[mb].nmol*
86 mtop->moltype[mtop->molblock[mb].type].cgs.nr;
89 return ncg;
92 void gmx_mtop_atomnr_to_atom(const gmx_mtop_t *mtop,int atnr_global,
93 t_atom **atom)
95 int mb,a_start,a_end,atnr_mol;
97 if (atnr_global < 0 || atnr_global >= mtop->natoms)
99 gmx_fatal(FARGS,"gmx_mtop_atomnr_to_atom was called with atnr_global=%d which is not in the atom range of this system (%d-%d)",
100 atnr_global,0,mtop->natoms-1);
103 mb = -1;
104 a_end = 0;
107 mb++;
108 a_start = a_end;
109 a_end = a_start + mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
111 while (atnr_global >= a_end);
113 atnr_mol = (atnr_global - a_start) % mtop->molblock[mb].natoms_mol;
115 *atom = &mtop->moltype[mtop->molblock[mb].type].atoms.atom[atnr_mol];
118 void gmx_mtop_atomnr_to_ilist(const gmx_mtop_t *mtop,int atnr_global,
119 t_ilist **ilist_mol,int *atnr_offset)
121 int mb,a_start,a_end,atnr_local;
123 if (atnr_global < 0 || atnr_global >= mtop->natoms)
125 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)",
126 atnr_global,0,mtop->natoms-1);
129 mb = -1;
130 a_end = 0;
133 mb++;
134 a_start = a_end;
135 a_end = a_start + mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
137 while (atnr_global >= a_end);
139 *ilist_mol = mtop->moltype[mtop->molblock[mb].type].ilist;
141 atnr_local = (atnr_global - a_start) % mtop->molblock[mb].natoms_mol;
143 *atnr_offset = atnr_global - atnr_local;
146 void gmx_mtop_atomnr_to_molblock_ind(const gmx_mtop_t *mtop,int atnr_global,
147 int *molb,int *molnr,int *atnr_mol)
149 int mb,a_start,a_end;
150 t_atoms *atoms;
152 if (atnr_global < 0 || atnr_global >= mtop->natoms)
154 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)",
155 atnr_global,0,mtop->natoms-1);
158 mb = -1;
159 a_end = 0;
162 mb++;
163 a_start = a_end;
164 a_end = a_start + mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
166 while (atnr_global >= a_end);
168 *molb = mb;
169 *molnr = (atnr_global - a_start) / mtop->molblock[mb].natoms_mol;
170 *atnr_mol = atnr_global - a_start - (*molnr)*mtop->molblock[mb].natoms_mol;
173 void gmx_mtop_atominfo_global(const gmx_mtop_t *mtop,int atnr_global,
174 char **atomname,int *resnr,char **resname)
176 int mb,a_start,a_end,maxresnr,at_loc;
177 gmx_molblock_t *molb;
178 t_atoms *atoms=NULL;
180 if (atnr_global < 0 || atnr_global >= mtop->natoms)
182 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)",
183 atnr_global,0,mtop->natoms-1);
186 mb = -1;
187 a_end = 0;
188 maxresnr = mtop->maxresnr;
191 if (mb >= 0)
193 if (atoms->nres <= mtop->maxres_renum)
195 /* Single residue molecule, keep counting */
196 maxresnr += mtop->molblock[mb].nmol*atoms->nres;
199 mb++;
200 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
201 a_start = a_end;
202 a_end = a_start + mtop->molblock[mb].nmol*atoms->nr;
204 while (atnr_global >= a_end);
206 at_loc = (atnr_global - a_start) % atoms->nr;
207 *atomname = *(atoms->atomname[at_loc]);
208 if (atoms->nres > mtop->maxres_renum)
210 *resnr = atoms->resinfo[atoms->atom[at_loc].resind].nr;
212 else
214 /* Single residue molecule, keep counting */
215 *resnr = maxresnr + 1 + at_loc/atoms->nr*atoms->nres + atoms->atom[at_loc].resind;
217 *resname = *(atoms->resinfo[atoms->atom[at_loc].resind].name);
220 typedef struct gmx_mtop_atomloop_all
222 const gmx_mtop_t *mtop;
223 int mblock;
224 t_atoms *atoms;
225 int mol;
226 int maxresnr;
227 int at_local;
228 int at_global;
229 } t_gmx_mtop_atomloop_all;
231 gmx_mtop_atomloop_all_t
232 gmx_mtop_atomloop_all_init(const gmx_mtop_t *mtop)
234 struct gmx_mtop_atomloop_all *aloop;
236 snew(aloop,1);
238 aloop->mtop = mtop;
239 aloop->mblock = 0;
240 aloop->atoms =
241 &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
242 aloop->mol = 0;
243 aloop->maxresnr = mtop->maxresnr;
244 aloop->at_local = -1;
245 aloop->at_global = -1;
247 return aloop;
250 static void gmx_mtop_atomloop_all_destroy(gmx_mtop_atomloop_all_t aloop)
252 sfree(aloop);
255 bool gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop,
256 int *at_global,t_atom **atom)
258 if (aloop == NULL)
260 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
263 aloop->at_local++;
264 aloop->at_global++;
266 if (aloop->at_local >= aloop->atoms->nr)
268 if (aloop->atoms->nres <= aloop->mtop->maxres_renum)
270 /* Single residue molecule, increase the count with one */
271 aloop->maxresnr += aloop->atoms->nres;
273 aloop->mol++;
274 aloop->at_local = 0;
275 if (aloop->mol >= aloop->mtop->molblock[aloop->mblock].nmol)
277 aloop->mblock++;
278 if (aloop->mblock >= aloop->mtop->nmolblock)
280 gmx_mtop_atomloop_all_destroy(aloop);
281 return FALSE;
283 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
284 aloop->mol = 0;
288 *at_global = aloop->at_global;
289 *atom = &aloop->atoms->atom[aloop->at_local];
291 return TRUE;
294 void gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop,
295 char **atomname,int *resnr,char **resname)
297 int resind_mol;
299 *atomname = *(aloop->atoms->atomname[aloop->at_local]);
300 resind_mol = aloop->atoms->atom[aloop->at_local].resind;
301 *resnr = aloop->atoms->resinfo[resind_mol].nr;
302 if (aloop->atoms->nres <= aloop->mtop->maxres_renum)
304 *resnr = aloop->maxresnr + 1 + resind_mol;
306 *resname = *(aloop->atoms->resinfo[resind_mol].name);
309 void gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop,
310 gmx_moltype_t **moltype,int *at_mol)
312 *moltype = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type];
313 *at_mol = aloop->at_local;
316 typedef struct gmx_mtop_atomloop_block
318 const gmx_mtop_t *mtop;
319 int mblock;
320 t_atoms *atoms;
321 int at_local;
322 } t_gmx_mtop_atomloop_block;
324 gmx_mtop_atomloop_block_t
325 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop)
327 struct gmx_mtop_atomloop_block *aloop;
329 snew(aloop,1);
331 aloop->mtop = mtop;
332 aloop->mblock = 0;
333 aloop->atoms = &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
334 aloop->at_local = -1;
336 return aloop;
339 static void gmx_mtop_atomloop_block_destroy(gmx_mtop_atomloop_block_t aloop)
341 sfree(aloop);
344 bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
345 t_atom **atom,int *nmol)
347 if (aloop == NULL)
349 gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
352 aloop->at_local++;
354 if (aloop->at_local >= aloop->atoms->nr)
356 aloop->mblock++;
357 if (aloop->mblock >= aloop->mtop->nmolblock)
359 gmx_mtop_atomloop_block_destroy(aloop);
360 return FALSE;
362 aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
363 aloop->at_local = 0;
366 *atom = &aloop->atoms->atom[aloop->at_local];
367 *nmol = aloop->mtop->molblock[aloop->mblock].nmol;
369 return TRUE;
372 typedef struct gmx_mtop_ilistloop
374 const gmx_mtop_t *mtop;
375 int mblock;
376 } t_gmx_mtop_ilist;
378 gmx_mtop_ilistloop_t
379 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop)
381 struct gmx_mtop_ilistloop *iloop;
383 snew(iloop,1);
385 iloop->mtop = mtop;
386 iloop->mblock = -1;
388 return iloop;
391 static void gmx_mtop_ilistloop_destroy(gmx_mtop_ilistloop_t iloop)
393 sfree(iloop);
396 bool gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
397 t_ilist **ilist_mol,int *nmol)
399 if (iloop == NULL)
401 gmx_incons("gmx_mtop_ilistloop_next called without calling gmx_mtop_ilistloop_init");
404 iloop->mblock++;
405 if (iloop->mblock == iloop->mtop->nmolblock)
407 gmx_mtop_ilistloop_destroy(iloop);
408 return FALSE;
411 *ilist_mol =
412 iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
414 *nmol = iloop->mtop->molblock[iloop->mblock].nmol;
416 return TRUE;
418 typedef struct gmx_mtop_ilistloop_all
420 const gmx_mtop_t *mtop;
421 int mblock;
422 int mol;
423 int a_offset;
424 } t_gmx_mtop_ilist_all;
426 gmx_mtop_ilistloop_all_t
427 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop)
429 struct gmx_mtop_ilistloop_all *iloop;
431 snew(iloop,1);
433 iloop->mtop = mtop;
434 iloop->mblock = 0;
435 iloop->mol = -1;
436 iloop->a_offset = 0;
438 return iloop;
441 static void gmx_mtop_ilistloop_all_destroy(gmx_mtop_ilistloop_all_t iloop)
443 sfree(iloop);
446 bool gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
447 t_ilist **ilist_mol,int *atnr_offset)
449 gmx_molblock_t *molb;
451 if (iloop == NULL)
453 gmx_incons("gmx_mtop_ilistloop_all_next called without calling gmx_mtop_ilistloop_all_init");
456 if (iloop->mol >= 0)
458 iloop->a_offset += iloop->mtop->molblock[iloop->mblock].natoms_mol;
461 iloop->mol++;
463 if (iloop->mol >= iloop->mtop->molblock[iloop->mblock].nmol) {
464 iloop->mblock++;
465 iloop->mol = 0;
466 if (iloop->mblock == iloop->mtop->nmolblock)
468 gmx_mtop_ilistloop_all_destroy(iloop);
469 return FALSE;
473 *ilist_mol =
474 iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
476 *atnr_offset = iloop->a_offset;
478 return TRUE;
481 int gmx_mtop_ftype_count(const gmx_mtop_t *mtop,int ftype)
483 gmx_mtop_ilistloop_t iloop;
484 t_ilist *il;
485 int n,nmol;
487 n = 0;
489 iloop = gmx_mtop_ilistloop_init(mtop);
490 while (gmx_mtop_ilistloop_next(iloop,&il,&nmol))
492 n += nmol*il[ftype].nr/(1+NRAL(ftype));
495 return n;
498 t_block gmx_mtop_global_cgs(const gmx_mtop_t *mtop)
500 t_block cgs_gl,*cgs_mol;
501 int mb,mol,cg;
502 gmx_molblock_t *molb;
503 t_atoms *atoms;
505 /* In most cases this is too much, but we realloc at the end */
506 snew(cgs_gl.index,mtop->natoms+1);
508 cgs_gl.nr = 0;
509 cgs_gl.index[0] = 0;
510 for(mb=0; mb<mtop->nmolblock; mb++)
512 molb = &mtop->molblock[mb];
513 cgs_mol = &mtop->moltype[molb->type].cgs;
514 for(mol=0; mol<molb->nmol; mol++)
516 for(cg=0; cg<cgs_mol->nr; cg++)
518 cgs_gl.index[cgs_gl.nr+1] =
519 cgs_gl.index[cgs_gl.nr] +
520 cgs_mol->index[cg+1] - cgs_mol->index[cg];
521 cgs_gl.nr++;
525 cgs_gl.nalloc_index = cgs_gl.nr + 1;
526 srenew(cgs_gl.index,cgs_gl.nalloc_index);
528 return cgs_gl;
531 static void atomcat(t_atoms *dest, t_atoms *src, int copies,
532 int maxres_renum, int *maxresnr)
534 int i,j,l,size;
535 int srcnr=src->nr;
536 int destnr=dest->nr;
538 if (srcnr)
540 size=destnr+copies*srcnr;
541 srenew(dest->atom,size);
542 srenew(dest->atomname,size);
543 srenew(dest->atomtype,size);
544 srenew(dest->atomtypeB,size);
546 if (src->nres)
548 size=dest->nres+copies*src->nres;
549 srenew(dest->resinfo,size);
552 /* residue information */
553 for (l=dest->nres,j=0; (j<copies); j++,l+=src->nres)
555 memcpy((char *) &(dest->resinfo[l]),(char *) &(src->resinfo[0]),
556 (size_t)(src->nres*sizeof(src->resinfo[0])));
559 for (l=destnr,j=0; (j<copies); j++,l+=srcnr)
561 memcpy((char *) &(dest->atomname[l]),(char *) &(src->atomname[0]),
562 (size_t)(srcnr*sizeof(src->atomname[0])));
563 memcpy((char *) &(dest->atomtype[l]),(char *) &(src->atomtype[0]),
564 (size_t)(srcnr*sizeof(src->atomtype[0])));
565 memcpy((char *) &(dest->atomtypeB[l]),(char *) &(src->atomtypeB[0]),
566 (size_t)(srcnr*sizeof(src->atomtypeB[0])));
567 memcpy((char *) &(dest->atom[l]),(char *) &(src->atom[0]),
568 (size_t)(srcnr*sizeof(src->atom[0])));
571 /* Increment residue indices */
572 for (l=destnr,j=0; (j<copies); j++)
574 for (i=0; (i<srcnr); i++,l++)
576 dest->atom[l].resind = dest->nres+j*src->nres+src->atom[i].resind;
580 if (src->nres <= maxres_renum)
582 /* Single residue molecule, continue counting residues */
583 for (j=0; (j<copies); j++)
585 for (l=0; l<src->nres; l++)
587 (*maxresnr)++;
588 dest->resinfo[dest->nres+j*src->nres+l].nr = *maxresnr;
593 dest->nres += copies*src->nres;
594 dest->nr += copies*src->nr;
597 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t *mtop)
599 t_atoms atoms;
600 int maxresnr,mb;
601 gmx_molblock_t *molb;
603 init_t_atoms(&atoms,0,FALSE);
605 maxresnr = mtop->maxresnr;
606 for(mb=0; mb<mtop->nmolblock; mb++)
608 molb = &mtop->molblock[mb];
609 atomcat(&atoms,&mtop->moltype[molb->type].atoms,molb->nmol,
610 mtop->maxres_renum,&maxresnr);
613 return atoms;
616 void gmx_mtop_make_atomic_charge_groups(gmx_mtop_t *mtop,
617 bool bKeepSingleMolCG)
619 int mb,cg;
620 t_block *cgs_mol;
622 for(mb=0; mb<mtop->nmolblock; mb++)
624 cgs_mol = &mtop->moltype[mtop->molblock[mb].type].cgs;
625 if (!(bKeepSingleMolCG && cgs_mol->nr == 1))
627 cgs_mol->nr = mtop->molblock[mb].natoms_mol;
628 cgs_mol->nalloc_index = cgs_mol->nr + 1;
629 srenew(cgs_mol->index,cgs_mol->nalloc_index);
630 for(cg=0; cg<cgs_mol->nr+1; cg++)
632 cgs_mol->index[cg] = cg;
639 * The cat routines below are old code from src/kernel/topcat.c
642 static void blockcat(t_block *dest,t_block *src,int copies,
643 int dnum,int snum)
645 int i,j,l,nra,size;
647 if (src->nr)
649 size=(dest->nr+copies*src->nr+1);
650 srenew(dest->index,size);
653 nra = dest->index[dest->nr];
654 for (l=dest->nr,j=0; (j<copies); j++)
656 for (i=0; (i<src->nr); i++)
658 dest->index[l++] = nra + src->index[i];
660 nra += src->index[src->nr];
662 dest->nr += copies*src->nr;
663 dest->index[dest->nr] = nra;
666 static void blockacat(t_blocka *dest,t_blocka *src,int copies,
667 int dnum,int snum)
669 int i,j,l,size;
670 int destnr = dest->nr;
671 int destnra = dest->nra;
673 if (src->nr)
675 size=(dest->nr+copies*src->nr+1);
676 srenew(dest->index,size);
678 if (src->nra)
680 size=(dest->nra+copies*src->nra);
681 srenew(dest->a,size);
684 for (l=destnr,j=0; (j<copies); j++)
686 for (i=0; (i<src->nr); i++)
688 dest->index[l++] = dest->nra+src->index[i];
690 dest->nra += src->nra;
692 for (l=destnra,j=0; (j<copies); j++)
694 for (i=0; (i<src->nra); i++)
696 dest->a[l++] = dnum+src->a[i];
698 dnum+=snum;
699 dest->nr += src->nr;
701 dest->index[dest->nr] = dest->nra;
704 static void ilistcat(int ftype,t_ilist *dest,t_ilist *src,int copies,
705 int dnum,int snum)
707 int nral,c,i,a;
709 nral = NRAL(ftype);
711 dest->nalloc = dest->nr + copies*src->nr;
712 srenew(dest->iatoms,dest->nalloc);
714 for(c=0; c<copies; c++)
716 for(i=0; i<src->nr; )
718 dest->iatoms[dest->nr++] = src->iatoms[i++];
719 for(a=0; a<nral; a++)
721 dest->iatoms[dest->nr++] = dnum + src->iatoms[i++];
724 dnum += snum;
728 static void set_posres_params(t_idef *idef,gmx_molblock_t *molb,
729 int i0,int a_offset)
731 t_ilist *il;
732 int i1,i,a_molb;
733 t_iparams *ip;
735 il = &idef->il[F_POSRES];
736 i1 = il->nr/2;
737 idef->iparams_posres_nalloc = i1;
738 srenew(idef->iparams_posres,idef->iparams_posres_nalloc);
739 for(i=i0; i<i1; i++)
741 ip = &idef->iparams_posres[i];
742 /* Copy the force constants */
743 *ip = idef->iparams[il->iatoms[i*2]];
744 a_molb = il->iatoms[i*2+1] - a_offset;
745 if (molb->nposres_xA == 0)
747 gmx_incons("Position restraint coordinates are missing");
749 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
750 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
751 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
752 if (molb->nposres_xB > 0)
754 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
755 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
756 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
758 else
760 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
761 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
762 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
764 /* Set the parameter index for idef->iparams_posre */
765 il->iatoms[i*2] = i;
769 static void gen_local_top(const gmx_mtop_t *mtop,const t_inputrec *ir,
770 bool bMergeConstr,
771 gmx_localtop_t *top)
773 int mb,srcnr,destnr,ftype,ftype_dest,mt,natoms,mol,nposre_old;
774 gmx_molblock_t *molb;
775 gmx_moltype_t *molt;
776 const gmx_ffparams_t *ffp;
777 t_idef *idef;
778 real *qA,*qB;
779 gmx_mtop_atomloop_all_t aloop;
780 int ag;
781 t_atom *atom;
783 top->atomtypes = mtop->atomtypes;
785 ffp = &mtop->ffparams;
787 idef = &top->idef;
788 idef->ntypes = ffp->ntypes;
789 idef->atnr = ffp->atnr;
790 idef->functype = ffp->functype;
791 idef->iparams = ffp->iparams;
792 idef->iparams_posres = NULL;
793 idef->iparams_posres_nalloc = 0;
794 idef->fudgeQQ = ffp->fudgeQQ;
795 idef->cmap_grid = ffp->cmap_grid;
796 idef->ilsort = ilsortUNKNOWN;
798 init_block(&top->cgs);
799 init_blocka(&top->excls);
800 for(ftype=0; ftype<F_NRE; ftype++)
802 idef->il[ftype].nr = 0;
803 idef->il[ftype].nalloc = 0;
804 idef->il[ftype].iatoms = NULL;
807 natoms = 0;
808 for(mb=0; mb<mtop->nmolblock; mb++)
810 molb = &mtop->molblock[mb];
811 molt = &mtop->moltype[molb->type];
813 srcnr = molt->atoms.nr;
814 destnr = natoms;
816 blockcat(&top->cgs,&molt->cgs,molb->nmol,destnr,srcnr);
818 blockacat(&top->excls,&molt->excls,molb->nmol,destnr,srcnr);
820 nposre_old = idef->il[F_POSRES].nr;
821 for(ftype=0; ftype<F_NRE; ftype++)
823 if (bMergeConstr &&
824 ftype == F_CONSTR && molt->ilist[F_CONSTRNC].nr > 0)
826 /* Merge all constrains into one ilist.
827 * This simplifies the constraint code.
829 for(mol=0; mol<molb->nmol; mol++) {
830 ilistcat(ftype,&idef->il[F_CONSTR],&molt->ilist[F_CONSTR],
831 1,destnr+mol*srcnr,srcnr);
832 ilistcat(ftype,&idef->il[F_CONSTR],&molt->ilist[F_CONSTRNC],
833 1,destnr+mol*srcnr,srcnr);
836 else if (!(bMergeConstr && ftype == F_CONSTRNC))
838 ilistcat(ftype,&idef->il[ftype],&molt->ilist[ftype],
839 molb->nmol,destnr,srcnr);
842 if (idef->il[F_POSRES].nr > nposre_old)
844 set_posres_params(idef,molb,nposre_old/2,natoms);
847 natoms += molb->nmol*srcnr;
850 if (ir == NULL)
852 top->idef.ilsort = ilsortUNKNOWN;
854 else
856 if (ir->efep != efepNO && gmx_mtop_bondeds_free_energy(mtop))
858 snew(qA,mtop->natoms);
859 snew(qB,mtop->natoms);
860 aloop = gmx_mtop_atomloop_all_init(mtop);
861 while (gmx_mtop_atomloop_all_next(aloop,&ag,&atom))
863 qA[ag] = atom->q;
864 qB[ag] = atom->qB;
866 gmx_sort_ilist_fe(&top->idef,qA,qB);
867 sfree(qA);
868 sfree(qB);
870 else
872 top->idef.ilsort = ilsortNO_FE;
877 gmx_localtop_t *gmx_mtop_generate_local_top(const gmx_mtop_t *mtop,
878 const t_inputrec *ir)
880 gmx_localtop_t *top;
882 snew(top,1);
884 gen_local_top(mtop,ir,TRUE,top);
886 return top;
889 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop)
891 int mt,mb;
892 gmx_localtop_t ltop;
893 t_topology top;
895 gen_local_top(mtop,NULL,FALSE,&ltop);
897 open_symtab(&top.symtab);
899 open_symtab(&top.symtab);
901 top.name = mtop->name;
902 top.idef = ltop.idef;
903 top.atomtypes = ltop.atomtypes;
904 top.cgs = ltop.cgs;
905 top.excls = ltop.excls;
906 top.atoms = gmx_mtop_global_atoms(mtop);
907 top.mols = mtop->mols;
909 /* We only need to free the moltype and molblock data,
910 * all other pointers have been copied to top.
912 for(mt=0; mt<mtop->nmoltype; mt++)
914 done_moltype(&mtop->moltype[mt]);
916 sfree(mtop->moltype);
918 for(mb=0; mb<mtop->nmolblock; mb++)
920 done_molblock(&mtop->molblock[mb]);
922 sfree(mtop->molblock);
924 return top;