Use ArrayRef in do_ewald
authorejjordan <ejjordan@kth.se>
Mon, 22 Mar 2021 11:49:54 +0000 (12:49 +0100)
committerejjordan <ejjordan@kth.se>
Mon, 22 Mar 2021 11:49:58 +0000 (12:49 +0100)
Also removed a function with no implementation.

src/gromacs/ewald/ewald.cpp
src/gromacs/ewald/ewald.h
src/gromacs/mdlib/force.cpp

index 26185de74a0b59323509449e1a5a54e65f51246a..3a62f306d1f906baf7ea398209fc52d02a0357e4 100644 (file)
@@ -68,6 +68,7 @@
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/interaction_const.h"
 #include "gromacs/mdtypes/md_enums.h"
+#include "gromacs/utility/arrayref.h"
 #include "gromacs/utility/fatalerror.h"
 #include "gromacs/utility/smalloc.h"
 
@@ -97,7 +98,7 @@ static void calc_lll(const rvec box, rvec lll)
 }
 
 //! Make tables for the structure factor parts
-static void tabulateStructureFactors(int natom, const rvec x[], int kmax, cvec** eir, const rvec lll)
+static void tabulateStructureFactors(int natom, gmx::ArrayRef<const gmx::RVec> x, int kmax, cvec** eir, const rvec lll)
 {
     int i, j, m;
 
@@ -130,28 +131,27 @@ static void tabulateStructureFactors(int natom, const rvec x[], int kmax, cvec**
     }
 }
 
-real do_ewald(const t_inputrec& ir,
-              const rvec        x[],
-              rvec              f[],
-              const real        chargeA[],
-              const real        chargeB[],
-              const matrix      box,
-              const t_commrec*  cr,
-              int               natoms,
-              matrix            lrvir,
-              real              ewaldcoeff,
-              real              lambda,
-              real*             dvdlambda,
-              gmx_ewald_tab_t*  et)
+real do_ewald(const t_inputrec&              ir,
+              gmx::ArrayRef<const gmx::RVec> x,
+              gmx::ArrayRef<gmx::RVec>       f,
+              gmx::ArrayRef<const real>      chargeA,
+              gmx::ArrayRef<const real>      chargeB,
+              const matrix                   box,
+              const t_commrec*               cr,
+              int                            natoms,
+              matrix                         lrvir,
+              real                           ewaldcoeff,
+              real                           lambda,
+              real*                          dvdlambda,
+              gmx_ewald_tab_t*               et)
 {
-    real        factor = -1.0 / (4 * ewaldcoeff * ewaldcoeff);
-    const real* charge;
-    real        energy_AB[2], energy;
-    rvec        lll;
-    int         lowiy, lowiz, ix, iy, iz, n, q;
-    real        tmp, cs, ss, ak, akv, mx, my, mz, m2, scale;
-    gmx_bool    bFreeEnergy;
-    cvec**      eir;
+    real   factor = -1.0 / (4 * ewaldcoeff * ewaldcoeff);
+    real   energy_AB[2], energy;
+    rvec   lll;
+    int    lowiy, lowiz, ix, iy, iz, n, q;
+    real   tmp, cs, ss, ak, akv, mx, my, mz, m2, scale;
+    cvec** eir;
+    bool   bFreeEnergy;
 
     if (cr != nullptr)
     {
@@ -191,6 +191,7 @@ real do_ewald(const t_inputrec& ir,
     calc_lll(boxDiag, lll);
     tabulateStructureFactors(natoms, x, et->kmax, eir, lll);
 
+    gmx::ArrayRef<const real> charge;
     for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
     {
         if (!bFreeEnergy)
index 9bf4b684bb43f0a7f3be41b8fcb9c9271eed05d3..2d4d20bc98c8ded83c0253c6693eb2963d27c4b0 100644 (file)
@@ -78,6 +78,12 @@ struct t_forcerec;
 struct t_inputrec;
 struct t_complex;
 
+namespace gmx
+{
+template<typename>
+class ArrayRef;
+}
+
 struct gmx_ewald_tab_t
 {
     gmx_ewald_tab_t(const t_inputrec& ir, FILE* fp);
@@ -93,23 +99,20 @@ struct gmx_ewald_tab_t
     std::vector<t_complex> tab_qxyz;
 };
 
-/*! \brief Initialize the tables used in the Ewald long-ranged part */
-void init_ewald_tab(struct gmx_ewald_tab_t** et, const t_inputrec& ir, FILE* fp);
-
 /*! \brief Do the long-ranged part of an Ewald calculation */
-real do_ewald(const t_inputrec& ir,
-              const rvec        x[],
-              rvec              f[],
-              const real        chargeA[],
-              const real        chargeB[],
-              const matrix      box,
-              const t_commrec*  cr,
-              int               natoms,
-              matrix            lrvir,
-              real              ewaldcoeff,
-              real              lambda,
-              real*             dvdlambda,
-              gmx_ewald_tab_t*  et);
+real do_ewald(const t_inputrec&              ir,
+              gmx::ArrayRef<const gmx::RVec> x,
+              gmx::ArrayRef<gmx::RVec>       f,
+              gmx::ArrayRef<const real>      chargeA,
+              gmx::ArrayRef<const real>      chargeB,
+              const matrix                   box,
+              const t_commrec*               cr,
+              int                            natoms,
+              matrix                         lrvir,
+              real                           ewaldcoeff,
+              real                           lambda,
+              real*                          dvdlambda,
+              gmx_ewald_tab_t*               et);
 
 /*! \brief Calculate the correction to the Ewald sum, due to a net system
  * charge.
index a598288f34af6692df77bd33b2f95f1d0ff60cdb..8c964b9945262683eeebcfbf276bade5a22beca3 100644 (file)
@@ -268,12 +268,11 @@ void calculateLongRangeNonbondeds(t_forcerec*                    fr,
 
         if (fr->ic->eeltype == CoulombInteractionType::Ewald)
         {
-            const rvec* x = as_rvec_array(coordinates.data());
-            Vlr_q         = do_ewald(ir,
-                             x,
-                             as_rvec_array(forceWithVirial->force_.data()),
-                             md->chargeA,
-                             md->chargeB,
+            Vlr_q = do_ewald(ir,
+                             coordinates,
+                             forceWithVirial->force_,
+                             gmx::arrayRefFromArray(md->chargeA, md->nr),
+                             gmx::arrayRefFromArray(md->chargeB, md->nr),
                              box,
                              cr,
                              md->homenr,