Change pointers to InteractionList to references
authorBerk Hess <hess@kth.se>
Sat, 15 Sep 2018 21:02:13 +0000 (23:02 +0200)
committerBerk Hess <hess@kth.se>
Mon, 17 Sep 2018 07:01:40 +0000 (09:01 +0200)
Change-Id: Ia2431137e7326a4fa68083b0e50a82455b426c86

src/gromacs/domdec/domdec_topology.cpp
src/gromacs/gmxpreprocess/grompp.cpp
src/gromacs/mdlib/calc_verletbuf.cpp
src/gromacs/topology/idef.cpp
src/gromacs/topology/idef.h
src/gromacs/topology/topology.cpp

index 58d1a80b2d3f55a474158ae6b347897214f72a06..4cd038040a0ec30644923157f87249f6dee83b45 100644 (file)
@@ -2458,25 +2458,25 @@ static void bonded_cg_distance_mol(const gmx_moltype_t *molt,
     {
         if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
         {
-            const auto    *il   = &molt->ilist[ftype];
+            const auto    &il   = molt->ilist[ftype];
             int            nral = NRAL(ftype);
             if (nral > 1)
             {
-                for (int i = 0; i < il->size(); i += 1+nral)
+                for (int i = 0; i < il.size(); i += 1+nral)
                 {
                     for (int ai = 0; ai < nral; ai++)
                     {
-                        int cgi = at2cg[il->iatoms[i+1+ai]];
+                        int cgi = at2cg[il.iatoms[i+1+ai]];
                         for (int aj = ai + 1; aj < nral; aj++)
                         {
-                            int cgj = at2cg[il->iatoms[i+1+aj]];
+                            int cgj = at2cg[il.iatoms[i+1+aj]];
                             if (cgi != cgj)
                             {
                                 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
 
                                 update_max_bonded_distance(rij2, ftype,
-                                                           il->iatoms[i+1+ai],
-                                                           il->iatoms[i+1+aj],
+                                                           il.iatoms[i+1+ai],
+                                                           il.iatoms[i+1+aj],
                                                            (nral == 2) ? bd_2b : bd_mb);
                             }
                         }
@@ -2521,24 +2521,24 @@ static void bonded_distance_intermol(const InteractionLists &ilists_intermol,
     {
         if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
         {
-            const auto    *il   = &ilists_intermol[ftype];
+            const auto    &il   = ilists_intermol[ftype];
             int            nral = NRAL(ftype);
 
             /* No nral>1 check here, since intermol interactions always
              * have nral>=2 (and the code is also correct for nral=1).
              */
-            for (int i = 0; i < il->size(); i += 1+nral)
+            for (int i = 0; i < il.size(); i += 1+nral)
             {
                 for (int ai = 0; ai < nral; ai++)
                 {
-                    int atom_i = il->iatoms[i + 1 + ai];
+                    int atom_i = il.iatoms[i + 1 + ai];
 
                     for (int aj = ai + 1; aj < nral; aj++)
                     {
                         rvec dx;
                         real rij2;
 
-                        int  atom_j = il->iatoms[i + 1 + aj];
+                        int  atom_j = il.iatoms[i + 1 + aj];
 
                         pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
 
index 2d6c9edae9b18f839c1ce6184665dc9432911c5a..1088cd591b90b60aa57fa0726200ef1d457d11fa 100644 (file)
@@ -266,8 +266,8 @@ static void check_bonds_timestep(const gmx_mtop_t *mtop, double dt, warninp_t wi
     {
         const t_atom           *atom  = moltype.atoms.atom;
         const InteractionLists &ilist = moltype.ilist;
-        const InteractionList  *ilc   = &ilist[F_CONSTR];
-        const InteractionList  *ils   = &ilist[F_SETTLE];
+        const InteractionList  &ilc   = ilist[F_CONSTR];
+        const InteractionList  &ils   = ilist[F_SETTLE];
         for (ftype = 0; ftype < F_NRE; ftype++)
         {
             if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
@@ -275,18 +275,18 @@ static void check_bonds_timestep(const gmx_mtop_t *mtop, double dt, warninp_t wi
                 continue;
             }
 
-            const InteractionList *ilb = &ilist[ftype];
-            for (i = 0; i < ilb->size(); i += 3)
+            const InteractionList &ilb = ilist[ftype];
+            for (i = 0; i < ilb.size(); i += 3)
             {
-                fc = ip[ilb->iatoms[i]].harmonic.krA;
-                re = ip[ilb->iatoms[i]].harmonic.rA;
+                fc = ip[ilb.iatoms[i]].harmonic.krA;
+                re = ip[ilb.iatoms[i]].harmonic.rA;
                 if (ftype == F_G96BONDS)
                 {
                     /* Convert squared sqaure fc to harmonic fc */
                     fc = 2*fc*re;
                 }
-                a1 = ilb->iatoms[i+1];
-                a2 = ilb->iatoms[i+2];
+                a1 = ilb.iatoms[i+1];
+                a2 = ilb.iatoms[i+2];
                 m1 = atom[a1].m;
                 m2 = atom[a2].m;
                 if (fc > 0 && m1 > 0 && m2 > 0)
@@ -305,18 +305,18 @@ static void check_bonds_timestep(const gmx_mtop_t *mtop, double dt, warninp_t wi
                 if (period2 < limit2)
                 {
                     bFound = false;
-                    for (j = 0; j < ilc->size(); j += 3)
+                    for (j = 0; j < ilc.size(); j += 3)
                     {
-                        if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
-                            (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
+                        if ((ilc.iatoms[j+1] == a1 && ilc.iatoms[j+2] == a2) ||
+                            (ilc.iatoms[j+1] == a2 && ilc.iatoms[j+2] == a1))
                         {
                             bFound = true;
                         }
                     }
-                    for (j = 0; j < ils->size(); j += 4)
+                    for (j = 0; j < ils.size(); j += 4)
                     {
-                        if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
-                            (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
+                        if ((a1 == ils.iatoms[j+1] || a1 == ils.iatoms[j+2] || a1 == ils.iatoms[j+3]) &&
+                            (a2 == ils.iatoms[j+1] || a2 == ils.iatoms[j+2] || a2 == ils.iatoms[j+3]))
                         {
                             bFound = true;
                         }
@@ -1318,14 +1318,14 @@ static void checkForUnboundAtoms(const gmx_moltype_t     *molt,
             (interaction_function[ftype].flags & IF_CONSTRAINT) ||
             ftype == F_SETTLE)
         {
-            const InteractionList *il   = &molt->ilist[ftype];
+            const InteractionList &il   = molt->ilist[ftype];
             const int              nral = NRAL(ftype);
 
-            for (int i = 0; i < il->size(); i += 1 + nral)
+            for (int i = 0; i < il.size(); i += 1 + nral)
             {
                 for (int j = 0; j < nral; j++)
                 {
-                    count[il->iatoms[i + 1 + j]]++;
+                    count[il.iatoms[i + 1 + j]]++;
                 }
             }
         }
@@ -1395,8 +1395,8 @@ static bool haveDecoupledModeInMol(const gmx_moltype_t *molt,
         if (interaction_function[ftype].flags & IF_ATYPE)
         {
             const int              nral = NRAL(ftype);
-            const InteractionList *il   = &molt->ilist[ftype];
-            for (int i = 0; i < il->size(); i += 1 + nral)
+            const InteractionList &il   = molt->ilist[ftype];
+            for (int i = 0; i < il.size(); i += 1 + nral)
             {
                 /* Here we check for the mass difference between the atoms
                  * at both ends of the angle, that the atoms at the ends have
@@ -1412,9 +1412,9 @@ static bool haveDecoupledModeInMol(const gmx_moltype_t *molt,
                  * occur in "normal" molecules and it doesn't hurt running
                  * those with higher accuracy settings as well.
                  */
-                int a0 = il->iatoms[1 + i];
-                int a1 = il->iatoms[1 + i + 1];
-                int a2 = il->iatoms[1 + i + 2];
+                int a0 = il.iatoms[1 + i];
+                int a1 = il.iatoms[1 + i + 1];
+                int a2 = il.iatoms[1 + i + 2];
                 if ((atom[a0].m > atom[a2].m*massFactorThreshold ||
                      atom[a2].m > atom[a0].m*massFactorThreshold) &&
                     atomToConstraints.index[a0 + 1] - atomToConstraints.index[a0] == 1 &&
index 80a7f19f7b2f1ab69da819694b0fb28a1d08a042..69debfc58bb813e2e5b297792132027266ef8652 100644 (file)
@@ -234,17 +234,17 @@ static void get_vsite_masses(const gmx_moltype_t  *moltype,
     {
         if (IS_VSITE(ft))
         {
-            const InteractionList *il = &moltype->ilist[ft];
+            const InteractionList &il = moltype->ilist[ft];
 
-            for (i = 0; i < il->size(); i += 1+NRAL(ft))
+            for (i = 0; i < il.size(); i += 1+NRAL(ft))
             {
                 const t_iparams *ip;
                 real             inv_mass, coeff, m_aj;
                 int              a1;
 
-                ip = &ffparams->iparams[il->iatoms[i]];
+                ip = &ffparams->iparams[il.iatoms[i]];
 
-                a1 = il->iatoms[i+1];
+                a1 = il.iatoms[i+1];
 
                 if (ft != F_VSITEN)
                 {
@@ -258,7 +258,7 @@ static void get_vsite_masses(const gmx_moltype_t  *moltype,
                     assert(maxj <= 5);
                     for (j = 1; j < maxj; j++)
                     {
-                        int aj = il->iatoms[i + 1 + j];
+                        int aj = il.iatoms[i + 1 + j];
                         cam[j] = getMass(moltype->atoms, aj, setMassesToOne);
                         if (cam[j] == 0)
                         {
@@ -314,10 +314,10 @@ static void get_vsite_masses(const gmx_moltype_t  *moltype,
 
                     /* Exact */
                     inv_mass = 0;
-                    for (j = 0; j < 3*ffparams->iparams[il->iatoms[i]].vsiten.n; j += 3)
+                    for (j = 0; j < 3*ffparams->iparams[il.iatoms[i]].vsiten.n; j += 3)
                     {
-                        int aj = il->iatoms[i + j + 2];
-                        coeff  = ffparams->iparams[il->iatoms[i+j]].vsiten.a;
+                        int aj = il.iatoms[i + j + 2];
+                        coeff  = ffparams->iparams[il.iatoms[i+j]].vsiten.a;
                         if (moltype->atoms.atom[aj].ptype == eptVSite)
                         {
                             m_aj = vsite_m[aj];
@@ -383,13 +383,13 @@ static void get_verlet_buffer_atomtypes(const gmx_mtop_t      *mtop,
 
         for (ft = F_CONSTR; ft <= F_CONSTRNC; ft++)
         {
-            const InteractionList *il = &moltype.ilist[ft];
+            const InteractionList &il = moltype.ilist[ft];
 
-            for (i = 0; i < il->size(); i += 1+NRAL(ft))
+            for (i = 0; i < il.size(); i += 1+NRAL(ft))
             {
-                ip         = &mtop->ffparams.iparams[il->iatoms[i]];
-                a1         = il->iatoms[i+1];
-                a2         = il->iatoms[i+2];
+                ip         = &mtop->ffparams.iparams[il.iatoms[i]];
+                a1         = il.iatoms[i+1];
+                a2         = il.iatoms[i+2];
                 real mass1 = getMass(*atoms, a1, setMassesToOne);
                 real mass2 = getMass(*atoms, a2, setMassesToOne);
                 if (mass2 > prop[a1].con_mass)
@@ -405,14 +405,14 @@ static void get_verlet_buffer_atomtypes(const gmx_mtop_t      *mtop,
             }
         }
 
-        const InteractionList *il = &moltype.ilist[F_SETTLE];
+        const InteractionList &il = moltype.ilist[F_SETTLE];
 
-        for (i = 0; i < il->size(); i += 1+NRAL(F_SETTLE))
+        for (i = 0; i < il.size(); i += 1+NRAL(F_SETTLE))
         {
-            ip         = &mtop->ffparams.iparams[il->iatoms[i]];
-            a1         = il->iatoms[i+1];
-            a2         = il->iatoms[i+2];
-            a3         = il->iatoms[i+3];
+            ip         = &mtop->ffparams.iparams[il.iatoms[i]];
+            a1         = il.iatoms[i+1];
+            a2         = il.iatoms[i+2];
+            a3         = il.iatoms[i+3];
             /* Usually the mass of a1 (usually oxygen) is larger than a2/a3.
              * If this is not the case, we overestimate the displacement,
              * which leads to a larger buffer (ok since this is an exotic case).
index 5f201efb8877883e7c26c6294689a8d9c3b1d9af..2274ef87fcdbef1510add66e35d6a01a62bd2fd6 100644 (file)
@@ -309,50 +309,47 @@ void pr_iparams(FILE *fp, t_functype ftype, const t_iparams *iparams)
 template <typename T>
 static void
 printIlist(FILE *fp, int indent, const char *title,
-           const t_functype *functype, const T *ilist,
+           const t_functype *functype, const T &ilist,
            gmx_bool bShowNumbers,
            gmx_bool bShowParameters, const t_iparams *iparams)
 {
     int      i, j, k, type, ftype;
 
-    if (available(fp, ilist, indent, title) && ilist->size() > 0)
+    indent = pr_title(fp, indent, title);
+    pr_indent(fp, indent);
+    fprintf(fp, "nr: %d\n", ilist.size());
+    if (ilist.size() > 0)
     {
-        indent = pr_title(fp, indent, title);
         pr_indent(fp, indent);
-        fprintf(fp, "nr: %d\n", ilist->size());
-        if (ilist->size() > 0)
+        fprintf(fp, "iatoms:\n");
+        for (i = j = 0; i < ilist.size(); )
         {
-            pr_indent(fp, indent);
-            fprintf(fp, "iatoms:\n");
-            for (i = j = 0; i < ilist->size(); )
+            pr_indent(fp, indent+INDENT);
+            type  = ilist.iatoms[i];
+            ftype = functype[type];
+            if (bShowNumbers)
             {
-                pr_indent(fp, indent+INDENT);
-                type  = ilist->iatoms[i];
-                ftype = functype[type];
-                if (bShowNumbers)
-                {
-                    fprintf(fp, "%d type=%d ", j, type);
-                }
-                j++;
-                printf("(%s)", interaction_function[ftype].name);
-                for (k = 0; k < interaction_function[ftype].nratoms; k++)
-                {
-                    fprintf(fp, " %3d", ilist->iatoms[i + 1 + k]);
-                }
-                if (bShowParameters)
-                {
-                    fprintf(fp, "  ");
-                    pr_iparams(fp, ftype,  &iparams[type]);
-                }
-                fprintf(fp, "\n");
-                i += 1+interaction_function[ftype].nratoms;
+                fprintf(fp, "%d type=%d ", j, type);
+            }
+            j++;
+            printf("(%s)", interaction_function[ftype].name);
+            for (k = 0; k < interaction_function[ftype].nratoms; k++)
+            {
+                fprintf(fp, " %3d", ilist.iatoms[i + 1 + k]);
             }
+            if (bShowParameters)
+            {
+                fprintf(fp, "  ");
+                pr_iparams(fp, ftype,  &iparams[type]);
+            }
+            fprintf(fp, "\n");
+            i += 1+interaction_function[ftype].nratoms;
         }
     }
 }
 
 void pr_ilist(FILE *fp, int indent, const char *title,
-              const t_functype *functype, const InteractionList *ilist,
+              const t_functype *functype, const InteractionList &ilist,
               gmx_bool bShowNumbers,
               gmx_bool bShowParameters, const t_iparams *iparams)
 {
@@ -448,7 +445,7 @@ void pr_idef(FILE *fp, int indent, const char *title, const t_idef *idef,
         for (j = 0; (j < F_NRE); j++)
         {
             printIlist(fp, indent, interaction_function[j].longname,
-                       idef->functype, &idef->il[j], bShowNumbers,
+                       idef->functype, idef->il[j], bShowNumbers,
                        bShowParameters, idef->iparams);
         }
     }
index 6f580f97aa0070e59187f092830d734ee506ba10..f4d4f95614b2ba00d7e3393793ccc4250e10cd50 100644 (file)
@@ -427,7 +427,7 @@ typedef struct t_idef
 
 void pr_iparams(FILE *fp, t_functype ftype, const t_iparams *iparams);
 void pr_ilist(FILE *fp, int indent, const char *title,
-              const t_functype *functype, const InteractionList *ilist,
+              const t_functype *functype, const InteractionList &ilist,
               gmx_bool bShowNumbers,
               gmx_bool bShowParameters, const t_iparams *iparams);
 void pr_ffparams(FILE *fp, int indent, const char *title,
index e26c64d2211df3cabc04d396cb0f63627eccc1b1..df6e857c57a5eda44afa86c76ec6936d5835f99a 100644 (file)
@@ -384,7 +384,7 @@ static void pr_moltype(FILE *fp, int indent, const char *title,
     for (j = 0; (j < F_NRE); j++)
     {
         pr_ilist(fp, indent, interaction_function[j].longname,
-                 ffparams->functype, &molt->ilist[j],
+                 ffparams->functype, molt->ilist[j],
                  bShowNumbers, bShowParameters, ffparams->iparams);
     }
 }
@@ -432,7 +432,7 @@ void pr_mtop(FILE *fp, int indent, const char *title, const gmx_mtop_t *mtop,
             {
                 pr_ilist(fp, indent, interaction_function[j].longname,
                          mtop->ffparams.functype,
-                         &(*mtop->intermolecular_ilist)[j],
+                         (*mtop->intermolecular_ilist)[j],
                          bShowNumbers, bShowParameters, mtop->ffparams.iparams);
             }
         }