ad33c0a41b55a1ac817694367b3e2243338d9faa
[alexxy/gromacs.git] / src / gromacs / topology / topology.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "topology.h"
40
41 #include <cstdio>
42
43 #include <algorithm>
44
45 #include "gromacs/math/vecdump.h"
46 #include "gromacs/topology/atoms.h"
47 #include "gromacs/topology/idef.h"
48 #include "gromacs/topology/ifunc.h"
49 #include "gromacs/topology/symtab.h"
50 #include "gromacs/utility/compare.h"
51 #include "gromacs/utility/gmxassert.h"
52 #include "gromacs/utility/smalloc.h"
53 #include "gromacs/utility/strconvert.h"
54 #include "gromacs/utility/txtdump.h"
55
56 const char *gtypes[egcNR+1] = {
57     "T-Coupling", "Energy Mon.", "Acceleration", "Freeze",
58     "User1", "User2", "VCM", "Compressed X", "Or. Res. Fit", "QMMM", nullptr
59 };
60
61 static void init_groups(gmx_groups_t *groups)
62 {
63     groups->ngrpname = 0;
64     groups->grpname  = nullptr;
65     for (int g = 0; g < egcNR; g++)
66     {
67         groups->grps[g].nr     = 0;
68         groups->grps[g].nm_ind = nullptr;
69         groups->ngrpnr[g]      = 0;
70         groups->grpnr[g]       = nullptr;
71     }
72
73 }
74
75 void init_mtop(gmx_mtop_t *mtop)
76 {
77     mtop->name = nullptr;
78
79     // TODO: Move to ffparams when that is converted to C++
80     mtop->ffparams.functype.clear();
81     mtop->ffparams.iparams.clear();
82     mtop->ffparams.cmap_grid.ngrid        = 0;
83     mtop->ffparams.cmap_grid.grid_spacing = 0;
84     mtop->ffparams.cmap_grid.cmapdata     = nullptr;
85
86     mtop->moltype.clear();
87     mtop->molblock.clear();
88     mtop->bIntermolecularInteractions = FALSE;
89     mtop->intermolecular_ilist        = nullptr;
90
91     mtop->natoms       = 0;
92     mtop->maxres_renum = 0;
93     mtop->maxresnr     = -1;
94     init_atomtypes(&mtop->atomtypes);
95     init_groups(&mtop->groups);
96     open_symtab(&mtop->symtab);
97 }
98
99 void init_top(t_topology *top)
100 {
101     top->name = nullptr;
102     init_idef(&top->idef);
103     init_atom(&(top->atoms));
104     init_atomtypes(&(top->atomtypes));
105     init_block(&top->cgs);
106     init_block(&top->mols);
107     init_blocka(&top->excls);
108     open_symtab(&top->symtab);
109 }
110
111
112 gmx_moltype_t::gmx_moltype_t() :
113     name(nullptr),
114     cgs(),
115     excls()
116 {
117     init_t_atoms(&atoms, 0, FALSE);
118 }
119
120 gmx_moltype_t::~gmx_moltype_t()
121 {
122     done_atom(&atoms);
123     done_block(&cgs);
124     done_blocka(&excls);
125 }
126
127 void done_gmx_groups_t(gmx_groups_t *g)
128 {
129     int i;
130
131     for (i = 0; (i < egcNR); i++)
132     {
133         if (nullptr != g->grps[i].nm_ind)
134         {
135             sfree(g->grps[i].nm_ind);
136             g->grps[i].nm_ind = nullptr;
137         }
138         if (nullptr != g->grpnr[i])
139         {
140             sfree(g->grpnr[i]);
141             g->grpnr[i] = nullptr;
142         }
143     }
144     /* The contents of this array is in symtab, don't free it here */
145     sfree(g->grpname);
146 }
147
148 gmx_mtop_t::gmx_mtop_t()
149 {
150     init_mtop(this);
151 }
152
153 gmx_mtop_t::~gmx_mtop_t()
154 {
155     done_symtab(&symtab);
156
157     for (int i = 0; i < ffparams.cmap_grid.ngrid; i++)
158     {
159         sfree(ffparams.cmap_grid.cmapdata[i].cmap);
160     }
161     sfree(ffparams.cmap_grid.cmapdata);
162
163     moltype.clear();
164     molblock.clear();
165     done_atomtypes(&atomtypes);
166     done_gmx_groups_t(&groups);
167 }
168
169 void done_top(t_topology *top)
170 {
171     done_idef(&top->idef);
172     done_atom(&(top->atoms));
173
174     /* For GB */
175     done_atomtypes(&(top->atomtypes));
176
177     done_symtab(&(top->symtab));
178     done_block(&(top->cgs));
179     done_block(&(top->mols));
180     done_blocka(&(top->excls));
181 }
182
183 void done_top_mtop(t_topology *top, gmx_mtop_t *mtop)
184 {
185     if (mtop != nullptr)
186     {
187         if (top != nullptr)
188         {
189             done_idef(&top->idef);
190             done_atom(&top->atoms);
191             done_block(&top->cgs);
192             done_blocka(&top->excls);
193             done_block(&top->mols);
194             done_symtab(&top->symtab);
195             open_symtab(&mtop->symtab);
196             done_atomtypes(&(top->atomtypes));
197         }
198
199         // Note that the rest of mtop will be freed by the destructor
200     }
201 }
202
203 void init_localtop(gmx_localtop_t *top)
204 {
205     init_block(&top->cgs);
206     init_blocka(&top->excls);
207     init_idef(&top->idef);
208     init_atomtypes(&top->atomtypes);
209 }
210
211 void done_localtop(gmx_localtop_t *top)
212 {
213     if (top == nullptr)
214     {
215         return;
216     }
217     done_idef(&top->idef);
218     done_block(&top->cgs);
219     done_blocka(&top->excls);
220     done_atomtypes(&top->atomtypes);
221 }
222
223 void done_and_sfree_localtop(gmx_localtop_t *top)
224 {
225     done_localtop(top);
226     sfree(top);
227 }
228
229 bool gmx_mtop_has_masses(const gmx_mtop_t *mtop)
230 {
231     if (mtop == nullptr)
232     {
233         return false;
234     }
235     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveMass;
236 }
237
238 bool gmx_mtop_has_charges(const gmx_mtop_t *mtop)
239 {
240     if (mtop == nullptr)
241     {
242         return false;
243     }
244     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveCharge;
245 }
246
247 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t &mtop)
248 {
249     for (const gmx_moltype_t &moltype : mtop.moltype)
250     {
251         const t_atoms &atoms = moltype.atoms;
252         if (atoms.haveBState)
253         {
254             for (int a = 0; a < atoms.nr; a++)
255             {
256                 if (atoms.atom[a].q != atoms.atom[a].qB)
257                 {
258                     return true;
259                 }
260             }
261         }
262     }
263     return false;
264 }
265
266 bool gmx_mtop_has_atomtypes(const gmx_mtop_t *mtop)
267 {
268     if (mtop == nullptr)
269     {
270         return false;
271     }
272     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveType;
273 }
274
275 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t *mtop)
276 {
277     if (mtop == nullptr)
278     {
279         return false;
280     }
281     return mtop->moltype.empty() || mtop->moltype[0].atoms.havePdbInfo;
282 }
283
284 static void pr_grps(FILE *fp, const char *title, const t_grps grps[], char **grpname[])
285 {
286     int i, j;
287
288     for (i = 0; (i < egcNR); i++)
289     {
290         fprintf(fp, "%s[%-12s] nr=%d, name=[", title, gtypes[i], grps[i].nr);
291         for (j = 0; (j < grps[i].nr); j++)
292         {
293             fprintf(fp, " %s", *(grpname[grps[i].nm_ind[j]]));
294         }
295         fprintf(fp, "]\n");
296     }
297 }
298
299 static void pr_groups(FILE *fp, int indent,
300                       const gmx_groups_t *groups,
301                       gmx_bool bShowNumbers)
302 {
303     int nat_max, i, g;
304
305     pr_grps(fp, "grp", groups->grps, groups->grpname);
306     pr_strings(fp, indent, "grpname", groups->grpname, groups->ngrpname, bShowNumbers);
307
308     pr_indent(fp, indent);
309     fprintf(fp, "groups          ");
310     for (g = 0; g < egcNR; g++)
311     {
312         printf(" %5.5s", gtypes[g]);
313     }
314     printf("\n");
315
316     pr_indent(fp, indent);
317     fprintf(fp, "allocated       ");
318     nat_max = 0;
319     for (g = 0; g < egcNR; g++)
320     {
321         printf(" %5d", groups->ngrpnr[g]);
322         nat_max = std::max(nat_max, groups->ngrpnr[g]);
323     }
324     printf("\n");
325
326     if (nat_max == 0)
327     {
328         pr_indent(fp, indent);
329         fprintf(fp, "groupnr[%5s] =", "*");
330         for (g = 0; g < egcNR; g++)
331         {
332             fprintf(fp, "  %3d ", 0);
333         }
334         fprintf(fp, "\n");
335     }
336     else
337     {
338         for (i = 0; i < nat_max; i++)
339         {
340             pr_indent(fp, indent);
341             fprintf(fp, "groupnr[%5d] =", i);
342             for (g = 0; g < egcNR; g++)
343             {
344                 fprintf(fp, "  %3d ",
345                         groups->grpnr[g] ? groups->grpnr[g][i] : 0);
346             }
347             fprintf(fp, "\n");
348         }
349     }
350 }
351
352 static void pr_moltype(FILE *fp, int indent, const char *title,
353                        const gmx_moltype_t *molt, int n,
354                        const gmx_ffparams_t *ffparams,
355                        gmx_bool bShowNumbers, gmx_bool bShowParameters)
356 {
357     int j;
358
359     indent = pr_title_n(fp, indent, title, n);
360     pr_indent(fp, indent);
361     fprintf(fp, "name=\"%s\"\n", *(molt->name));
362     pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
363     pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
364     pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
365     for (j = 0; (j < F_NRE); j++)
366     {
367         pr_ilist(fp, indent, interaction_function[j].longname,
368                  ffparams->functype.data(), molt->ilist[j],
369                  bShowNumbers, bShowParameters, ffparams->iparams.data());
370     }
371 }
372
373 static void pr_molblock(FILE *fp, int indent, const char *title,
374                         const gmx_molblock_t *molb, int n,
375                         const std::vector<gmx_moltype_t> &molt)
376 {
377     indent = pr_title_n(fp, indent, title, n);
378     pr_indent(fp, indent);
379     fprintf(fp, "%-20s = %d \"%s\"\n",
380             "moltype", molb->type, *(molt[molb->type].name));
381     pr_int(fp, indent, "#molecules", molb->nmol);
382     pr_int(fp, indent, "#posres_xA", molb->posres_xA.size());
383     if (!molb->posres_xA.empty())
384     {
385         pr_rvecs(fp, indent, "posres_xA", as_rvec_array(molb->posres_xA.data()), molb->posres_xA.size());
386     }
387     pr_int(fp, indent, "#posres_xB", molb->posres_xB.size());
388     if (!molb->posres_xB.empty())
389     {
390         pr_rvecs(fp, indent, "posres_xB", as_rvec_array(molb->posres_xB.data()), molb->posres_xB.size());
391     }
392 }
393
394 void pr_mtop(FILE *fp, int indent, const char *title, const gmx_mtop_t *mtop,
395              gmx_bool bShowNumbers, gmx_bool bShowParameters)
396 {
397     if (available(fp, mtop, indent, title))
398     {
399         indent = pr_title(fp, indent, title);
400         pr_indent(fp, indent);
401         fprintf(fp, "name=\"%s\"\n", *(mtop->name));
402         pr_int(fp, indent, "#atoms", mtop->natoms);
403         pr_int(fp, indent, "#molblock", mtop->molblock.size());
404         for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
405         {
406             pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
407         }
408         pr_str(fp, indent, "bIntermolecularInteractions",
409                gmx::boolToString(mtop->bIntermolecularInteractions));
410         if (mtop->bIntermolecularInteractions)
411         {
412             for (int j = 0; j < F_NRE; j++)
413             {
414                 pr_ilist(fp, indent, interaction_function[j].longname,
415                          mtop->ffparams.functype.data(),
416                          (*mtop->intermolecular_ilist)[j],
417                          bShowNumbers, bShowParameters, mtop->ffparams.iparams.data());
418             }
419         }
420         pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
421         pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
422         for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
423         {
424             pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
425                        &mtop->ffparams, bShowNumbers, bShowParameters);
426         }
427         pr_groups(fp, indent, &mtop->groups, bShowNumbers);
428     }
429 }
430
431 void pr_top(FILE *fp, int indent, const char *title, const t_topology *top,
432             gmx_bool bShowNumbers, gmx_bool bShowParameters)
433 {
434     if (available(fp, top, indent, title))
435     {
436         indent = pr_title(fp, indent, title);
437         pr_indent(fp, indent);
438         fprintf(fp, "name=\"%s\"\n", *(top->name));
439         pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
440         pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
441         pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
442         pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
443         pr_str(fp, indent, "bIntermolecularInteractions",
444                gmx::boolToString(top->bIntermolecularInteractions));
445         pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
446         pr_idef(fp, indent, "idef", &top->idef, bShowNumbers, bShowParameters);
447     }
448 }
449
450 static void cmp_ilist(FILE *fp, int ftype, const t_ilist *il1, const t_ilist *il2)
451 {
452     int  i;
453     char buf[256];
454
455     fprintf(fp, "comparing ilist %s\n", interaction_function[ftype].name);
456     sprintf(buf, "%s->nr", interaction_function[ftype].name);
457     cmp_int(fp, buf, -1, il1->nr, il2->nr);
458     sprintf(buf, "%s->iatoms", interaction_function[ftype].name);
459     if (((il1->nr > 0) && (!il1->iatoms)) ||
460         ((il2->nr > 0) && (!il2->iatoms)) ||
461         ((il1->nr != il2->nr)))
462     {
463         fprintf(fp, "Comparing radically different topologies - %s is different\n",
464                 buf);
465     }
466     else
467     {
468         for (i = 0; (i < il1->nr); i++)
469         {
470             cmp_int(fp, buf, i, il1->iatoms[i], il2->iatoms[i]);
471         }
472     }
473 }
474
475 static void cmp_iparm(FILE *fp, const char *s, t_functype ft,
476                       const t_iparams &ip1, const t_iparams &ip2, real ftol, real abstol)
477 {
478     int      i;
479     gmx_bool bDiff;
480
481     bDiff = FALSE;
482     for (i = 0; i < MAXFORCEPARAM && !bDiff; i++)
483     {
484         bDiff = !equal_real(ip1.generic.buf[i], ip2.generic.buf[i], ftol, abstol);
485     }
486     if (bDiff)
487     {
488         fprintf(fp, "%s1: ", s);
489         pr_iparams(fp, ft, &ip1);
490         fprintf(fp, "%s2: ", s);
491         pr_iparams(fp, ft, &ip2);
492     }
493 }
494
495 static void cmp_iparm_AB(FILE *fp, const char *s, t_functype ft,
496                          const t_iparams &ip1, real ftol, real abstol)
497 {
498     int      nrfpA, nrfpB, p0, i;
499     gmx_bool bDiff;
500
501     /* Normally the first parameter is perturbable */
502     p0    = 0;
503     nrfpA = interaction_function[ft].nrfpA;
504     nrfpB = interaction_function[ft].nrfpB;
505     if (ft == F_PDIHS)
506     {
507         nrfpB = 2;
508     }
509     else if (interaction_function[ft].flags & IF_TABULATED)
510     {
511         /* For tabulated interactions only the second parameter is perturbable */
512         p0    = 1;
513         nrfpB = 1;
514     }
515     bDiff = FALSE;
516     for (i = 0; i < nrfpB && !bDiff; i++)
517     {
518         bDiff = !equal_real(ip1.generic.buf[p0+i], ip1.generic.buf[nrfpA+i], ftol, abstol);
519     }
520     if (bDiff)
521     {
522         fprintf(fp, "%s: ", s);
523         pr_iparams(fp, ft, &ip1);
524     }
525 }
526
527 static void cmp_cmap(FILE *fp, const gmx_cmap_t *cmap1, const gmx_cmap_t *cmap2, real ftol, real abstol)
528 {
529     cmp_int(fp, "cmap ngrid", -1, cmap1->ngrid, cmap2->ngrid);
530     cmp_int(fp, "cmap grid_spacing", -1, cmap1->grid_spacing, cmap2->grid_spacing);
531     if (cmap1->ngrid == cmap2->ngrid &&
532         cmap1->grid_spacing == cmap2->grid_spacing)
533     {
534         int g;
535
536         for (g = 0; g < cmap1->ngrid; g++)
537         {
538             int i;
539
540             fprintf(fp, "comparing cmap %d\n", g);
541
542             for (i = 0; i < 4*cmap1->grid_spacing*cmap1->grid_spacing; i++)
543             {
544                 cmp_real(fp, "", i, cmap1->cmapdata[g].cmap[i], cmap2->cmapdata[g].cmap[i], ftol, abstol);
545             }
546         }
547     }
548 }
549
550 static void cmp_idef(FILE *fp, const t_idef *id1, const t_idef *id2, real ftol, real abstol)
551 {
552     int  i;
553     char buf1[64], buf2[64];
554
555     fprintf(fp, "comparing idef\n");
556     if (id2)
557     {
558         cmp_int(fp, "idef->ntypes", -1, id1->ntypes, id2->ntypes);
559         cmp_int(fp, "idef->atnr",  -1, id1->atnr, id2->atnr);
560         for (i = 0; (i < std::min(id1->ntypes, id2->ntypes)); i++)
561         {
562             sprintf(buf1, "idef->functype[%d]", i);
563             sprintf(buf2, "idef->iparam[%d]", i);
564             cmp_int(fp, buf1, i, static_cast<int>(id1->functype[i]), static_cast<int>(id2->functype[i]));
565             cmp_iparm(fp, buf2, id1->functype[i],
566                       id1->iparams[i], id2->iparams[i], ftol, abstol);
567         }
568         cmp_real(fp, "fudgeQQ", -1, id1->fudgeQQ, id2->fudgeQQ, ftol, abstol);
569         cmp_cmap(fp, &id1->cmap_grid, &id2->cmap_grid, ftol, abstol);
570         for (i = 0; (i < F_NRE); i++)
571         {
572             cmp_ilist(fp, i, &(id1->il[i]), &(id2->il[i]));
573         }
574     }
575     else
576     {
577         for (i = 0; (i < id1->ntypes); i++)
578         {
579             cmp_iparm_AB(fp, "idef->iparam", id1->functype[i], id1->iparams[i], ftol, abstol);
580         }
581     }
582 }
583
584 static void cmp_block(FILE *fp, const t_block *b1, const t_block *b2, const char *s)
585 {
586     char buf[32];
587
588     fprintf(fp, "comparing block %s\n", s);
589     sprintf(buf, "%s.nr", s);
590     cmp_int(fp, buf, -1, b1->nr, b2->nr);
591 }
592
593 static void cmp_blocka(FILE *fp, const t_blocka *b1, const t_blocka *b2, const char *s)
594 {
595     char buf[32];
596
597     fprintf(fp, "comparing blocka %s\n", s);
598     sprintf(buf, "%s.nr", s);
599     cmp_int(fp, buf, -1, b1->nr, b2->nr);
600     sprintf(buf, "%s.nra", s);
601     cmp_int(fp, buf, -1, b1->nra, b2->nra);
602 }
603
604 void cmp_top(FILE *fp, const t_topology *t1, const t_topology *t2, real ftol, real abstol)
605 {
606     fprintf(fp, "comparing top\n");
607     if (t2)
608     {
609         cmp_idef(fp, &(t1->idef), &(t2->idef), ftol, abstol);
610         cmp_atoms(fp, &(t1->atoms), &(t2->atoms), ftol, abstol);
611         cmp_block(fp, &t1->cgs, &t2->cgs, "cgs");
612         cmp_block(fp, &t1->mols, &t2->mols, "mols");
613         cmp_bool(fp, "bIntermolecularInteractions", -1, t1->bIntermolecularInteractions, t2->bIntermolecularInteractions);
614         cmp_blocka(fp, &t1->excls, &t2->excls, "excls");
615     }
616     else
617     {
618         cmp_idef(fp, &(t1->idef), nullptr, ftol, abstol);
619         cmp_atoms(fp, &(t1->atoms), nullptr, ftol, abstol);
620     }
621 }
622
623 void cmp_groups(FILE *fp, const gmx_groups_t *g0, const gmx_groups_t *g1,
624                 int natoms0, int natoms1)
625 {
626     char buf[32];
627
628     fprintf(fp, "comparing groups\n");
629
630     for (int i = 0; i < egcNR; i++)
631     {
632         sprintf(buf, "grps[%d].nr", i);
633         cmp_int(fp, buf, -1, g0->grps[i].nr, g1->grps[i].nr);
634         if (g0->grps[i].nr == g1->grps[i].nr)
635         {
636             for (int j = 0; j < g0->grps[i].nr; j++)
637             {
638                 sprintf(buf, "grps[%d].name[%d]", i, j);
639                 cmp_str(fp, buf, -1,
640                         *g0->grpname[g0->grps[i].nm_ind[j]],
641                         *g1->grpname[g1->grps[i].nm_ind[j]]);
642             }
643         }
644         cmp_int(fp, "ngrpnr", i, g0->ngrpnr[i], g1->ngrpnr[i]);
645         if (g0->ngrpnr[i] == g1->ngrpnr[i] && natoms0 == natoms1 &&
646             (g0->grpnr[i] != nullptr || g1->grpnr[i] != nullptr))
647         {
648             for (int j = 0; j < natoms0; j++)
649             {
650                 cmp_int(fp, gtypes[i], j, getGroupType(g0, i, j), getGroupType(g1, i, j));
651             }
652         }
653     }
654     /* We have compared the names in the groups lists,
655      * so we can skip the grpname list comparison.
656      */
657 }
658
659 int getGroupType(const gmx_groups_t *group, int type, int atom)
660 {
661     return (group->grpnr[type] ? group->grpnr[type][atom] : 0);
662 }