568dbee391cad1f6c2e703dfcca4fd259a1a662b
[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               = nullptr;
81     mtop->ffparams.iparams                = nullptr;
82     mtop->ffparams.cmap_grid.ngrid        = 0;
83     mtop->ffparams.cmap_grid.grid_spacing = 0;
84     mtop->ffparams.cmap_grid.cmapdata     = nullptr;
85     mtop->ffparams.ntypes                 = 0;
86
87     mtop->moltype.clear();
88     mtop->molblock.clear();
89     mtop->bIntermolecularInteractions = FALSE;
90     mtop->intermolecular_ilist        = nullptr;
91
92     mtop->natoms       = 0;
93     mtop->maxres_renum = 0;
94     mtop->maxresnr     = -1;
95     init_atomtypes(&mtop->atomtypes);
96     init_groups(&mtop->groups);
97     open_symtab(&mtop->symtab);
98 }
99
100 void init_top(t_topology *top)
101 {
102     top->name = nullptr;
103     init_idef(&top->idef);
104     init_atom(&(top->atoms));
105     init_atomtypes(&(top->atomtypes));
106     init_block(&top->cgs);
107     init_block(&top->mols);
108     init_blocka(&top->excls);
109     open_symtab(&top->symtab);
110 }
111
112
113 gmx_moltype_t::gmx_moltype_t() :
114     name(nullptr),
115     cgs(),
116     excls()
117 {
118     init_t_atoms(&atoms, 0, FALSE);
119 }
120
121 gmx_moltype_t::~gmx_moltype_t()
122 {
123     done_atom(&atoms);
124     done_block(&cgs);
125     done_blocka(&excls);
126 }
127
128 void done_gmx_groups_t(gmx_groups_t *g)
129 {
130     int i;
131
132     for (i = 0; (i < egcNR); i++)
133     {
134         if (nullptr != g->grps[i].nm_ind)
135         {
136             sfree(g->grps[i].nm_ind);
137             g->grps[i].nm_ind = nullptr;
138         }
139         if (nullptr != g->grpnr[i])
140         {
141             sfree(g->grpnr[i]);
142             g->grpnr[i] = nullptr;
143         }
144     }
145     /* The contents of this array is in symtab, don't free it here */
146     sfree(g->grpname);
147 }
148
149 gmx_mtop_t::gmx_mtop_t()
150 {
151     init_mtop(this);
152 }
153
154 gmx_mtop_t::~gmx_mtop_t()
155 {
156     done_symtab(&symtab);
157
158     sfree(ffparams.functype);
159     sfree(ffparams.iparams);
160     for (int i = 0; i < ffparams.cmap_grid.ngrid; i++)
161     {
162         sfree(ffparams.cmap_grid.cmapdata[i].cmap);
163     }
164     sfree(ffparams.cmap_grid.cmapdata);
165
166     moltype.clear();
167     molblock.clear();
168     done_atomtypes(&atomtypes);
169     done_gmx_groups_t(&groups);
170 }
171
172 void done_top(t_topology *top)
173 {
174     done_idef(&top->idef);
175     done_atom(&(top->atoms));
176
177     /* For GB */
178     done_atomtypes(&(top->atomtypes));
179
180     done_symtab(&(top->symtab));
181     done_block(&(top->cgs));
182     done_block(&(top->mols));
183     done_blocka(&(top->excls));
184 }
185
186 void done_top_mtop(t_topology *top, gmx_mtop_t *mtop)
187 {
188     if (mtop != nullptr)
189     {
190         if (top != nullptr)
191         {
192             done_idef(&top->idef);
193             done_atom(&top->atoms);
194             done_block(&top->cgs);
195             done_blocka(&top->excls);
196             done_block(&top->mols);
197             done_symtab(&top->symtab);
198             open_symtab(&mtop->symtab);
199             done_atomtypes(&(top->atomtypes));
200         }
201
202         // Note that the rest of mtop will be freed by the destructor
203     }
204 }
205
206 void init_localtop(gmx_localtop_t *top)
207 {
208     init_block(&top->cgs);
209     init_blocka(&top->excls);
210     init_idef(&top->idef);
211     init_atomtypes(&top->atomtypes);
212 }
213
214 void done_localtop(gmx_localtop_t *top)
215 {
216     if (top == nullptr)
217     {
218         return;
219     }
220     done_idef(&top->idef);
221     done_block(&top->cgs);
222     done_blocka(&top->excls);
223     done_atomtypes(&top->atomtypes);
224 }
225
226 void done_and_sfree_localtop(gmx_localtop_t *top)
227 {
228     done_localtop(top);
229     sfree(top);
230 }
231
232 bool gmx_mtop_has_masses(const gmx_mtop_t *mtop)
233 {
234     if (mtop == nullptr)
235     {
236         return false;
237     }
238     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveMass;
239 }
240
241 bool gmx_mtop_has_charges(const gmx_mtop_t *mtop)
242 {
243     if (mtop == nullptr)
244     {
245         return false;
246     }
247     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveCharge;
248 }
249
250 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t &mtop)
251 {
252     for (const gmx_moltype_t &moltype : mtop.moltype)
253     {
254         const t_atoms &atoms = moltype.atoms;
255         if (atoms.haveBState)
256         {
257             for (int a = 0; a < atoms.nr; a++)
258             {
259                 if (atoms.atom[a].q != atoms.atom[a].qB)
260                 {
261                     return true;
262                 }
263             }
264         }
265     }
266     return false;
267 }
268
269 bool gmx_mtop_has_atomtypes(const gmx_mtop_t *mtop)
270 {
271     if (mtop == nullptr)
272     {
273         return false;
274     }
275     return mtop->moltype.empty() || mtop->moltype[0].atoms.haveType;
276 }
277
278 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t *mtop)
279 {
280     if (mtop == nullptr)
281     {
282         return false;
283     }
284     return mtop->moltype.empty() || mtop->moltype[0].atoms.havePdbInfo;
285 }
286
287 static void pr_grps(FILE *fp, const char *title, const t_grps grps[], char **grpname[])
288 {
289     int i, j;
290
291     for (i = 0; (i < egcNR); i++)
292     {
293         fprintf(fp, "%s[%-12s] nr=%d, name=[", title, gtypes[i], grps[i].nr);
294         for (j = 0; (j < grps[i].nr); j++)
295         {
296             fprintf(fp, " %s", *(grpname[grps[i].nm_ind[j]]));
297         }
298         fprintf(fp, "]\n");
299     }
300 }
301
302 static void pr_groups(FILE *fp, int indent,
303                       const gmx_groups_t *groups,
304                       gmx_bool bShowNumbers)
305 {
306     int nat_max, i, g;
307
308     pr_grps(fp, "grp", groups->grps, groups->grpname);
309     pr_strings(fp, indent, "grpname", groups->grpname, groups->ngrpname, bShowNumbers);
310
311     pr_indent(fp, indent);
312     fprintf(fp, "groups          ");
313     for (g = 0; g < egcNR; g++)
314     {
315         printf(" %5.5s", gtypes[g]);
316     }
317     printf("\n");
318
319     pr_indent(fp, indent);
320     fprintf(fp, "allocated       ");
321     nat_max = 0;
322     for (g = 0; g < egcNR; g++)
323     {
324         printf(" %5d", groups->ngrpnr[g]);
325         nat_max = std::max(nat_max, groups->ngrpnr[g]);
326     }
327     printf("\n");
328
329     if (nat_max == 0)
330     {
331         pr_indent(fp, indent);
332         fprintf(fp, "groupnr[%5s] =", "*");
333         for (g = 0; g < egcNR; g++)
334         {
335             fprintf(fp, "  %3d ", 0);
336         }
337         fprintf(fp, "\n");
338     }
339     else
340     {
341         for (i = 0; i < nat_max; i++)
342         {
343             pr_indent(fp, indent);
344             fprintf(fp, "groupnr[%5d] =", i);
345             for (g = 0; g < egcNR; g++)
346             {
347                 fprintf(fp, "  %3d ",
348                         groups->grpnr[g] ? groups->grpnr[g][i] : 0);
349             }
350             fprintf(fp, "\n");
351         }
352     }
353 }
354
355 static void pr_moltype(FILE *fp, int indent, const char *title,
356                        const gmx_moltype_t *molt, int n,
357                        const gmx_ffparams_t *ffparams,
358                        gmx_bool bShowNumbers, gmx_bool bShowParameters)
359 {
360     int j;
361
362     indent = pr_title_n(fp, indent, title, n);
363     pr_indent(fp, indent);
364     fprintf(fp, "name=\"%s\"\n", *(molt->name));
365     pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
366     pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
367     pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
368     for (j = 0; (j < F_NRE); j++)
369     {
370         pr_ilist(fp, indent, interaction_function[j].longname,
371                  ffparams->functype, molt->ilist[j],
372                  bShowNumbers, bShowParameters, ffparams->iparams);
373     }
374 }
375
376 static void pr_molblock(FILE *fp, int indent, const char *title,
377                         const gmx_molblock_t *molb, int n,
378                         const std::vector<gmx_moltype_t> &molt)
379 {
380     indent = pr_title_n(fp, indent, title, n);
381     pr_indent(fp, indent);
382     fprintf(fp, "%-20s = %d \"%s\"\n",
383             "moltype", molb->type, *(molt[molb->type].name));
384     pr_int(fp, indent, "#molecules", molb->nmol);
385     pr_int(fp, indent, "#posres_xA", molb->posres_xA.size());
386     if (!molb->posres_xA.empty())
387     {
388         pr_rvecs(fp, indent, "posres_xA", as_rvec_array(molb->posres_xA.data()), molb->posres_xA.size());
389     }
390     pr_int(fp, indent, "#posres_xB", molb->posres_xB.size());
391     if (!molb->posres_xB.empty())
392     {
393         pr_rvecs(fp, indent, "posres_xB", as_rvec_array(molb->posres_xB.data()), molb->posres_xB.size());
394     }
395 }
396
397 void pr_mtop(FILE *fp, int indent, const char *title, const gmx_mtop_t *mtop,
398              gmx_bool bShowNumbers, gmx_bool bShowParameters)
399 {
400     if (available(fp, mtop, indent, title))
401     {
402         indent = pr_title(fp, indent, title);
403         pr_indent(fp, indent);
404         fprintf(fp, "name=\"%s\"\n", *(mtop->name));
405         pr_int(fp, indent, "#atoms", mtop->natoms);
406         pr_int(fp, indent, "#molblock", mtop->molblock.size());
407         for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
408         {
409             pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
410         }
411         pr_str(fp, indent, "bIntermolecularInteractions",
412                gmx::boolToString(mtop->bIntermolecularInteractions));
413         if (mtop->bIntermolecularInteractions)
414         {
415             for (int j = 0; j < F_NRE; j++)
416             {
417                 pr_ilist(fp, indent, interaction_function[j].longname,
418                          mtop->ffparams.functype,
419                          (*mtop->intermolecular_ilist)[j],
420                          bShowNumbers, bShowParameters, mtop->ffparams.iparams);
421             }
422         }
423         pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
424         pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
425         for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
426         {
427             pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
428                        &mtop->ffparams, bShowNumbers, bShowParameters);
429         }
430         pr_groups(fp, indent, &mtop->groups, bShowNumbers);
431     }
432 }
433
434 void pr_top(FILE *fp, int indent, const char *title, const t_topology *top,
435             gmx_bool bShowNumbers, gmx_bool bShowParameters)
436 {
437     if (available(fp, top, indent, title))
438     {
439         indent = pr_title(fp, indent, title);
440         pr_indent(fp, indent);
441         fprintf(fp, "name=\"%s\"\n", *(top->name));
442         pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
443         pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
444         pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
445         pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
446         pr_str(fp, indent, "bIntermolecularInteractions",
447                gmx::boolToString(top->bIntermolecularInteractions));
448         pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
449         pr_idef(fp, indent, "idef", &top->idef, bShowNumbers, bShowParameters);
450     }
451 }
452
453 static void cmp_ilist(FILE *fp, int ftype, const t_ilist *il1, const t_ilist *il2)
454 {
455     int  i;
456     char buf[256];
457
458     fprintf(fp, "comparing ilist %s\n", interaction_function[ftype].name);
459     sprintf(buf, "%s->nr", interaction_function[ftype].name);
460     cmp_int(fp, buf, -1, il1->nr, il2->nr);
461     sprintf(buf, "%s->iatoms", interaction_function[ftype].name);
462     if (((il1->nr > 0) && (!il1->iatoms)) ||
463         ((il2->nr > 0) && (!il2->iatoms)) ||
464         ((il1->nr != il2->nr)))
465     {
466         fprintf(fp, "Comparing radically different topologies - %s is different\n",
467                 buf);
468     }
469     else
470     {
471         for (i = 0; (i < il1->nr); i++)
472         {
473             cmp_int(fp, buf, i, il1->iatoms[i], il2->iatoms[i]);
474         }
475     }
476 }
477
478 static void cmp_iparm(FILE *fp, const char *s, t_functype ft,
479                       const t_iparams &ip1, const t_iparams &ip2, real ftol, real abstol)
480 {
481     int      i;
482     gmx_bool bDiff;
483
484     bDiff = FALSE;
485     for (i = 0; i < MAXFORCEPARAM && !bDiff; i++)
486     {
487         bDiff = !equal_real(ip1.generic.buf[i], ip2.generic.buf[i], ftol, abstol);
488     }
489     if (bDiff)
490     {
491         fprintf(fp, "%s1: ", s);
492         pr_iparams(fp, ft, &ip1);
493         fprintf(fp, "%s2: ", s);
494         pr_iparams(fp, ft, &ip2);
495     }
496 }
497
498 static void cmp_iparm_AB(FILE *fp, const char *s, t_functype ft,
499                          const t_iparams &ip1, real ftol, real abstol)
500 {
501     int      nrfpA, nrfpB, p0, i;
502     gmx_bool bDiff;
503
504     /* Normally the first parameter is perturbable */
505     p0    = 0;
506     nrfpA = interaction_function[ft].nrfpA;
507     nrfpB = interaction_function[ft].nrfpB;
508     if (ft == F_PDIHS)
509     {
510         nrfpB = 2;
511     }
512     else if (interaction_function[ft].flags & IF_TABULATED)
513     {
514         /* For tabulated interactions only the second parameter is perturbable */
515         p0    = 1;
516         nrfpB = 1;
517     }
518     bDiff = FALSE;
519     for (i = 0; i < nrfpB && !bDiff; i++)
520     {
521         bDiff = !equal_real(ip1.generic.buf[p0+i], ip1.generic.buf[nrfpA+i], ftol, abstol);
522     }
523     if (bDiff)
524     {
525         fprintf(fp, "%s: ", s);
526         pr_iparams(fp, ft, &ip1);
527     }
528 }
529
530 static void cmp_cmap(FILE *fp, const gmx_cmap_t *cmap1, const gmx_cmap_t *cmap2, real ftol, real abstol)
531 {
532     cmp_int(fp, "cmap ngrid", -1, cmap1->ngrid, cmap2->ngrid);
533     cmp_int(fp, "cmap grid_spacing", -1, cmap1->grid_spacing, cmap2->grid_spacing);
534     if (cmap1->ngrid == cmap2->ngrid &&
535         cmap1->grid_spacing == cmap2->grid_spacing)
536     {
537         int g;
538
539         for (g = 0; g < cmap1->ngrid; g++)
540         {
541             int i;
542
543             fprintf(fp, "comparing cmap %d\n", g);
544
545             for (i = 0; i < 4*cmap1->grid_spacing*cmap1->grid_spacing; i++)
546             {
547                 cmp_real(fp, "", i, cmap1->cmapdata[g].cmap[i], cmap2->cmapdata[g].cmap[i], ftol, abstol);
548             }
549         }
550     }
551 }
552
553 static void cmp_idef(FILE *fp, const t_idef *id1, const t_idef *id2, real ftol, real abstol)
554 {
555     int  i;
556     char buf1[64], buf2[64];
557
558     fprintf(fp, "comparing idef\n");
559     if (id2)
560     {
561         cmp_int(fp, "idef->ntypes", -1, id1->ntypes, id2->ntypes);
562         cmp_int(fp, "idef->atnr",  -1, id1->atnr, id2->atnr);
563         for (i = 0; (i < std::min(id1->ntypes, id2->ntypes)); i++)
564         {
565             sprintf(buf1, "idef->functype[%d]", i);
566             sprintf(buf2, "idef->iparam[%d]", i);
567             cmp_int(fp, buf1, i, static_cast<int>(id1->functype[i]), static_cast<int>(id2->functype[i]));
568             cmp_iparm(fp, buf2, id1->functype[i],
569                       id1->iparams[i], id2->iparams[i], ftol, abstol);
570         }
571         cmp_real(fp, "fudgeQQ", -1, id1->fudgeQQ, id2->fudgeQQ, ftol, abstol);
572         cmp_cmap(fp, &id1->cmap_grid, &id2->cmap_grid, ftol, abstol);
573         for (i = 0; (i < F_NRE); i++)
574         {
575             cmp_ilist(fp, i, &(id1->il[i]), &(id2->il[i]));
576         }
577     }
578     else
579     {
580         for (i = 0; (i < id1->ntypes); i++)
581         {
582             cmp_iparm_AB(fp, "idef->iparam", id1->functype[i], id1->iparams[i], ftol, abstol);
583         }
584     }
585 }
586
587 static void cmp_block(FILE *fp, const t_block *b1, const t_block *b2, const char *s)
588 {
589     char buf[32];
590
591     fprintf(fp, "comparing block %s\n", s);
592     sprintf(buf, "%s.nr", s);
593     cmp_int(fp, buf, -1, b1->nr, b2->nr);
594 }
595
596 static void cmp_blocka(FILE *fp, const t_blocka *b1, const t_blocka *b2, const char *s)
597 {
598     char buf[32];
599
600     fprintf(fp, "comparing blocka %s\n", s);
601     sprintf(buf, "%s.nr", s);
602     cmp_int(fp, buf, -1, b1->nr, b2->nr);
603     sprintf(buf, "%s.nra", s);
604     cmp_int(fp, buf, -1, b1->nra, b2->nra);
605 }
606
607 void cmp_top(FILE *fp, const t_topology *t1, const t_topology *t2, real ftol, real abstol)
608 {
609     fprintf(fp, "comparing top\n");
610     if (t2)
611     {
612         cmp_idef(fp, &(t1->idef), &(t2->idef), ftol, abstol);
613         cmp_atoms(fp, &(t1->atoms), &(t2->atoms), ftol, abstol);
614         cmp_block(fp, &t1->cgs, &t2->cgs, "cgs");
615         cmp_block(fp, &t1->mols, &t2->mols, "mols");
616         cmp_bool(fp, "bIntermolecularInteractions", -1, t1->bIntermolecularInteractions, t2->bIntermolecularInteractions);
617         cmp_blocka(fp, &t1->excls, &t2->excls, "excls");
618     }
619     else
620     {
621         cmp_idef(fp, &(t1->idef), nullptr, ftol, abstol);
622         cmp_atoms(fp, &(t1->atoms), nullptr, ftol, abstol);
623     }
624 }
625
626 void cmp_groups(FILE *fp, const gmx_groups_t *g0, const gmx_groups_t *g1,
627                 int natoms0, int natoms1)
628 {
629     char buf[32];
630
631     fprintf(fp, "comparing groups\n");
632
633     for (int i = 0; i < egcNR; i++)
634     {
635         sprintf(buf, "grps[%d].nr", i);
636         cmp_int(fp, buf, -1, g0->grps[i].nr, g1->grps[i].nr);
637         if (g0->grps[i].nr == g1->grps[i].nr)
638         {
639             for (int j = 0; j < g0->grps[i].nr; j++)
640             {
641                 sprintf(buf, "grps[%d].name[%d]", i, j);
642                 cmp_str(fp, buf, -1,
643                         *g0->grpname[g0->grps[i].nm_ind[j]],
644                         *g1->grpname[g1->grps[i].nm_ind[j]]);
645             }
646         }
647         cmp_int(fp, "ngrpnr", i, g0->ngrpnr[i], g1->ngrpnr[i]);
648         if (g0->ngrpnr[i] == g1->ngrpnr[i] && natoms0 == natoms1 &&
649             (g0->grpnr[i] != nullptr || g1->grpnr[i] != nullptr))
650         {
651             for (int j = 0; j < natoms0; j++)
652             {
653                 cmp_int(fp, gtypes[i], j, getGroupType(g0, i, j), getGroupType(g1, i, j));
654             }
655         }
656     }
657     /* We have compared the names in the groups lists,
658      * so we can skip the grpname list comparison.
659      */
660 }
661
662 int getGroupType(const gmx_groups_t *group, int type, int atom)
663 {
664     return (group->grpnr[type] ? group->grpnr[type][atom] : 0);
665 }