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