Use more arrayref in listed forces sigatures
authorJoe Jordan <ejjordan12@gmail.com>
Fri, 19 Mar 2021 11:45:23 +0000 (11:45 +0000)
committerJoe Jordan <ejjordan12@gmail.com>
Fri, 19 Mar 2021 11:45:23 +0000 (11:45 +0000)
Pass members of t_mdatoms directly in many places in listed_forces.
Part of work on making refactor of t_mdatoms easier.

src/gromacs/gmxana/gmx_disre.cpp
src/gromacs/listed_forces/bonded.cpp
src/gromacs/listed_forces/bonded.h
src/gromacs/listed_forces/disre.cpp
src/gromacs/listed_forces/disre.h
src/gromacs/listed_forces/listed_forces.cpp
src/gromacs/listed_forces/orires.cpp
src/gromacs/listed_forces/orires.h
src/gromacs/listed_forces/pairs.cpp
src/gromacs/listed_forces/pairs.h
src/gromacs/listed_forces/tests/bonded.cpp

index b2b2728b1dd165a24cb22f5475d4f928697fa02c..5f32405b245bf9ab134edf680966acfd041b46c2 100644 (file)
@@ -237,7 +237,7 @@ static void check_viol(FILE*                          log,
         dr[clust_id].aver_6[ndr] += disresdata->Rt_6[label];
 
         snew(fshift, SHIFTS);
-        ta_disres(n, &forceatoms[i], forceparams.data(), x, f, fshift, pbc, lam, &dvdl, nullptr, nullptr, disresdata, nullptr, nullptr);
+        ta_disres(n, &forceatoms[i], forceparams.data(), x, f, fshift, pbc, lam, &dvdl, {}, nullptr, disresdata, nullptr, nullptr);
         sfree(fshift);
         viol = disresdata->sumviol;
 
index 5460db9b6a64e330afb4753ec3ab23ffbf8b3637..fbddaf5621daaef04316b65283f532f36cd0d54d 100644 (file)
@@ -70,6 +70,7 @@
 #include "gromacs/simd/simd.h"
 #include "gromacs/simd/simd_math.h"
 #include "gromacs/simd/vector_operations.h"
+#include "gromacs/utility/arrayref.h"
 #include "gromacs/utility/basedefinitions.h"
 #include "gromacs/utility/enumerationhelpers.h"
 #include "gromacs/utility/fatalerror.h"
@@ -91,20 +92,20 @@ namespace
 {
 
 //! Type of CPU function to compute a bonded interaction.
-using BondedFunction = real (*)(int              nbonds,
-                                const t_iatom    iatoms[],
-                                const t_iparams  iparams[],
-                                const rvec       x[],
-                                rvec4            f[],
-                                rvec             fshift[],
-                                const t_pbc*     pbc,
-                                real             lambda,
-                                real*            dvdlambda,
-                                const t_mdatoms* md,
-                                t_fcdata*        fcd,
-                                t_disresdata*    disresdata,
-                                t_oriresdata*    oriresdata,
-                                int*             ddgatindex);
+using BondedFunction = real (*)(int                       nbonds,
+                                const t_iatom             iatoms[],
+                                const t_iparams           iparams[],
+                                const rvec                x[],
+                                rvec4                     f[],
+                                rvec                      fshift[],
+                                const t_pbc*              pbc,
+                                real                      lambda,
+                                real*                     dvdlambda,
+                                gmx::ArrayRef<const real> charge,
+                                t_fcdata*                 fcd,
+                                t_disresdata*             disresdata,
+                                t_oriresdata*             oriresdata,
+                                int*                      ddgatindex);
 
 /*! \brief Mysterious CMAP coefficient matrix */
 const int cmap_coeff_matrix[] = {
@@ -256,7 +257,7 @@ real morse_bonds(int             nbonds,
                  const t_pbc*    pbc,
                  real            lambda,
                  real*           dvdlambda,
-                 const t_mdatoms gmx_unused* md,
+                 gmx::ArrayRef<const real> /*charge*/,
                  t_fcdata gmx_unused* fcd,
                  t_disresdata gmx_unused* disresdata,
                  t_oriresdata gmx_unused* oriresdata,
@@ -326,7 +327,7 @@ real cubic_bonds(int             nbonds,
                  const t_pbc*    pbc,
                  real gmx_unused lambda,
                  real gmx_unused* dvdlambda,
-                 const t_mdatoms gmx_unused* md,
+                 gmx::ArrayRef<const real> /*charge*/,
                  t_fcdata gmx_unused* fcd,
                  t_disresdata gmx_unused* disresdata,
                  t_oriresdata gmx_unused* oriresdata,
@@ -383,7 +384,7 @@ real FENE_bonds(int             nbonds,
                 const t_pbc*    pbc,
                 real gmx_unused lambda,
                 real gmx_unused* dvdlambda,
-                const t_mdatoms gmx_unused* md,
+                gmx::ArrayRef<const real> /*charge*/,
                 t_fcdata gmx_unused* fcd,
                 t_disresdata gmx_unused* disresdata,
                 t_oriresdata gmx_unused* oriresdata,
@@ -474,7 +475,7 @@ bonds(int             nbonds,
       const t_pbc*    pbc,
       real            lambda,
       real*           dvdlambda,
-      const t_mdatoms gmx_unused* md,
+      gmx::ArrayRef<const real> /*charge*/,
       t_fcdata gmx_unused* fcd,
       t_disresdata gmx_unused* disresdata,
       t_oriresdata gmx_unused* oriresdata,
@@ -536,7 +537,7 @@ bonds(int             nbonds,
       const t_pbc*    pbc,
       real gmx_unused lambda,
       real gmx_unused* dvdlambda,
-      const t_mdatoms gmx_unused* md,
+      gmx::ArrayRef<const real> /*charge*/,
       t_fcdata gmx_unused* fcd,
       t_disresdata gmx_unused* disresdata,
       t_oriresdata gmx_unused* oriresdata,
@@ -628,7 +629,7 @@ real restraint_bonds(int             nbonds,
                      const t_pbc*    pbc,
                      real            lambda,
                      real*           dvdlambda,
-                     const t_mdatoms gmx_unused* md,
+                     gmx::ArrayRef<const real> /*charge*/,
                      t_fcdata gmx_unused* fcd,
                      t_disresdata gmx_unused* disresdata,
                      t_oriresdata gmx_unused* oriresdata,
@@ -709,16 +710,16 @@ real restraint_bonds(int             nbonds,
 }
 
 template<BondedKernelFlavor flavor>
-real polarize(int              nbonds,
-              const t_iatom    forceatoms[],
-              const t_iparams  forceparams[],
-              const rvec       x[],
-              rvec4            f[],
-              rvec             fshift[],
-              const t_pbc*     pbc,
-              real             lambda,
-              real*            dvdlambda,
-              const t_mdatoms* md,
+real polarize(int                       nbonds,
+              const t_iatom             forceatoms[],
+              const t_iparams           forceparams[],
+              const rvec                x[],
+              rvec4                     f[],
+              rvec                      fshift[],
+              const t_pbc*              pbc,
+              real                      lambda,
+              real*                     dvdlambda,
+              gmx::ArrayRef<const real> charge,
               t_fcdata gmx_unused* fcd,
               t_disresdata gmx_unused* disresdata,
               t_oriresdata gmx_unused* oriresdata,
@@ -734,7 +735,7 @@ real polarize(int              nbonds,
         type = forceatoms[i++];
         ai   = forceatoms[i++];
         aj   = forceatoms[i++];
-        ksh  = gmx::square(md->chargeA[aj]) * gmx::c_one4PiEps0 / forceparams[type].polarize.alpha;
+        ksh  = gmx::square(charge[aj]) * gmx::c_one4PiEps0 / forceparams[type].polarize.alpha;
 
         ki  = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /*   3      */
         dr2 = iprod(dx, dx);                       /*   5              */
@@ -756,16 +757,16 @@ real polarize(int              nbonds,
 }
 
 template<BondedKernelFlavor flavor>
-real anharm_polarize(int              nbonds,
-                     const t_iatom    forceatoms[],
-                     const t_iparams  forceparams[],
-                     const rvec       x[],
-                     rvec4            f[],
-                     rvec             fshift[],
-                     const t_pbc*     pbc,
-                     real             lambda,
-                     real*            dvdlambda,
-                     const t_mdatoms* md,
+real anharm_polarize(int                       nbonds,
+                     const t_iatom             forceatoms[],
+                     const t_iparams           forceparams[],
+                     const rvec                x[],
+                     rvec4                     f[],
+                     rvec                      fshift[],
+                     const t_pbc*              pbc,
+                     real                      lambda,
+                     real*                     dvdlambda,
+                     gmx::ArrayRef<const real> charge,
                      t_fcdata gmx_unused* fcd,
                      t_disresdata gmx_unused* disresdata,
                      t_oriresdata gmx_unused* oriresdata,
@@ -781,7 +782,7 @@ real anharm_polarize(int              nbonds,
         type = forceatoms[i++];
         ai   = forceatoms[i++];
         aj   = forceatoms[i++];
-        ksh = gmx::square(md->chargeA[aj]) * gmx::c_one4PiEps0 / forceparams[type].anharm_polarize.alpha; /* 7*/
+        ksh = gmx::square(charge[aj]) * gmx::c_one4PiEps0 / forceparams[type].anharm_polarize.alpha; /* 7*/
         khyp  = forceparams[type].anharm_polarize.khyp;
         drcut = forceparams[type].anharm_polarize.drcut;
 
@@ -820,8 +821,8 @@ real water_pol(int             nbonds,
                rvec gmx_unused fshift[],
                const t_pbc gmx_unused* pbc,
                real gmx_unused lambda,
-               real gmx_unused* dvdlambda,
-               const t_mdatoms gmx_unused* md,
+               real gmx_unused*          dvdlambda,
+               gmx::ArrayRef<const real> charge,
                t_fcdata gmx_unused* fcd,
                t_disresdata gmx_unused* disresdata,
                t_oriresdata gmx_unused* oriresdata,
@@ -840,7 +841,7 @@ real water_pol(int             nbonds,
     {
         type0  = forceatoms[0];
         aS     = forceatoms[5];
-        qS     = md->chargeA[aS];
+        qS     = charge[aS];
         kk[XX] = gmx::square(qS) * gmx::c_one4PiEps0 / forceparams[type0].wpol.al_x;
         kk[YY] = gmx::square(qS) * gmx::c_one4PiEps0 / forceparams[type0].wpol.al_y;
         kk[ZZ] = gmx::square(qS) * gmx::c_one4PiEps0 / forceparams[type0].wpol.al_z;
@@ -965,8 +966,8 @@ real thole_pol(int             nbonds,
                rvec            fshift[],
                const t_pbc*    pbc,
                real gmx_unused lambda,
-               real gmx_unused* dvdlambda,
-               const t_mdatoms* md,
+               real gmx_unused*          dvdlambda,
+               gmx::ArrayRef<const real> charge,
                t_fcdata gmx_unused* fcd,
                t_disresdata gmx_unused* disresdata,
                t_oriresdata gmx_unused* oriresdata,
@@ -984,8 +985,8 @@ real thole_pol(int             nbonds,
         da1  = forceatoms[i++];
         a2   = forceatoms[i++];
         da2  = forceatoms[i++];
-        q1   = md->chargeA[da1];
-        q2   = md->chargeA[da2];
+        q1   = charge[da1];
+        q2   = charge[da2];
         a    = forceparams[type].thole.a;
         al1  = forceparams[type].thole.alpha1;
         al2  = forceparams[type].thole.alpha2;
@@ -1019,7 +1020,7 @@ angles(int             nbonds,
        const t_pbc*    pbc,
        real            lambda,
        real*           dvdlambda,
-       const t_mdatoms gmx_unused* md,
+       gmx::ArrayRef<const real> /*charge*/,
        t_fcdata gmx_unused* fcd,
        t_disresdata gmx_unused* disresdata,
        t_oriresdata gmx_unused* oriresdata,
@@ -1113,7 +1114,7 @@ angles(int             nbonds,
        const t_pbc*    pbc,
        real gmx_unused lambda,
        real gmx_unused* dvdlambda,
-       const t_mdatoms gmx_unused* md,
+       gmx::ArrayRef<const real> /*charge*/,
        t_fcdata gmx_unused* fcd,
        t_disresdata gmx_unused* disresdata,
        t_oriresdata gmx_unused* oriresdata,
@@ -1273,7 +1274,7 @@ real linear_angles(int             nbonds,
                    const t_pbc*    pbc,
                    real            lambda,
                    real*           dvdlambda,
-                   const t_mdatoms gmx_unused* md,
+                   gmx::ArrayRef<const real> /*charge*/,
                    t_fcdata gmx_unused* fcd,
                    t_disresdata gmx_unused* disresdata,
                    t_oriresdata gmx_unused* oriresdata,
@@ -1345,7 +1346,7 @@ urey_bradley(int             nbonds,
              const t_pbc*    pbc,
              real            lambda,
              real*           dvdlambda,
-             const t_mdatoms gmx_unused* md,
+             gmx::ArrayRef<const real> /*charge*/,
              t_fcdata gmx_unused* fcd,
              t_disresdata gmx_unused* disresdata,
              t_oriresdata gmx_unused* oriresdata,
@@ -1457,7 +1458,7 @@ urey_bradley(int             nbonds,
              const t_pbc*    pbc,
              real gmx_unused lambda,
              real gmx_unused* dvdlambda,
-             const t_mdatoms gmx_unused* md,
+             gmx::ArrayRef<const real> /*charge*/,
              t_fcdata gmx_unused* fcd,
              t_disresdata gmx_unused* disresdata,
              t_oriresdata gmx_unused* oriresdata,
@@ -1601,7 +1602,7 @@ real quartic_angles(int             nbonds,
                     const t_pbc*    pbc,
                     real gmx_unused lambda,
                     real gmx_unused* dvdlambda,
-                    const t_mdatoms gmx_unused* md,
+                    gmx::ArrayRef<const real> /*charge*/,
                     t_fcdata gmx_unused* fcd,
                     t_disresdata gmx_unused* disresdata,
                     t_oriresdata gmx_unused* oriresdata,
@@ -1968,7 +1969,7 @@ pdihs(int             nbonds,
       const t_pbc*    pbc,
       real            lambda,
       real*           dvdlambda,
-      const t_mdatoms gmx_unused* md,
+      gmx::ArrayRef<const real> /*charge*/,
       t_fcdata gmx_unused* fcd,
       t_disresdata gmx_unused* disresdata,
       t_oriresdata gmx_unused* oriresdata,
@@ -2031,7 +2032,7 @@ pdihs(int             nbonds,
       const t_pbc*    pbc,
       real gmx_unused lambda,
       real gmx_unused* dvdlambda,
-      const t_mdatoms gmx_unused* md,
+      gmx::ArrayRef<const real> /*charge*/,
       t_fcdata gmx_unused* fcd,
       t_disresdata gmx_unused* disresdata,
       t_oriresdata gmx_unused* oriresdata,
@@ -2147,7 +2148,7 @@ rbdihs(int             nbonds,
        const t_pbc*    pbc,
        real gmx_unused lambda,
        real gmx_unused* dvdlambda,
-       const t_mdatoms gmx_unused* md,
+       gmx::ArrayRef<const real> /*charge*/,
        t_fcdata gmx_unused* fcd,
        t_disresdata gmx_unused* disresdata,
        t_oriresdata gmx_unused* oriresdata,
@@ -2275,7 +2276,7 @@ real idihs(int             nbonds,
            const t_pbc*    pbc,
            real            lambda,
            real*           dvdlambda,
-           const t_mdatoms gmx_unused* md,
+           gmx::ArrayRef<const real> /*charge*/,
            t_fcdata gmx_unused* fcd,
            t_disresdata gmx_unused* disresdata,
            t_oriresdata gmx_unused* oriresdata,
@@ -2443,7 +2444,7 @@ real angres(int             nbonds,
             const t_pbc*    pbc,
             real            lambda,
             real*           dvdlambda,
-            const t_mdatoms gmx_unused* md,
+            gmx::ArrayRef<const real> /*charge*/,
             t_fcdata gmx_unused* fcd,
             t_disresdata gmx_unused* disresdata,
             t_oriresdata gmx_unused* oriresdata,
@@ -2462,7 +2463,7 @@ real angresz(int             nbonds,
              const t_pbc*    pbc,
              real            lambda,
              real*           dvdlambda,
-             const t_mdatoms gmx_unused* md,
+             gmx::ArrayRef<const real> /*charge*/,
              t_fcdata gmx_unused* fcd,
              t_disresdata gmx_unused* disresdata,
              t_oriresdata gmx_unused* oriresdata,
@@ -2481,7 +2482,7 @@ real dihres(int             nbonds,
             const t_pbc*    pbc,
             real            lambda,
             real*           dvdlambda,
-            const t_mdatoms gmx_unused* md,
+            gmx::ArrayRef<const real> /*charge*/,
             t_fcdata gmx_unused* fcd,
             t_disresdata gmx_unused* disresdata,
             t_oriresdata gmx_unused* oriresdata,
@@ -2576,7 +2577,7 @@ real unimplemented(int gmx_unused nbonds,
                    const t_pbc gmx_unused* pbc,
                    real gmx_unused lambda,
                    real gmx_unused* dvdlambda,
-                   const t_mdatoms gmx_unused* md,
+                   gmx::ArrayRef<const real> /*charge*/,
                    t_fcdata gmx_unused* fcd,
                    t_disresdata gmx_unused* disresdata,
                    t_oriresdata gmx_unused* oriresdata,
@@ -2595,7 +2596,7 @@ real restrangles(int             nbonds,
                  const t_pbc*    pbc,
                  real gmx_unused lambda,
                  real gmx_unused* dvdlambda,
-                 const t_mdatoms gmx_unused* md,
+                 gmx::ArrayRef<const real> /*charge*/,
                  t_fcdata gmx_unused* fcd,
                  t_disresdata gmx_unused* disresdata,
                  t_oriresdata gmx_unused* oriresdata,
@@ -2699,7 +2700,7 @@ real restrdihs(int             nbonds,
                const t_pbc*    pbc,
                real gmx_unused lambda,
                real gmx_unused* dvlambda,
-               const t_mdatoms gmx_unused* md,
+               gmx::ArrayRef<const real> /*charge*/,
                t_fcdata gmx_unused* fcd,
                t_disresdata gmx_unused* disresdata,
                t_oriresdata gmx_unused* oriresdata,
@@ -2823,7 +2824,7 @@ real cbtdihs(int             nbonds,
              const t_pbc*    pbc,
              real gmx_unused lambda,
              real gmx_unused* dvdlambda,
-             const t_mdatoms gmx_unused* md,
+             gmx::ArrayRef<const real> /*charge*/,
              t_fcdata gmx_unused* fcd,
              t_disresdata gmx_unused* disresdata,
              t_oriresdata gmx_unused* oriresdata,
@@ -2940,7 +2941,7 @@ rbdihs(int             nbonds,
        const t_pbc*    pbc,
        real            lambda,
        real*           dvdlambda,
-       const t_mdatoms gmx_unused* md,
+       gmx::ArrayRef<const real> /*charge*/,
        t_fcdata gmx_unused* fcd,
        t_disresdata gmx_unused* disresdata,
        t_oriresdata gmx_unused* oriresdata,
@@ -3092,7 +3093,7 @@ real cmap_dihs(int                 nbonds,
                const struct t_pbc* pbc,
                real gmx_unused lambda,
                real gmx_unused* dvdlambda,
-               const t_mdatoms gmx_unused* md,
+               gmx::ArrayRef<const real> /*charge*/,
                t_fcdata gmx_unused* fcd,
                t_disresdata gmx_unused* disresdata,
                t_oriresdata gmx_unused* oriresdata,
@@ -3498,7 +3499,7 @@ real g96bonds(int             nbonds,
               const t_pbc*    pbc,
               real            lambda,
               real*           dvdlambda,
-              const t_mdatoms gmx_unused* md,
+              gmx::ArrayRef<const real> /*charge*/,
               t_fcdata gmx_unused* fcd,
               t_disresdata gmx_unused* disresdata,
               t_oriresdata gmx_unused* oriresdata,
@@ -3557,7 +3558,7 @@ real g96angles(int             nbonds,
                const t_pbc*    pbc,
                real            lambda,
                real*           dvdlambda,
-               const t_mdatoms gmx_unused* md,
+               gmx::ArrayRef<const real> /*charge*/,
                t_fcdata gmx_unused* fcd,
                t_disresdata gmx_unused* disresdata,
                t_oriresdata gmx_unused* oriresdata,
@@ -3626,7 +3627,7 @@ real cross_bond_bond(int             nbonds,
                      const t_pbc*    pbc,
                      real gmx_unused lambda,
                      real gmx_unused* dvdlambda,
-                     const t_mdatoms gmx_unused* md,
+                     gmx::ArrayRef<const real> /*charge*/,
                      t_fcdata gmx_unused* fcd,
                      t_disresdata gmx_unused* disresdata,
                      t_oriresdata gmx_unused* oriresdata,
@@ -3700,7 +3701,7 @@ real cross_bond_angle(int             nbonds,
                       const t_pbc*    pbc,
                       real gmx_unused lambda,
                       real gmx_unused* dvdlambda,
-                      const t_mdatoms gmx_unused* md,
+                      gmx::ArrayRef<const real> /*charge*/,
                       t_fcdata gmx_unused* fcd,
                       t_disresdata gmx_unused* disresdata,
                       t_oriresdata gmx_unused* oriresdata,
@@ -3838,8 +3839,8 @@ real tab_bonds(int             nbonds,
                const t_pbc*    pbc,
                real            lambda,
                real*           dvdlambda,
-               const t_mdatoms gmx_unused* md,
-               t_fcdata*                   fcd,
+               gmx::ArrayRef<const real> /*charge*/,
+               t_fcdata*    fcd,
                t_disresdata gmx_unused* disresdata,
                t_oriresdata gmx_unused* oriresdata,
                int gmx_unused* global_atom_index)
@@ -3895,8 +3896,8 @@ real tab_angles(int             nbonds,
                 const t_pbc*    pbc,
                 real            lambda,
                 real*           dvdlambda,
-                const t_mdatoms gmx_unused* md,
-                t_fcdata*                   fcd,
+                gmx::ArrayRef<const real> /*charge*/,
+                t_fcdata*    fcd,
                 t_disresdata gmx_unused* disresdata,
                 t_oriresdata gmx_unused* oriresdata,
                 int gmx_unused* global_atom_index)
@@ -3977,8 +3978,8 @@ real tab_dihs(int             nbonds,
               const t_pbc*    pbc,
               real            lambda,
               real*           dvdlambda,
-              const t_mdatoms gmx_unused* md,
-              t_fcdata*                   fcd,
+              gmx::ArrayRef<const real> /*charge*/,
+              t_fcdata*    fcd,
               t_disresdata gmx_unused* disresdata,
               t_oriresdata gmx_unused* oriresdata,
               int gmx_unused* global_atom_index)
@@ -4145,27 +4146,27 @@ gmx::EnumerationArray<BondedKernelFlavor, std::array<BondedInteractions, F_NRE>>
 
 } // namespace
 
-real calculateSimpleBond(const int           ftype,
-                         const int           numForceatoms,
-                         const t_iatom       forceatoms[],
-                         const t_iparams     forceparams[],
-                         const rvec          x[],
-                         rvec4               f[],
-                         rvec                fshift[],
-                         const struct t_pbc* pbc,
-                         const real          lambda,
-                         real*               dvdlambda,
-                         const t_mdatoms*    md,
-                         t_fcdata*           fcd,
-                         t_disresdata*       disresdata,
-                         t_oriresdata*       oriresdata,
+real calculateSimpleBond(const int                 ftype,
+                         const int                 numForceatoms,
+                         const t_iatom             forceatoms[],
+                         const t_iparams           forceparams[],
+                         const rvec                x[],
+                         rvec4                     f[],
+                         rvec                      fshift[],
+                         const struct t_pbc*       pbc,
+                         const real                lambda,
+                         real*                     dvdlambda,
+                         gmx::ArrayRef<const real> charge,
+                         t_fcdata*                 fcd,
+                         t_disresdata*             disresdata,
+                         t_oriresdata*             oriresdata,
                          int gmx_unused*          global_atom_index,
                          const BondedKernelFlavor bondedKernelFlavor)
 {
     const BondedInteractions& bonded = c_bondedInteractionFunctionsPerFlavor[bondedKernelFlavor][ftype];
 
     real v = bonded.function(
-            numForceatoms, forceatoms, forceparams, x, f, fshift, pbc, lambda, dvdlambda, md, fcd, disresdata, oriresdata, global_atom_index);
+            numForceatoms, forceatoms, forceparams, x, f, fshift, pbc, lambda, dvdlambda, charge, fcd, disresdata, oriresdata, global_atom_index);
 
     return v;
 }
index d3c006145e1070282bc78e585f6961391eed46e5..7a2b6282b91b064588fbe2441aae6f73ea89cc48 100644 (file)
@@ -56,7 +56,6 @@
 
 struct gmx_cmap_t;
 struct t_fcdata;
-struct t_mdatom;
 struct t_nrnb;
 struct t_pbc;
 struct t_disresdata;
@@ -66,6 +65,8 @@ namespace gmx
 {
 template<typename EnumType, typename DataType, EnumType ArraySize>
 struct EnumerationArray;
+template<typename>
+class ArrayRef;
 } // namespace gmx
 
 /*! \brief Calculate bond-angle. No PBC is taken into account (use mol-shift) */
@@ -127,7 +128,7 @@ real cmap_dihs(int                 nbonds,
                const struct t_pbc* pbc,
                real gmx_unused lambda,
                real gmx_unused* dvdlambda,
-               const t_mdatoms gmx_unused* md,
+               gmx::ArrayRef<const real> /*charge*/,
                t_fcdata gmx_unused* fcd,
                t_disresdata gmx_unused* disresdata,
                t_oriresdata gmx_unused* oriresdata,
@@ -172,20 +173,20 @@ static constexpr inline bool computeEnergyOrVirial(const BondedKernelFlavor flav
  * All pointers should be non-null, except for pbc and g which can be nullptr.
  * \returns the energy or 0 when \p bondedKernelFlavor did not request the energy.
  */
-real calculateSimpleBond(int                 ftype,
-                         int                 numForceatoms,
-                         const t_iatom       forceatoms[],
-                         const t_iparams     forceparams[],
-                         const rvec          x[],
-                         rvec4               f[],
-                         rvec                fshift[],
-                         const struct t_pbc* pbc,
-                         real                lambda,
-                         real*               dvdlambda,
-                         const t_mdatoms*    md,
-                         t_fcdata*           fcd,
-                         t_disresdata*       disresdata,
-                         t_oriresdata*       oriresdata,
+real calculateSimpleBond(int                       ftype,
+                         int                       numForceatoms,
+                         const t_iatom             forceatoms[],
+                         const t_iparams           forceparams[],
+                         const rvec                x[],
+                         rvec4                     f[],
+                         rvec                      fshift[],
+                         const struct t_pbc*       pbc,
+                         real                      lambda,
+                         real*                     dvdlambda,
+                         gmx::ArrayRef<const real> charge,
+                         t_fcdata*                 fcd,
+                         t_disresdata*             disresdata,
+                         t_oriresdata*             oriresdata,
                          int gmx_unused*    global_atom_index,
                          BondedKernelFlavor bondedKernelFlavor);
 
index 89c48d588c637dd2c0aaedfd2a06d85852a05e6d..b5630b4aeb1cd49782dd0eb8c6af80fb7244fb0c 100644 (file)
 #include "gromacs/mdtypes/commrec.h"
 #include "gromacs/mdtypes/fcdata.h"
 #include "gromacs/mdtypes/inputrec.h"
+#include "gromacs/mdtypes/mdatom.h"
 #include "gromacs/mdtypes/md_enums.h"
 #include "gromacs/mdtypes/state.h"
 #include "gromacs/pbcutil/ishift.h"
 #include "gromacs/pbcutil/pbc.h"
 #include "gromacs/topology/mtop_util.h"
 #include "gromacs/topology/topology.h"
+#include "gromacs/utility/arrayref.h"
 #include "gromacs/utility/fatalerror.h"
 #include "gromacs/utility/futil.h"
 #include "gromacs/utility/gmxassert.h"
@@ -424,7 +426,7 @@ real ta_disres(int              nfa,
                const t_pbc*     pbc,
                real gmx_unused lambda,
                real gmx_unused* dvdlambda,
-               const t_mdatoms gmx_unused* md,
+               gmx::ArrayRef<const real> /*charge*/,
                t_fcdata gmx_unused* fcd,
                t_disresdata*        disresdata,
                t_oriresdata gmx_unused* oriresdata,
index 756cce341c207800fa99a80f7abf9badb25670cb..7d885ab944d7ae810ec9fed17a9cd059f1600221 100644 (file)
@@ -64,6 +64,12 @@ class t_state;
 enum class DDRole;
 enum class NumRanks;
 
+namespace gmx
+{
+template<typename>
+class ArrayRef;
+} // namespace gmx
+
 //! Whether distance restraints are called from mdrun or from an analysis tool
 enum class DisResRunMode
 {
@@ -107,16 +113,16 @@ void calc_disres_R_6(const t_commrec*      cr,
                      const history_t*      hist);
 
 //! Calculates the distance restraint forces, return the potential.
-real ta_disres(int              nfa,
-               const t_iatom*   forceatoms,
-               const t_iparams* ip,
-               const rvec*      x,
-               rvec4*           f,
-               rvec*            fshift,
-               const t_pbc*     pbc,
-               real             lambda,
-               real*            dvdlambda,
-               const t_mdatoms* md,
+real ta_disres(int                       nfa,
+               const t_iatom*            forceatoms,
+               const t_iparams*          ip,
+               const rvec*               x,
+               rvec4*                    f,
+               rvec*                     fshift,
+               const t_pbc*              pbc,
+               real                      lambda,
+               real*                     dvdlambda,
+               gmx::ArrayRef<const real> charge,
                t_fcdata gmx_unused* fcd,
                t_disresdata*        disresdata,
                t_oriresdata gmx_unused* oriresdata,
index e8b92fc82b7dfa408f5df69f74cad165d23d643d..82b93fa3e2d256c89bee174e3a0792d876ff7dee 100644 (file)
@@ -66,6 +66,7 @@
 #include "gromacs/mdlib/force.h"
 #include "gromacs/mdlib/gmx_omp_nthreads.h"
 #include "gromacs/mdtypes/commrec.h"
+#include "gromacs/mdtypes/mdatom.h"
 #include "gromacs/mdtypes/fcdata.h"
 #include "gromacs/mdtypes/forceoutput.h"
 #include "gromacs/mdtypes/forcerec.h"
@@ -470,7 +471,7 @@ real calc_one_bond(int                           thread,
                           pbc,
                           lambda[static_cast<int>(efptFTYPE)],
                           &(dvdl[static_cast<int>(efptFTYPE)]),
-                          md,
+                          gmx::arrayRefFromArray(md->chargeA, md->nr),
                           fcd,
                           nullptr,
                           nullptr,
@@ -488,7 +489,7 @@ real calc_one_bond(int                           thread,
                                     pbc,
                                     lambda[static_cast<int>(efptFTYPE)],
                                     &(dvdl[static_cast<int>(efptFTYPE)]),
-                                    md,
+                                    gmx::arrayRefFromArray(md->chargeA, md->nr),
                                     fcd,
                                     fcd->disres,
                                     fcd->orires,
index 0ad51f201e7998e1e3d9c3d78cc35c8c13669bd7..be62d1139b3b5bff96fec672073f000a521977e8 100644 (file)
@@ -446,14 +446,16 @@ real calc_orires_dev(const gmx_multisim_t* ms,
     }
 
     clear_rvec(com);
-    mtot  = 0;
-    int j = 0;
+    mtot       = 0;
+    int  j     = 0;
+    auto massT = md->massT;
+    auto cORF  = md->cORF;
     for (int i = 0; i < md->nr; i++)
     {
-        if (md->cORF[i] == 0)
+        if (cORF[i] == 0)
         {
             copy_rvec(xWholeMolecules[i], xtmp[j]);
-            mref[j] = md->massT[i];
+            mref[j] = massT[i];
             for (int d = 0; d < DIM; d++)
             {
                 com[d] += mref[j] * xtmp[j][d];
@@ -657,7 +659,7 @@ real orires(int             nfa,
             const t_pbc*    pbc,
             real gmx_unused lambda,
             real gmx_unused* dvdlambda,
-            const t_mdatoms gmx_unused* md,
+            gmx::ArrayRef<const real> /*charge*/,
             t_fcdata gmx_unused* fcd,
             t_disresdata gmx_unused* disresdata,
             t_oriresdata*            oriresdata,
index 7c269b02174db8960938e81e46d013d99dd8e64c..bff4c3102fc9606513c0e42cc1b27e0c115d3d86 100644 (file)
@@ -59,6 +59,7 @@ struct t_oriresdata;
 struct t_disresdata;
 struct t_fcdata;
 class t_state;
+struct t_mdatoms;
 
 namespace gmx
 {
@@ -110,20 +111,20 @@ void diagonalize_orires_tensors(t_oriresdata* od);
 void print_orires_log(FILE* log, t_oriresdata* od);
 
 //! Calculates the orientation restraint forces.
-real orires(int              nfa,
-            const t_iatom    forceatoms[],
-            const t_iparams  ip[],
-            const rvec       x[],
-            rvec4            f[],
-            rvec             fshift[],
-            const t_pbc*     pbc,
-            real             lambda,
-            real*            dvdlambda,
-            const t_mdatoms* md,
-            t_fcdata*        fcd,
-            t_disresdata*    disresdata,
-            t_oriresdata*    oriresdata,
-            int*             global_atom_index);
+real orires(int                       nfa,
+            const t_iatom             forceatoms[],
+            const t_iparams           ip[],
+            const rvec                x[],
+            rvec4                     f[],
+            rvec                      fshift[],
+            const t_pbc*              pbc,
+            real                      lambda,
+            real*                     dvdlambda,
+            gmx::ArrayRef<const real> charge,
+            t_fcdata*                 fcd,
+            t_disresdata*             disresdata,
+            t_oriresdata*             oriresdata,
+            int*                      global_atom_index);
 
 //! Copies the new time averages that have been calculated in calc_orires_dev.
 void update_orires_history(const t_oriresdata& oriresdata, history_t* hist);
index f84fe890abbe35410078b102d070f45a1a8844e0..43d292b6fc2bb166652e94b80fb4e6a9a4a7294d 100644 (file)
@@ -54,8 +54,8 @@
 #include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/group.h"
 #include "gromacs/mdtypes/interaction_const.h"
-#include "gromacs/mdtypes/md_enums.h"
 #include "gromacs/mdtypes/mdatom.h"
+#include "gromacs/mdtypes/md_enums.h"
 #include "gromacs/mdtypes/nblist.h"
 #include "gromacs/mdtypes/simulation_workload.h"
 #include "gromacs/pbcutil/ishift.h"
@@ -427,23 +427,27 @@ static real do_pairs_general(int                 ftype,
 
     const real epsfac = fr->ic->epsfac;
 
-    bFreeEnergy = FALSE;
+    bFreeEnergy     = FALSE;
+    auto cENER      = md->cENER;
+    auto bPerturbed = md->bPerturbed;
+    auto chargeA    = md->chargeA;
+    auto chargeB    = md->chargeB;
     for (i = 0; (i < nbonds);)
     {
         itype = iatoms[i++];
         ai    = iatoms[i++];
         aj    = iatoms[i++];
-        gid   = GID(md->cENER[ai], md->cENER[aj], md->nenergrp);
+        gid   = GID(cENER[ai], cENER[aj], md->nenergrp);
 
         /* Get parameters */
         switch (ftype)
         {
             case F_LJ14:
                 bFreeEnergy = (fr->efep != FreeEnergyPerturbationType::No
-                               && (((md->nPerturbed != 0) && (md->bPerturbed[ai] || md->bPerturbed[aj]))
+                               && ((md->bPerturbed && (bPerturbed[ai] || bPerturbed[aj]))
                                    || iparams[itype].lj14.c6A != iparams[itype].lj14.c6B
                                    || iparams[itype].lj14.c12A != iparams[itype].lj14.c12B));
-                qq          = md->chargeA[ai] * md->chargeA[aj] * epsfac * fr->fudgeQQ;
+                qq          = chargeA[ai] * chargeA[aj] * epsfac * fr->fudgeQQ;
                 c6          = iparams[itype].lj14.c6A;
                 c12         = iparams[itype].lj14.c12A;
                 break;
@@ -499,7 +503,7 @@ static real do_pairs_general(int                 ftype,
         if (bFreeEnergy)
         {
             /* Currently free energy is only supported for F_LJ14, so no need to check for that if we got here */
-            qqB  = md->chargeB[ai] * md->chargeB[aj] * epsfac * fr->fudgeQQ;
+            qqB  = chargeB[ai] * chargeB[aj] * epsfac * fr->fudgeQQ;
             c6B  = iparams[itype].lj14.c6B * 6.0;
             c12B = iparams[itype].lj14.c12B * 12.0;
 
@@ -588,6 +592,7 @@ static void do_pairs_simple(int              nbonds,
     std::int32_t aj[pack_size];
     real         coeff[3 * pack_size];
 #endif
+    auto chargeA = md->chargeA;
 
     /* nbonds is #pairs*nfa1, here we step pack_size pairs */
     for (int i = 0; i < nbonds; i += pack_size * nfa1)
@@ -606,7 +611,7 @@ static void do_pairs_simple(int              nbonds,
             {
                 coeff[0 * pack_size + s] = iparams[itype].lj14.c6A;
                 coeff[1 * pack_size + s] = iparams[itype].lj14.c12A;
-                coeff[2 * pack_size + s] = md->chargeA[ai[s]] * md->chargeA[aj[s]];
+                coeff[2 * pack_size + s] = chargeA[ai[s]] * chargeA[aj[s]];
 
                 /* Avoid indexing the iatoms array out of bounds.
                  * We pad the coordinate indices with the last atom pair.
index af9320e014af230b454861088594b11455fe788a..982c3f03b6b1ab4f35e514dc4e5d3e1bb63d882d 100644 (file)
@@ -2,7 +2,7 @@
  * This file is part of the GROMACS molecular simulation package.
  *
  * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
- * Copyright (c) 2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
 struct gmx_grppairener_t;
 struct t_forcerec;
 struct t_pbc;
+struct t_mdatoms;
 
 namespace gmx
 {
 class StepWorkload;
-}
+} // namespace gmx
 
 /*! \brief Calculate VdW/charge listed pair interactions (usually 1-4
  * interactions).
index 2564a77f3c659629c157b50f21a3629a37aed5a0..6299f3384a49f28dc2ab907591478a93fced8cb2 100644 (file)
@@ -577,9 +577,7 @@ protected:
     {
         SCOPED_TRACE(std::string("Testing PBC type: ") + c_pbcTypeNames[pbcType_]);
         std::vector<int>  ddgatindex = { 0, 1, 2, 3 };
-        std::vector<real> chargeA    = { 1.5, -2.0, 1.5, -1.0 };
-        t_mdatoms         mdatoms    = { 0 };
-        mdatoms.chargeA              = chargeA.data();
+        std::vector<real> charge     = { 1.5, -2.0, 1.5, -1.0 };
         /* Here we run both the standard, plain-C force+shift-forces+energy+free-energy
          * kernel flavor and the potentially optimized, with SIMD and less output,
          * force only kernels. Note that we also run the optimized kernel for free-energy
@@ -604,7 +602,7 @@ protected:
                                                 &pbc_,
                                                 lambda,
                                                 &output.dvdlambda,
-                                                &mdatoms,
+                                                charge,
                                                 /* struct t_fcdata * */ nullptr,
                                                 nullptr,
                                                 nullptr,