flops table now reports CUDA analytical Ewald
authorBerk Hess <hess@kth.se>
Thu, 23 May 2013 13:21:21 +0000 (15:21 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Thu, 20 Jun 2013 14:20:56 +0000 (16:20 +0200)
Change-Id: I72dfdc3b3f1ecb651897922c2e7ffec784c3bf60

include/nbnxn_cuda_data_mgmt.h
src/mdlib/nbnxn_cuda/nbnxn_cuda_data_mgmt.cu
src/mdlib/sim_util.c

index e542ed978e6ae3ba5898cb2000e1d2de8c34ea57..1d29225cbaffa92fb954d7c12f4005ec043620dc 100644 (file)
@@ -128,6 +128,17 @@ int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
 }
 #endif
 
+/*! Returns if analytical Ewald CUDA kernels are used. */
+FUNC_QUALIFIER
+gmx_bool nbnxn_cuda_is_kernel_ewald_analytical(const nbnxn_cuda_ptr_t cu_nb)
+#ifdef GMX_GPU
+;
+#else
+{
+    return FALSE;
+}
+#endif
+
 #ifdef __cplusplus
 }
 #endif
index 882063167495fbd1775ad6658c0a52577df759d5..f42fa1e688107501cb0604509775bd28b1a38f7a 100644 (file)
@@ -964,3 +964,9 @@ int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
         gpu_min_ci_balanced_factor*cu_nb->dev_info->prop.multiProcessorCount : 0;
 
 }
+
+gmx_bool nbnxn_cuda_is_kernel_ewald_analytical(const nbnxn_cuda_ptr_t cu_nb)
+{
+    return ((cu_nb->nbparam->eeltype == eelCuEWALD_ANA) ||
+            (cu_nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
+}
index e0797cdc0026302c5b77b8c3e11d4e6ece1c73f3..f964ce4799ed00bb3a4d791fa945e9ee24165ba6 100644 (file)
@@ -618,6 +618,7 @@ static void do_nb_verlet(t_forcerec *fr,
     int                        nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
     char                      *env;
     nonbonded_verlet_group_t  *nbvg;
+    gmx_bool                  bCUDA;
 
     if (!(flags & GMX_FORCE_NONBONDED))
     {
@@ -633,7 +634,9 @@ static void do_nb_verlet(t_forcerec *fr,
         gmx_incons("Invalid cut-off scheme passed!");
     }
 
-    if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
+    bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
+
+    if (!bCUDA)
     {
         wallcycle_sub_start(wcycle, ewcsNONBONDED);
     }
@@ -701,7 +704,7 @@ static void do_nb_verlet(t_forcerec *fr,
             gmx_incons("Invalid nonbonded kernel type passed!");
 
     }
-    if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
+    if (!bCUDA)
     {
         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
     }
@@ -710,13 +713,14 @@ static void do_nb_verlet(t_forcerec *fr,
     {
         enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
     }
-    else if (nbvg->ewald_excl == ewaldexclTable)
+    else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
+             (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
     {
-        enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
+        enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
     }
     else
     {
-        enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
+        enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
     }
     enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
     if (flags & GMX_FORCE_ENERGY)