minor function argument cleanup
authorSzilard Pall <pszilard@cbr.su.se>
Mon, 27 May 2013 17:37:18 +0000 (19:37 +0200)
committerSzilard Pall <pszilard@cbr.su.se>
Mon, 27 May 2013 18:22:46 +0000 (20:22 +0200)
Removed an unused function argument and changed anoehter one to not pass
the entire nonbonded verlet data structure to the nbnxn_cuda module, but
only what's needed from it.

Change-Id: I89699548f69245eba6d8451683ba6dc36af75dc2

include/nbnxn_cuda_data_mgmt.h
src/kernel/runner.c
src/mdlib/forcerec.c
src/mdlib/nbnxn_cuda/nbnxn_cuda_data_mgmt.cu

index 0d8ce3d0c82606a3cc988e40f34c1e72bbdd543a..e542ed978e6ae3ba5898cb2000e1d2de8c34ea57 100644 (file)
@@ -66,9 +66,9 @@ void nbnxn_cuda_init(FILE *fplog,
 
 /*! Initializes simulation constant data. */
 FUNC_QUALIFIER
-void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t           p_cu_nb,
-                           const interaction_const_t *ic,
-                           const nonbonded_verlet_t  *nbv) FUNC_TERM
+void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t                cu_nb,
+                           const interaction_const_t      *ic,
+                           const nonbonded_verlet_group_t *nbv_group) FUNC_TERM
 
 /*! Initializes pair-list data for GPU, called at every pair search step. */
 FUNC_QUALIFIER
index d9bfae0482ed39702cf10ea0d99d8af24414da7f..ff98cfba63f8effde389e3452767a7313b4b7d96 100644 (file)
@@ -630,7 +630,6 @@ static void increase_nstlist(FILE *fp, t_commrec *cr,
 static void prepare_verlet_scheme(FILE             *fplog,
                                   gmx_hw_info_t    *hwinfo,
                                   t_commrec        *cr,
-                                  gmx_hw_opt_t     *hw_opt,
                                   const char       *nbpu_opt,
                                   t_inputrec       *ir,
                                   const gmx_mtop_t *mtop,
@@ -987,7 +986,7 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
 
         if (inputrec->cutoff_scheme == ecutsVERLET)
         {
-            prepare_verlet_scheme(fplog, hwinfo, cr, hw_opt, nbpu_opt,
+            prepare_verlet_scheme(fplog, hwinfo, cr, nbpu_opt,
                                   inputrec, mtop, state->box,
                                   &minf.bUseGPU);
         }
index 2253eb105106bf3568d65dbccd8c4c9ea7a14cb8..2c2e4da07bd7f4c4cd2bd345feb436e913d6ec5c 100644 (file)
@@ -1863,7 +1863,7 @@ void init_interaction_const(FILE                 *fp,
 
     if (fr->nbv != NULL && fr->nbv->bUseGPU)
     {
-        nbnxn_cuda_init_const(fr->nbv->cu_nbv, ic, fr->nbv);
+        nbnxn_cuda_init_const(fr->nbv->cu_nbv, ic, fr->nbv->grp);
     }
 
     bUsesSimpleTables = uses_simple_tables(fr->cutoff_scheme, fr->nbv, -1);
index f2d427b18cbaf0628c77049fe46768774c81d71c..a7ef47ce2e301396c1ed2086b9fdce22c12d7df6 100644 (file)
@@ -237,13 +237,13 @@ static int pick_ewald_kernel_type(bool                   bTwinCut,
 /*! Initializes the nonbonded parameter data structure. */
 static void init_nbparam(cu_nbparam_t *nbp,
                          const interaction_const_t *ic,
-                         const nonbonded_verlet_t *nbv,
+                         const nbnxn_atomdata_t *nbat,
                          const cuda_dev_info_t *dev_info)
 {
     cudaError_t stat;
     int         ntypes, nnbfp;
 
-    ntypes  = nbv->grp[0].nbat->ntype;
+    ntypes  = nbat->ntype;
 
     nbp->ewald_beta = ic->ewaldcoeff;
     nbp->sh_ewald   = ic->sh_ewald;
@@ -284,7 +284,7 @@ static void init_nbparam(cu_nbparam_t *nbp,
     nnbfp = 2*ntypes*ntypes;
     stat = cudaMalloc((void **)&nbp->nbfp, nnbfp*sizeof(*nbp->nbfp));
     CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp");
-    cu_copy_H2D(nbp->nbfp, nbv->grp[0].nbat->nbfp, nnbfp*sizeof(*nbp->nbfp));
+    cu_copy_H2D(nbp->nbfp, nbat->nbfp, nnbfp*sizeof(*nbp->nbfp));
 
     cudaChannelFormatDesc cd   = cudaCreateChannelDesc<float>();
     stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_texref(),
@@ -642,12 +642,12 @@ void nbnxn_cuda_init(FILE *fplog,
     }
 }
 
-void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t cu_nb,
-                           const interaction_const_t *ic,
-                           const nonbonded_verlet_t *nbv)
+void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t                cu_nb,
+                           const interaction_const_t      *ic,
+                           const nonbonded_verlet_group_t *nbv_group)
 {
-    init_atomdata_first(cu_nb->atdat, nbv->grp[0].nbat->ntype);
-    init_nbparam(cu_nb->nbparam, ic, nbv, cu_nb->dev_info);
+    init_atomdata_first(cu_nb->atdat, nbv_group[0].nbat->ntype);
+    init_nbparam(cu_nb->nbparam, ic, nbv_group[0].nbat, cu_nb->dev_info);
 
     /* clear energy and shift force outputs */
     nbnxn_cuda_clear_e_fshift(cu_nb);