Remove t_mdatoms from Nbnxm setAtomProperties()
authorBerk Hess <hess@kth.se>
Wed, 20 May 2020 14:56:32 +0000 (16:56 +0200)
committerBerk Hess <hess@kth.se>
Wed, 20 May 2020 15:22:18 +0000 (17:22 +0200)
src/gromacs/mdlib/sim_util.cpp
src/gromacs/mdrun/tpi.cpp
src/gromacs/nbnxm/atomdata.cpp
src/gromacs/nbnxm/atomdata.h
src/gromacs/nbnxm/benchmark/bench_setup.cpp
src/gromacs/nbnxm/nbnxm.cpp
src/gromacs/nbnxm/nbnxm.h

index 31b0955b5b421dc503cefe7fba855ea0723d453f..ca42600eda8170b33f94e5dd331d1c0e06d1fe9d 100644 (file)
@@ -1213,7 +1213,8 @@ void do_force(FILE*                               fplog,
             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
         }
 
-        nbv->setAtomProperties(*mdatoms, fr->cginfo);
+        nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
+                               gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr), fr->cginfo);
 
         wallcycle_stop(wcycle, ewcNS);
 
index ac2340ce17ebef9a3764f1d50ece9a97cde7ae4f..280d1f22cb2b352ee0308d38b5d25eacc93f76b0 100644 (file)
@@ -677,7 +677,9 @@ void LegacySimulator::do_tpi()
                                   -1, fr->cginfo, x, 0, nullptr);
 
                 /* TODO: Avoid updating all atoms at every bNS step */
-                fr->nbv->setAtomProperties(*mdatoms, fr->cginfo);
+                fr->nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
+                                           gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
+                                           fr->cginfo);
 
                 fr->nbv->constructPairlist(InteractionLocality::Local, top.excls, step, nrnb);
 
index a64fe1a3812d1037afb4bc38464e55738fdb19c8..35901f40f122e766b13747d0d3879589e5e98531 100644 (file)
@@ -52,7 +52,6 @@
 #include "gromacs/math/vec.h"
 #include "gromacs/mdlib/gmx_omp_nthreads.h"
 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
-#include "gromacs/mdtypes/mdatom.h"
 #include "gromacs/nbnxm/nbnxm.h"
 #include "gromacs/pbcutil/ishift.h"
 #include "gromacs/simd/simd.h"
@@ -724,7 +723,7 @@ static void copy_lj_to_nbat_lj_comb(gmx::ArrayRef<const real> ljparam_type, cons
 /* Sets the atom type in nbnxn_atomdata_t */
 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t::Params* params,
                                          const Nbnxm::GridSet&     gridSet,
-                                         const int*                type)
+                                         ArrayRef<const int>       atomTypes)
 {
     params->type.resize(gridSet.numGridAtomsTotal());
 
@@ -736,8 +735,9 @@ static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t::Params* params,
             const int numAtoms   = grid.paddedNumAtomsInColumn(i);
             const int atomOffset = grid.firstAtomInColumn(i);
 
-            copy_int_to_nbat_int(gridSet.atomIndices().data() + atomOffset, grid.numAtomsInColumn(i),
-                                 numAtoms, type, params->numTypes - 1, params->type.data() + atomOffset);
+            copy_int_to_nbat_int(gridSet.atomIndices().data() + atomOffset,
+                                 grid.numAtomsInColumn(i), numAtoms, atomTypes.data(),
+                                 params->numTypes - 1, params->type.data() + atomOffset);
         }
     }
 }
@@ -782,7 +782,9 @@ static void nbnxn_atomdata_set_ljcombparams(nbnxn_atomdata_t::Params* params,
 }
 
 /* Sets the charges in nbnxn_atomdata_t *nbat */
-static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t* nbat, const Nbnxm::GridSet& gridSet, const real* charge)
+static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t*     nbat,
+                                       const Nbnxm::GridSet& gridSet,
+                                       ArrayRef<const real>  charges)
 {
     if (nbat->XFormat != nbatXYZQ)
     {
@@ -804,7 +806,7 @@ static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t* nbat, const Nbnxm::Grid
                 int   i;
                 for (i = 0; i < numAtoms; i++)
                 {
-                    *q = charge[gridSet.atomIndices()[atomOffset + i]];
+                    *q = charges[gridSet.atomIndices()[atomOffset + i]];
                     q += STRIDE_XYZQ;
                 }
                 /* Complete the partially filled last cell with zeros */
@@ -820,7 +822,7 @@ static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t* nbat, const Nbnxm::Grid
                 int   i;
                 for (i = 0; i < numAtoms; i++)
                 {
-                    *q = charge[gridSet.atomIndices()[atomOffset + i]];
+                    *q = charges[gridSet.atomIndices()[atomOffset + i]];
                     q++;
                 }
                 /* Complete the partially filled last cell with zeros */
@@ -926,7 +928,7 @@ copy_egp_to_nbat_egps(const int* a, int na, int na_round, int na_c, int bit_shif
 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t::Params* params,
                                             const Nbnxm::GridSet&     gridSet,
-                                            const int*                atinfo)
+                                            ArrayRef<const int>       atomInfo)
 {
     if (params->nenergrp == 1)
     {
@@ -944,7 +946,7 @@ static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t::Params* params,
             const int atomOffset = grid.firstAtomInColumn(i);
 
             copy_egp_to_nbat_egps(gridSet.atomIndices().data() + atomOffset, grid.numAtomsInColumn(i),
-                                  numAtoms, c_nbnxnCpuIClusterSize, params->neg_2log, atinfo,
+                                  numAtoms, c_nbnxnCpuIClusterSize, params->neg_2log, atomInfo.data(),
                                   params->energrp.data() + grid.atomToCluster(atomOffset));
         }
     }
@@ -953,14 +955,15 @@ static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t::Params* params,
 /* Sets all required atom parameter data in nbnxn_atomdata_t */
 void nbnxn_atomdata_set(nbnxn_atomdata_t*     nbat,
                         const Nbnxm::GridSet& gridSet,
-                        const t_mdatoms*      mdatoms,
-                        const int*            atinfo)
+                        ArrayRef<const int>   atomTypes,
+                        ArrayRef<const real>  atomCharges,
+                        ArrayRef<const int>   atomInfo)
 {
     nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
 
-    nbnxn_atomdata_set_atomtypes(&params, gridSet, mdatoms->typeA);
+    nbnxn_atomdata_set_atomtypes(&params, gridSet, atomTypes);
 
-    nbnxn_atomdata_set_charges(nbat, gridSet, mdatoms->chargeA);
+    nbnxn_atomdata_set_charges(nbat, gridSet, atomCharges);
 
     if (gridSet.haveFep())
     {
@@ -970,7 +973,7 @@ void nbnxn_atomdata_set(nbnxn_atomdata_t*     nbat,
     /* This must be done after masking types for FEP */
     nbnxn_atomdata_set_ljcombparams(&params, nbat->XFormat, gridSet);
 
-    nbnxn_atomdata_set_energygroups(&params, gridSet, atinfo);
+    nbnxn_atomdata_set_energygroups(&params, gridSet, atomInfo);
 }
 
 /* Copies the shift vector array to nbnxn_atomdata_t */
index 22d115706e0a394f3f4f55e95676ff07675548bc..4b3a398562b87f3fd1c485d2c99326c5f48d761e 100644 (file)
@@ -64,7 +64,6 @@ class MDLogger;
 struct NbnxmGpu;
 struct nbnxn_atomdata_t;
 struct nonbonded_verlet_t;
-struct t_mdatoms;
 struct tMPI_Atomic;
 
 class GpuEventSynchronizer;
@@ -340,10 +339,11 @@ void nbnxn_atomdata_init(const gmx::MDLogger&      mdlog,
                          int                       nout);
 
 //! Sets the atomdata after pair search
-void nbnxn_atomdata_set(nbnxn_atomdata_t*     nbat,
-                        const Nbnxm::GridSet& gridSet,
-                        const t_mdatoms*      mdatoms,
-                        const int*            atinfo);
+void nbnxn_atomdata_set(nbnxn_atomdata_t*         nbat,
+                        const Nbnxm::GridSet&     gridSet,
+                        gmx::ArrayRef<const int>  atomTypes,
+                        gmx::ArrayRef<const real> atomCharges,
+                        gmx::ArrayRef<const int>  atomInfo);
 
 //! Copy the shift vectors to nbat
 void nbnxn_atomdata_copy_shiftvec(gmx_bool dynamic_box, rvec* shift_vec, nbnxn_atomdata_t* nbat);
index b933cee1b25afaf4902abf31b0458a75c33e93ab..d5d4ab9d22a217a3ba04efb7080154568ab72a53 100644 (file)
@@ -227,11 +227,7 @@ static std::unique_ptr<nonbonded_verlet_t> setupNbnxmForBenchInstance(const Kern
 
     nbv->constructPairlist(gmx::InteractionLocality::Local, system.excls, 0, &nrnb);
 
-    t_mdatoms mdatoms;
-    // We only use (read) the atom type and charge from mdatoms
-    mdatoms.typeA   = const_cast<int*>(system.atomTypes.data());
-    mdatoms.chargeA = const_cast<real*>(system.charges.data());
-    nbv->setAtomProperties(mdatoms, atomInfo);
+    nbv->setAtomProperties(system.atomTypes, system.charges, atomInfo);
 
     return nbv;
 }
index ab79a8cf5b8482d453348d108e48e59fa07a77e5..20495969e0a9aa42ee464e40df15d0822d504ff5 100644 (file)
@@ -119,9 +119,11 @@ void nonbonded_verlet_t::setLocalAtomOrder()
     pairSearch_->setLocalAtomOrder();
 }
 
-void nonbonded_verlet_t::setAtomProperties(const t_mdatoms& mdatoms, gmx::ArrayRef<const int> atomInfo)
+void nonbonded_verlet_t::setAtomProperties(gmx::ArrayRef<const int>  atomTypes,
+                                           gmx::ArrayRef<const real> atomCharges,
+                                           gmx::ArrayRef<const int>  atomInfo)
 {
-    nbnxn_atomdata_set(nbat.get(), pairSearch_->gridSet(), &mdatoms, atomInfo.data());
+    nbnxn_atomdata_set(nbat.get(), pairSearch_->gridSet(), atomTypes, atomCharges, atomInfo);
 }
 
 void nonbonded_verlet_t::convertCoordinates(const gmx::AtomLocality        locality,
index d5ced753f7cc2b2f2d8dac8d3e4508a486b39d28..f1d7bed568c33e9481702e9ea8415f4a8240df0d 100644 (file)
@@ -266,7 +266,9 @@ public:
                            t_nrnb*                      nrnb);
 
     //! Updates all the atom properties in Nbnxm
-    void setAtomProperties(const t_mdatoms& mdatoms, gmx::ArrayRef<const int> atomInfo);
+    void setAtomProperties(gmx::ArrayRef<const int>  atomTypes,
+                           gmx::ArrayRef<const real> atomCharges,
+                           gmx::ArrayRef<const int>  atomInfo);
 
     /*!\brief Convert the coordinates to NBNXM format for the given locality.
      *