Remove bDomDec from t_forcerec
authorMark Abraham <mark.j.abraham@gmail.com>
Tue, 22 May 2018 19:42:35 +0000 (21:42 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Thu, 24 May 2018 19:23:26 +0000 (21:23 +0200)
Whether DD is in use should be conveyed in a consistent manner, which
for now is best done with cr.

Change-Id: I1e9e9763b5beee2ad1540c404f57c91597abc919

src/gromacs/mdlib/force.cpp
src/gromacs/mdlib/force.h
src/gromacs/mdlib/forcerec.cpp
src/gromacs/mdlib/rf_util.cpp
src/gromacs/mdlib/sim_util.cpp
src/gromacs/mdtypes/forcerec.h

index 416c95c3cd5225fdc943b512c964f363d733ccc6..205b065f5312d438757993257098e12763b7a7e9 100644 (file)
@@ -549,7 +549,8 @@ void do_force_lowlevel(t_forcerec           *fr,
         {
             real dvdl_rf_excl      = 0;
             enerd->term[F_RF_EXCL] =
-                RF_excl_correction(fr, graph, md, excl, x, forceForUseWithShiftForces,
+                RF_excl_correction(fr, graph, md, excl, DOMAINDECOMP(cr),
+                                   x, forceForUseWithShiftForces,
                                    fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
 
             enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
index 34075da9ce248c390fbc195a53a3483b61348db3..793eab2b85d3b4bee5e80c2eedf29f437a91e231 100644 (file)
@@ -83,6 +83,7 @@ real RF_excl_correction(const t_forcerec *fr,
                         const t_graph    *g,
                         const t_mdatoms  *mdatoms,
                         const t_blocka   *excl,
+                        bool              usingDomainDecomposition,
                         rvec              x[],
                         rvec              f[],
                         rvec             *fshift,
index 480c012ca9587568b113c591d118d83f94dd15fe..e6a4c8e73aa90b897a35e26db13ca7f771004f52 100644 (file)
@@ -2337,8 +2337,6 @@ void init_forcerec(FILE                             *fp,
     /* By default we turn SIMD kernels on, but it might be turned off further down... */
     fr->use_simd_kernels = TRUE;
 
-    fr->bDomDec = DOMAINDECOMP(cr);
-
     if (check_box(ir->ePBC, box))
     {
         gmx_fatal(FARGS, check_box(ir->ePBC, box));
index fc554d6b2cbce26402b64b5ee01b90ceeefb4d95..dd99949e6bf86b3cd4919f0fddf632175993d962 100644 (file)
@@ -56,6 +56,7 @@ real RF_excl_correction(const t_forcerec *fr,
                         const t_graph    *g,
                         const t_mdatoms  *mdatoms,
                         const t_blocka   *excl,
+                        bool              usingDomainDecomposition,
                         rvec              x[],
                         rvec              f[],
                         rvec             *fshift,
@@ -92,7 +93,7 @@ real RF_excl_correction(const t_forcerec *fr,
     AA      = excl->a;
     ki      = CENTRAL;
 
-    if (fr->bDomDec)
+    if (usingDomainDecomposition)
     {
         niat = excl->nr;
     }
index 924cd4ec18e4647feb42ff179c08bf7d8ce1a192..e3cff865845d1b60f785fc4f4b8d2ca2227e0ee4 100644 (file)
@@ -1197,7 +1197,7 @@ static void do_force_cutsVERLET(FILE *fplog,
         box_diag[ZZ] = box[ZZ][ZZ];
 
         wallcycle_start(wcycle, ewcNS);
-        if (!fr->bDomDec)
+        if (!DOMAINDECOMP(cr))
         {
             wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
             nbnxn_put_on_grid(nbv->nbs.get(), fr->ePBC, box,
index b140dee8e07fdcb8bf89e1f157557ce47005aa24..52878e380dda62928ea77eb17caf96d72df66697 100644 (file)
@@ -157,11 +157,9 @@ struct ewald_corr_thread_t;
 struct t_forcerec {
     struct interaction_const_t *ic;
 
-    /* Domain Decomposition */
-    gmx_bool bDomDec;
-
     /* PBC stuff */
     int                         ePBC;
+    //! Whether PBC must be considered for bonded interactions.
     gmx_bool                    bMolPBC;
     int                         rc_scaling;
     rvec                        posres_com;