19fd0927dcb472479cad967a581139c2ada0ad9d
[alexxy/gromacs.git] / src / gromacs / topology / mtop_util.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "mtop_util.h"
38
39 #include <limits.h>
40 #include <stddef.h>
41 #include <stdio.h>
42 #include <stdlib.h>
43 #include <string.h>
44
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"
58
59 static int gmx_mtop_maxresnr(const gmx_mtop_t *mtop, int maxres_renum)
60 {
61     int            maxresnr, mt, r;
62     const t_atoms *atoms;
63
64     maxresnr = 0;
65
66     for (mt = 0; mt < mtop->nmoltype; mt++)
67     {
68         atoms = &mtop->moltype[mt].atoms;
69         if (atoms->nres > maxres_renum)
70         {
71             for (r = 0; r < atoms->nres; r++)
72             {
73                 if (atoms->resinfo[r].nr > maxresnr)
74                 {
75                     maxresnr = atoms->resinfo[r].nr;
76                 }
77             }
78         }
79     }
80
81     return maxresnr;
82 }
83
84 void gmx_mtop_finalize(gmx_mtop_t *mtop)
85 {
86     char *env;
87
88     mtop->maxres_renum = 1;
89
90     env = getenv("GMX_MAXRESRENUM");
91     if (env != NULL)
92     {
93         sscanf(env, "%d", &mtop->maxres_renum);
94     }
95     if (mtop->maxres_renum == -1)
96     {
97         /* -1 signals renumber residues in all molecules */
98         mtop->maxres_renum = INT_MAX;
99     }
100
101     mtop->maxresnr = gmx_mtop_maxresnr(mtop, mtop->maxres_renum);
102 }
103
104 void gmx_mtop_count_atomtypes(const gmx_mtop_t *mtop, int state, int typecount[])
105 {
106     int      i, mb, nmol, tpi;
107     t_atoms *atoms;
108
109     for (i = 0; i < mtop->ffparams.atnr; ++i)
110     {
111         typecount[i] = 0;
112     }
113     for (mb = 0; mb < mtop->nmolblock; ++mb)
114     {
115         nmol  = mtop->molblock[mb].nmol;
116         atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
117         for (i = 0; i < atoms->nr; ++i)
118         {
119             if (state == 0)
120             {
121                 tpi = atoms->atom[i].type;
122             }
123             else
124             {
125                 tpi = atoms->atom[i].typeB;
126             }
127             typecount[tpi] += nmol;
128         }
129     }
130 }
131
132 int ncg_mtop(const gmx_mtop_t *mtop)
133 {
134     int ncg;
135     int mb;
136
137     ncg = 0;
138     for (mb = 0; mb < mtop->nmolblock; mb++)
139     {
140         ncg +=
141             mtop->molblock[mb].nmol*
142             mtop->moltype[mtop->molblock[mb].type].cgs.nr;
143     }
144
145     return ncg;
146 }
147
148 void gmx_mtop_remove_chargegroups(gmx_mtop_t *mtop)
149 {
150     int      mt;
151     t_block *cgs;
152     int      i;
153
154     for (mt = 0; mt < mtop->nmoltype; mt++)
155     {
156         cgs = &mtop->moltype[mt].cgs;
157         if (cgs->nr < mtop->moltype[mt].atoms.nr)
158         {
159             cgs->nr = mtop->moltype[mt].atoms.nr;
160             srenew(cgs->index, cgs->nr+1);
161             for (i = 0; i < cgs->nr+1; i++)
162             {
163                 cgs->index[i] = i;
164             }
165         }
166     }
167 }
168
169
170 typedef struct
171 {
172     int a_start;
173     int a_end;
174     int na_mol;
175 } mb_at_t;
176
177 typedef struct gmx_mtop_atomlookup
178 {
179     const gmx_mtop_t *mtop;
180     int               nmb;
181     int               mb_start;
182     mb_at_t          *mba;
183 } t_gmx_mtop_atomlookup;
184
185
186 gmx_mtop_atomlookup_t
187 gmx_mtop_atomlookup_init(const gmx_mtop_t *mtop)
188 {
189     t_gmx_mtop_atomlookup *alook;
190     int                    mb;
191     int                    a_start, a_end, na, na_start = -1;
192
193     snew(alook, 1);
194
195     alook->mtop     = mtop;
196     alook->nmb      = mtop->nmolblock;
197     alook->mb_start = 0;
198     snew(alook->mba, alook->nmb);
199
200     a_start = 0;
201     for (mb = 0; mb < mtop->nmolblock; mb++)
202     {
203         na    = mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
204         a_end = a_start + na;
205
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;
209
210         /* We start the binary search with the largest block */
211         if (mb == 0 || na > na_start)
212         {
213             alook->mb_start = mb;
214             na_start        = na;
215         }
216
217         a_start = a_end;
218     }
219
220     return alook;
221 }
222
223 gmx_mtop_atomlookup_t
224 gmx_mtop_atomlookup_settle_init(const gmx_mtop_t *mtop)
225 {
226     t_gmx_mtop_atomlookup *alook;
227     int                    mb;
228     int                    na, na_start = -1;
229
230     alook = gmx_mtop_atomlookup_init(mtop);
231
232     /* Check if the starting molblock has settle */
233     if (mtop->moltype[mtop->molblock[alook->mb_start].type].ilist[F_SETTLE].nr  == 0)
234     {
235         /* Search the largest molblock with settle */
236         alook->mb_start = -1;
237         for (mb = 0; mb < mtop->nmolblock; mb++)
238         {
239             if (mtop->moltype[mtop->molblock[mb].type].ilist[F_SETTLE].nr > 0)
240             {
241                 na = alook->mba[mb].a_end - alook->mba[mb].a_start;
242                 if (alook->mb_start == -1 || na > na_start)
243                 {
244                     alook->mb_start = mb;
245                     na_start        = na;
246                 }
247             }
248         }
249
250         if (alook->mb_start == -1)
251         {
252             gmx_incons("gmx_mtop_atomlookup_settle_init called without settles");
253         }
254     }
255
256     return alook;
257 }
258
259 void
260 gmx_mtop_atomlookup_destroy(gmx_mtop_atomlookup_t alook)
261 {
262     sfree(alook->mba);
263     sfree(alook);
264 }
265
266 void gmx_mtop_atomnr_to_atom(const gmx_mtop_atomlookup_t alook,
267                              int                         atnr_global,
268                              t_atom                    **atom)
269 {
270     int mb0, mb1, mb;
271     int a_start, atnr_mol;
272
273 #ifdef DEBUG_MTOP
274     if (atnr_global < 0 || atnr_global >= mtop->natoms)
275     {
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);
278     }
279 #endif
280
281     mb0 = -1;
282     mb1 = alook->nmb;
283     mb  = alook->mb_start;
284
285     while (TRUE)
286     {
287         a_start = alook->mba[mb].a_start;
288         if (atnr_global < a_start)
289         {
290             mb1 = mb;
291         }
292         else if (atnr_global >= alook->mba[mb].a_end)
293         {
294             mb0 = mb;
295         }
296         else
297         {
298             break;
299         }
300         mb = ((mb0 + mb1 + 1)>>1);
301     }
302
303     atnr_mol = (atnr_global - a_start) % alook->mba[mb].na_mol;
304
305     *atom = &alook->mtop->moltype[alook->mtop->molblock[mb].type].atoms.atom[atnr_mol];
306 }
307
308 void gmx_mtop_atomnr_to_ilist(const gmx_mtop_atomlookup_t alook,
309                               int atnr_global,
310                               t_ilist **ilist_mol, int *atnr_offset)
311 {
312     int mb0, mb1, mb;
313     int a_start, atnr_local;
314
315 #ifdef DEBUG_MTOP
316     if (atnr_global < 0 || atnr_global >= mtop->natoms)
317     {
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);
320     }
321 #endif
322
323     mb0 = -1;
324     mb1 = alook->nmb;
325     mb  = alook->mb_start;
326
327     while (TRUE)
328     {
329         a_start = alook->mba[mb].a_start;
330         if (atnr_global < a_start)
331         {
332             mb1 = mb;
333         }
334         else if (atnr_global >= alook->mba[mb].a_end)
335         {
336             mb0 = mb;
337         }
338         else
339         {
340             break;
341         }
342         mb = ((mb0 + mb1 + 1)>>1);
343     }
344
345     *ilist_mol = alook->mtop->moltype[alook->mtop->molblock[mb].type].ilist;
346
347     atnr_local = (atnr_global - a_start) % alook->mba[mb].na_mol;
348
349     *atnr_offset = atnr_global - atnr_local;
350 }
351
352 void gmx_mtop_atomnr_to_molblock_ind(const gmx_mtop_atomlookup_t alook,
353                                      int atnr_global,
354                                      int *molb, int *molnr, int *atnr_mol)
355 {
356     int mb0, mb1, mb;
357     int a_start;
358
359 #ifdef DEBUG_MTOP
360     if (atnr_global < 0 || atnr_global >= mtop->natoms)
361     {
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);
364     }
365 #endif
366
367     mb0 = -1;
368     mb1 = alook->nmb;
369     mb  = alook->mb_start;
370
371     while (TRUE)
372     {
373         a_start = alook->mba[mb].a_start;
374         if (atnr_global < a_start)
375         {
376             mb1 = mb;
377         }
378         else if (atnr_global >= alook->mba[mb].a_end)
379         {
380             mb0 = mb;
381         }
382         else
383         {
384             break;
385         }
386         mb = ((mb0 + mb1 + 1)>>1);
387     }
388
389     *molb     = mb;
390     *molnr    = (atnr_global - a_start) / alook->mba[mb].na_mol;
391     *atnr_mol = atnr_global - a_start - (*molnr)*alook->mba[mb].na_mol;
392 }
393
394 void gmx_mtop_atominfo_global(const gmx_mtop_t *mtop, int atnr_global,
395                               char **atomname, int *resnr, char **resname)
396 {
397     int             mb, a_start, a_end, maxresnr, at_loc;
398     gmx_molblock_t *molb;
399     t_atoms        *atoms = NULL;
400
401     if (atnr_global < 0 || atnr_global >= mtop->natoms)
402     {
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);
405     }
406
407     mb       = -1;
408     a_end    = 0;
409     maxresnr = mtop->maxresnr;
410     do
411     {
412         if (mb >= 0)
413         {
414             if (atoms->nres <= mtop->maxres_renum)
415             {
416                 /* Single residue molecule, keep counting */
417                 maxresnr += mtop->molblock[mb].nmol*atoms->nres;
418             }
419         }
420         mb++;
421         atoms   = &mtop->moltype[mtop->molblock[mb].type].atoms;
422         a_start = a_end;
423         a_end   = a_start + mtop->molblock[mb].nmol*atoms->nr;
424     }
425     while (atnr_global >= a_end);
426
427     at_loc    = (atnr_global - a_start) % atoms->nr;
428     *atomname = *(atoms->atomname[at_loc]);
429     if (atoms->nres > mtop->maxres_renum)
430     {
431         *resnr = atoms->resinfo[atoms->atom[at_loc].resind].nr;
432     }
433     else
434     {
435         /* Single residue molecule, keep counting */
436         *resnr = maxresnr + 1 + (atnr_global - a_start)/atoms->nr*atoms->nres + atoms->atom[at_loc].resind;
437     }
438     *resname  = *(atoms->resinfo[atoms->atom[at_loc].resind].name);
439 }
440
441 typedef struct gmx_mtop_atomloop_all
442 {
443     const gmx_mtop_t *mtop;
444     int               mblock;
445     t_atoms          *atoms;
446     int               mol;
447     int               maxresnr;
448     int               at_local;
449     int               at_global;
450 } t_gmx_mtop_atomloop_all;
451
452 gmx_mtop_atomloop_all_t
453 gmx_mtop_atomloop_all_init(const gmx_mtop_t *mtop)
454 {
455     struct gmx_mtop_atomloop_all *aloop;
456
457     snew(aloop, 1);
458
459     aloop->mtop         = mtop;
460     aloop->mblock       = 0;
461     aloop->atoms        =
462         &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
463     aloop->mol          = 0;
464     aloop->maxresnr     = mtop->maxresnr;
465     aloop->at_local     = -1;
466     aloop->at_global    = -1;
467
468     return aloop;
469 }
470
471 static void gmx_mtop_atomloop_all_destroy(gmx_mtop_atomloop_all_t aloop)
472 {
473     sfree(aloop);
474 }
475
476 gmx_bool gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop,
477                                     int *at_global, t_atom **atom)
478 {
479     if (aloop == NULL)
480     {
481         gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
482     }
483
484     aloop->at_local++;
485     aloop->at_global++;
486
487     if (aloop->at_local >= aloop->atoms->nr)
488     {
489         if (aloop->atoms->nres <= aloop->mtop->maxres_renum)
490         {
491             /* Single residue molecule, increase the count with one */
492             aloop->maxresnr += aloop->atoms->nres;
493         }
494         aloop->mol++;
495         aloop->at_local = 0;
496         if (aloop->mol >= aloop->mtop->molblock[aloop->mblock].nmol)
497         {
498             aloop->mblock++;
499             if (aloop->mblock >= aloop->mtop->nmolblock)
500             {
501                 gmx_mtop_atomloop_all_destroy(aloop);
502                 return FALSE;
503             }
504             aloop->atoms = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
505             aloop->mol   = 0;
506         }
507     }
508
509     *at_global = aloop->at_global;
510     *atom      = &aloop->atoms->atom[aloop->at_local];
511
512     return TRUE;
513 }
514
515 void gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop,
516                                  char **atomname, int *resnr, char **resname)
517 {
518     int resind_mol;
519
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)
524     {
525         *resnr = aloop->maxresnr + 1 + resind_mol;
526     }
527     *resname  = *(aloop->atoms->resinfo[resind_mol].name);
528 }
529
530 void gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop,
531                                    gmx_moltype_t **moltype, int *at_mol)
532 {
533     *moltype = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type];
534     *at_mol  = aloop->at_local;
535 }
536
537 typedef struct gmx_mtop_atomloop_block
538 {
539     const gmx_mtop_t *mtop;
540     int               mblock;
541     t_atoms          *atoms;
542     int               at_local;
543 } t_gmx_mtop_atomloop_block;
544
545 gmx_mtop_atomloop_block_t
546 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop)
547 {
548     struct gmx_mtop_atomloop_block *aloop;
549
550     snew(aloop, 1);
551
552     aloop->mtop      = mtop;
553     aloop->mblock    = 0;
554     aloop->atoms     = &mtop->moltype[mtop->molblock[aloop->mblock].type].atoms;
555     aloop->at_local  = -1;
556
557     return aloop;
558 }
559
560 static void gmx_mtop_atomloop_block_destroy(gmx_mtop_atomloop_block_t aloop)
561 {
562     sfree(aloop);
563 }
564
565 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
566                                       t_atom **atom, int *nmol)
567 {
568     if (aloop == NULL)
569     {
570         gmx_incons("gmx_mtop_atomloop_all_next called without calling gmx_mtop_atomloop_all_init");
571     }
572
573     aloop->at_local++;
574
575     if (aloop->at_local >= aloop->atoms->nr)
576     {
577         aloop->mblock++;
578         if (aloop->mblock >= aloop->mtop->nmolblock)
579         {
580             gmx_mtop_atomloop_block_destroy(aloop);
581             return FALSE;
582         }
583         aloop->atoms    = &aloop->mtop->moltype[aloop->mtop->molblock[aloop->mblock].type].atoms;
584         aloop->at_local = 0;
585     }
586
587     *atom = &aloop->atoms->atom[aloop->at_local];
588     *nmol = aloop->mtop->molblock[aloop->mblock].nmol;
589
590     return TRUE;
591 }
592
593 typedef struct gmx_mtop_ilistloop
594 {
595     const gmx_mtop_t *mtop;
596     int               mblock;
597 } t_gmx_mtop_ilist;
598
599 gmx_mtop_ilistloop_t
600 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop)
601 {
602     struct gmx_mtop_ilistloop *iloop;
603
604     snew(iloop, 1);
605
606     iloop->mtop      = mtop;
607     iloop->mblock    = -1;
608
609     return iloop;
610 }
611
612 static void gmx_mtop_ilistloop_destroy(gmx_mtop_ilistloop_t iloop)
613 {
614     sfree(iloop);
615 }
616
617 gmx_bool gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
618                                  t_ilist **ilist_mol, int *nmol)
619 {
620     if (iloop == NULL)
621     {
622         gmx_incons("gmx_mtop_ilistloop_next called without calling gmx_mtop_ilistloop_init");
623     }
624
625     iloop->mblock++;
626     if (iloop->mblock == iloop->mtop->nmolblock)
627     {
628         gmx_mtop_ilistloop_destroy(iloop);
629         return FALSE;
630     }
631
632     *ilist_mol =
633         iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
634
635     *nmol = iloop->mtop->molblock[iloop->mblock].nmol;
636
637     return TRUE;
638 }
639 typedef struct gmx_mtop_ilistloop_all
640 {
641     const gmx_mtop_t *mtop;
642     int               mblock;
643     int               mol;
644     int               a_offset;
645 } t_gmx_mtop_ilist_all;
646
647 gmx_mtop_ilistloop_all_t
648 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop)
649 {
650     struct gmx_mtop_ilistloop_all *iloop;
651
652     snew(iloop, 1);
653
654     iloop->mtop      = mtop;
655     iloop->mblock    = 0;
656     iloop->mol       = -1;
657     iloop->a_offset  = 0;
658
659     return iloop;
660 }
661
662 static void gmx_mtop_ilistloop_all_destroy(gmx_mtop_ilistloop_all_t iloop)
663 {
664     sfree(iloop);
665 }
666
667 gmx_bool gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
668                                      t_ilist **ilist_mol, int *atnr_offset)
669 {
670     gmx_molblock_t *molb;
671
672     if (iloop == NULL)
673     {
674         gmx_incons("gmx_mtop_ilistloop_all_next called without calling gmx_mtop_ilistloop_all_init");
675     }
676
677     if (iloop->mol >= 0)
678     {
679         iloop->a_offset += iloop->mtop->molblock[iloop->mblock].natoms_mol;
680     }
681
682     iloop->mol++;
683
684     if (iloop->mol >= iloop->mtop->molblock[iloop->mblock].nmol)
685     {
686         iloop->mblock++;
687         iloop->mol = 0;
688         if (iloop->mblock == iloop->mtop->nmolblock)
689         {
690             gmx_mtop_ilistloop_all_destroy(iloop);
691             return FALSE;
692         }
693     }
694
695     *ilist_mol =
696         iloop->mtop->moltype[iloop->mtop->molblock[iloop->mblock].type].ilist;
697
698     *atnr_offset = iloop->a_offset;
699
700     return TRUE;
701 }
702
703 int gmx_mtop_ftype_count(const gmx_mtop_t *mtop, int ftype)
704 {
705     gmx_mtop_ilistloop_t iloop;
706     t_ilist             *il;
707     int                  n, nmol;
708
709     n = 0;
710
711     iloop = gmx_mtop_ilistloop_init(mtop);
712     while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
713     {
714         n += nmol*il[ftype].nr/(1+NRAL(ftype));
715     }
716
717     return n;
718 }
719
720 t_block gmx_mtop_global_cgs(const gmx_mtop_t *mtop)
721 {
722     t_block         cgs_gl, *cgs_mol;
723     int             mb, mol, cg;
724     gmx_molblock_t *molb;
725     t_atoms        *atoms;
726
727     /* In most cases this is too much, but we realloc at the end */
728     snew(cgs_gl.index, mtop->natoms+1);
729
730     cgs_gl.nr       = 0;
731     cgs_gl.index[0] = 0;
732     for (mb = 0; mb < mtop->nmolblock; mb++)
733     {
734         molb    = &mtop->molblock[mb];
735         cgs_mol = &mtop->moltype[molb->type].cgs;
736         for (mol = 0; mol < molb->nmol; mol++)
737         {
738             for (cg = 0; cg < cgs_mol->nr; cg++)
739             {
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];
743                 cgs_gl.nr++;
744             }
745         }
746     }
747     cgs_gl.nalloc_index = cgs_gl.nr + 1;
748     srenew(cgs_gl.index, cgs_gl.nalloc_index);
749
750     return cgs_gl;
751 }
752
753 static void atomcat(t_atoms *dest, t_atoms *src, int copies,
754                     int maxres_renum, int *maxresnr)
755 {
756     int i, j, l, size;
757     int srcnr  = src->nr;
758     int destnr = dest->nr;
759
760     if (srcnr)
761     {
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);
767     }
768     if (src->nres)
769     {
770         size = dest->nres+copies*src->nres;
771         srenew(dest->resinfo, size);
772     }
773
774     /* residue information */
775     for (l = dest->nres, j = 0; (j < copies); j++, l += src->nres)
776     {
777         memcpy((char *) &(dest->resinfo[l]), (char *) &(src->resinfo[0]),
778                (size_t)(src->nres*sizeof(src->resinfo[0])));
779     }
780
781     for (l = destnr, j = 0; (j < copies); j++, l += srcnr)
782     {
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])));
791     }
792
793     /* Increment residue indices */
794     for (l = destnr, j = 0; (j < copies); j++)
795     {
796         for (i = 0; (i < srcnr); i++, l++)
797         {
798             dest->atom[l].resind = dest->nres+j*src->nres+src->atom[i].resind;
799         }
800     }
801
802     if (src->nres <= maxres_renum)
803     {
804         /* Single residue molecule, continue counting residues */
805         for (j = 0; (j < copies); j++)
806         {
807             for (l = 0; l < src->nres; l++)
808             {
809                 (*maxresnr)++;
810                 dest->resinfo[dest->nres+j*src->nres+l].nr = *maxresnr;
811             }
812         }
813     }
814
815     dest->nres += copies*src->nres;
816     dest->nr   += copies*src->nr;
817 }
818
819 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t *mtop)
820 {
821     t_atoms         atoms;
822     int             maxresnr, mb;
823     gmx_molblock_t *molb;
824
825     init_t_atoms(&atoms, 0, FALSE);
826
827     maxresnr = mtop->maxresnr;
828     for (mb = 0; mb < mtop->nmolblock; mb++)
829     {
830         molb = &mtop->molblock[mb];
831         atomcat(&atoms, &mtop->moltype[molb->type].atoms, molb->nmol,
832                 mtop->maxres_renum, &maxresnr);
833     }
834
835     return atoms;
836 }
837
838 void gmx_mtop_make_atomic_charge_groups(gmx_mtop_t *mtop,
839                                         gmx_bool    bKeepSingleMolCG)
840 {
841     int      mb, cg;
842     t_block *cgs_mol;
843
844     for (mb = 0; mb < mtop->nmolblock; mb++)
845     {
846         cgs_mol = &mtop->moltype[mtop->molblock[mb].type].cgs;
847         if (!(bKeepSingleMolCG && cgs_mol->nr == 1))
848         {
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++)
853             {
854                 cgs_mol->index[cg] = cg;
855             }
856         }
857     }
858 }
859
860 /*
861  * The cat routines below are old code from src/kernel/topcat.c
862  */
863
864 static void blockcat(t_block *dest, t_block *src, int copies)
865 {
866     int i, j, l, nra, size;
867
868     if (src->nr)
869     {
870         size = (dest->nr+copies*src->nr+1);
871         srenew(dest->index, size);
872     }
873
874     nra = dest->index[dest->nr];
875     for (l = dest->nr, j = 0; (j < copies); j++)
876     {
877         for (i = 0; (i < src->nr); i++)
878         {
879             dest->index[l++] = nra + src->index[i];
880         }
881         nra += src->index[src->nr];
882     }
883     dest->nr             += copies*src->nr;
884     dest->index[dest->nr] = nra;
885 }
886
887 static void blockacat(t_blocka *dest, t_blocka *src, int copies,
888                       int dnum, int snum)
889 {
890     int i, j, l, size;
891     int destnr  = dest->nr;
892     int destnra = dest->nra;
893
894     if (src->nr)
895     {
896         size = (dest->nr+copies*src->nr+1);
897         srenew(dest->index, size);
898     }
899     if (src->nra)
900     {
901         size = (dest->nra+copies*src->nra);
902         srenew(dest->a, size);
903     }
904
905     for (l = destnr, j = 0; (j < copies); j++)
906     {
907         for (i = 0; (i < src->nr); i++)
908         {
909             dest->index[l++] = dest->nra+src->index[i];
910         }
911         dest->nra += src->nra;
912     }
913     for (l = destnra, j = 0; (j < copies); j++)
914     {
915         for (i = 0; (i < src->nra); i++)
916         {
917             dest->a[l++] = dnum+src->a[i];
918         }
919         dnum     += snum;
920         dest->nr += src->nr;
921     }
922     dest->index[dest->nr] = dest->nra;
923 }
924
925 static void ilistcat(int ftype, t_ilist *dest, t_ilist *src, int copies,
926                      int dnum, int snum)
927 {
928     int nral, c, i, a;
929
930     nral = NRAL(ftype);
931
932     dest->nalloc = dest->nr + copies*src->nr;
933     srenew(dest->iatoms, dest->nalloc);
934
935     for (c = 0; c < copies; c++)
936     {
937         for (i = 0; i < src->nr; )
938         {
939             dest->iatoms[dest->nr++] = src->iatoms[i++];
940             for (a = 0; a < nral; a++)
941             {
942                 dest->iatoms[dest->nr++] = dnum + src->iatoms[i++];
943             }
944         }
945         dnum += snum;
946     }
947 }
948
949 static void set_posres_params(t_idef *idef, gmx_molblock_t *molb,
950                               int i0, int a_offset)
951 {
952     t_ilist   *il;
953     int        i1, i, a_molb;
954     t_iparams *ip;
955
956     il = &idef->il[F_POSRES];
957     i1 = il->nr/2;
958     idef->iparams_posres_nalloc = i1;
959     srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
960     for (i = i0; i < i1; i++)
961     {
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)
967         {
968             gmx_incons("Position restraint coordinates are missing");
969         }
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)
974         {
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];
978         }
979         else
980         {
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];
984         }
985         /* Set the parameter index for idef->iparams_posre */
986         il->iatoms[i*2] = i;
987     }
988 }
989
990 static void set_fbposres_params(t_idef *idef, gmx_molblock_t *molb,
991                                 int i0, int a_offset)
992 {
993     t_ilist   *il;
994     int        i1, i, a_molb;
995     t_iparams *ip;
996
997     il = &idef->il[F_FBPOSRES];
998     i1 = il->nr/2;
999     idef->iparams_fbposres_nalloc = i1;
1000     srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
1001     for (i = i0; i < i1; i++)
1002     {
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)
1008         {
1009             gmx_incons("Position restraint coordinates are missing");
1010         }
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 */
1016
1017         /* Set the parameter index for idef->iparams_posre */
1018         il->iatoms[i*2] = i;
1019     }
1020 }
1021
1022 static void gen_local_top(const gmx_mtop_t *mtop, const t_inputrec *ir,
1023                           gmx_bool bMergeConstr,
1024                           gmx_localtop_t *top)
1025 {
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;
1030     t_idef                 *idef;
1031     real                   *qA, *qB;
1032     gmx_mtop_atomloop_all_t aloop;
1033     int                     ag;
1034     t_atom                 *atom;
1035
1036     top->atomtypes = mtop->atomtypes;
1037
1038     ffp = &mtop->ffparams;
1039
1040     idef                          = &top->idef;
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;
1052
1053     init_block(&top->cgs);
1054     init_blocka(&top->excls);
1055     for (ftype = 0; ftype < F_NRE; ftype++)
1056     {
1057         idef->il[ftype].nr     = 0;
1058         idef->il[ftype].nalloc = 0;
1059         idef->il[ftype].iatoms = NULL;
1060     }
1061
1062     natoms = 0;
1063     for (mb = 0; mb < mtop->nmolblock; mb++)
1064     {
1065         molb = &mtop->molblock[mb];
1066         molt = &mtop->moltype[molb->type];
1067
1068         srcnr  = molt->atoms.nr;
1069         destnr = natoms;
1070
1071         blockcat(&top->cgs, &molt->cgs, molb->nmol);
1072
1073         blockacat(&top->excls, &molt->excls, molb->nmol, destnr, srcnr);
1074
1075         nposre_old   = idef->il[F_POSRES].nr;
1076         nfbposre_old = idef->il[F_FBPOSRES].nr;
1077         for (ftype = 0; ftype < F_NRE; ftype++)
1078         {
1079             if (bMergeConstr &&
1080                 ftype == F_CONSTR && molt->ilist[F_CONSTRNC].nr > 0)
1081             {
1082                 /* Merge all constrains into one ilist.
1083                  * This simplifies the constraint code.
1084                  */
1085                 for (mol = 0; mol < molb->nmol; mol++)
1086                 {
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);
1091                 }
1092             }
1093             else if (!(bMergeConstr && ftype == F_CONSTRNC))
1094             {
1095                 ilistcat(ftype, &idef->il[ftype], &molt->ilist[ftype],
1096                          molb->nmol, destnr, srcnr);
1097             }
1098         }
1099         if (idef->il[F_POSRES].nr > nposre_old)
1100         {
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);
1104         }
1105         if (idef->il[F_FBPOSRES].nr > nfbposre_old)
1106         {
1107             set_fbposres_params(idef, molb, nfbposre_old/2, natoms);
1108         }
1109
1110         natoms += molb->nmol*srcnr;
1111     }
1112
1113     if (ir == NULL)
1114     {
1115         top->idef.ilsort = ilsortUNKNOWN;
1116     }
1117     else
1118     {
1119         if (ir->efep != efepNO && gmx_mtop_bondeds_free_energy(mtop))
1120         {
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))
1125             {
1126                 qA[ag] = atom->q;
1127                 qB[ag] = atom->qB;
1128             }
1129             gmx_sort_ilist_fe(&top->idef, qA, qB);
1130             sfree(qA);
1131             sfree(qB);
1132         }
1133         else
1134         {
1135             top->idef.ilsort = ilsortNO_FE;
1136         }
1137     }
1138 }
1139
1140 gmx_localtop_t *gmx_mtop_generate_local_top(const gmx_mtop_t *mtop,
1141                                             const t_inputrec *ir)
1142 {
1143     gmx_localtop_t *top;
1144
1145     snew(top, 1);
1146
1147     gen_local_top(mtop, ir, TRUE, top);
1148
1149     return top;
1150 }
1151
1152 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop)
1153 {
1154     int            mt, mb;
1155     gmx_localtop_t ltop;
1156     t_topology     top;
1157
1158     gen_local_top(mtop, NULL, FALSE, &ltop);
1159
1160     top.name      = mtop->name;
1161     top.idef      = ltop.idef;
1162     top.atomtypes = ltop.atomtypes;
1163     top.cgs       = ltop.cgs;
1164     top.excls     = ltop.excls;
1165     top.atoms     = gmx_mtop_global_atoms(mtop);
1166     top.mols      = mtop->mols;
1167     top.symtab    = mtop->symtab;
1168
1169     /* We only need to free the moltype and molblock data,
1170      * all other pointers have been copied to top.
1171      *
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.
1174      */
1175     for (mt = 0; mt < mtop->nmoltype; mt++)
1176     {
1177         done_moltype(&mtop->moltype[mt]);
1178     }
1179     sfree(mtop->moltype);
1180
1181     for (mb = 0; mb < mtop->nmolblock; mb++)
1182     {
1183         done_molblock(&mtop->molblock[mb]);
1184     }
1185     sfree(mtop->molblock);
1186
1187     return top;
1188 }