Disable AllvsAll kernels for very small systems (< 64 atoms).
authorRossen Apostolov <rossen@cbr.su.se>
Fri, 2 Jul 2010 11:13:48 +0000 (13:13 +0200)
committerRossen Apostolov <rossen@cbr.su.se>
Fri, 2 Jul 2010 11:13:48 +0000 (13:13 +0200)
This is only a temporary fix, and bug 416 should stay open until resolved.

include/force.h
src/kernel/runner.c
src/mdlib/forcerec.c
src/tools/gmx_membed.c

index df58728da64acc0384d107011f5a5df8164f9ce0..ce484a0dd45aaf7efedb9c98964b4f13b97824e2 100644 (file)
@@ -126,7 +126,7 @@ forcerec_set_ranges(t_forcerec *fr,
                    int natoms_force,int natoms_f_novirsum);
 /* Set the number of cg's and atoms for the force calculation */
 
-extern bool can_use_allvsall(const t_inputrec *ir,
+extern bool can_use_allvsall(const t_inputrec *ir, const gmx_mtop_t *mtop,
                              bool bPrintNote,t_commrec *cr,FILE *fp);
 /* Returns if we can use all-vs-all loops.
  * If bPrintNote==TRUE, prints a note, if necessary, to stderr
index 4f57fcd3f1442e903cd915059dcc007dc17af84d..2f763b890ca77f20b5497fc20984acaa279887eb 100644 (file)
@@ -421,7 +421,7 @@ int mdrunner(int nthreads_requested, FILE *fplog,t_commrec *cr,int nfile,
             );
     }
 
-    if (can_use_allvsall(inputrec,TRUE,cr,fplog))
+    if (can_use_allvsall(inputrec,mtop,TRUE,cr,fplog))
     {
         /* All-vs-all loops do not work with domain decomposition */
         Flags |= MD_PARTDEC;
index 50b1396f8627d897e0d79a59bdf500159abfddc5..0872800b5ceb78d05f79044a5bd9a5e98bf55a77 100644 (file)
@@ -1162,7 +1162,7 @@ static real cutoff_inf(real cutoff)
     return cutoff;
 }
 
-bool can_use_allvsall(const t_inputrec *ir,
+bool can_use_allvsall(const t_inputrec *ir, const gmx_mtop_t *mtop,
                       bool bPrintNote,t_commrec *cr,FILE *fp)
 {
     bool bAllvsAll;
@@ -1172,7 +1172,10 @@ bool can_use_allvsall(const t_inputrec *ir,
     bAllvsAll = FALSE;
 #else
     bAllvsAll =
-        (ir->rlist==0            &&
+        (
+         /* disable for very small systems (bug 416) */
+         mtop->natoms > 64       &&
+         ir->rlist==0            &&
          ir->rcoulomb==0         &&
          ir->rvdw==0             &&
          ir->ePBC==epbcNONE      &&
@@ -1276,7 +1279,7 @@ void init_forcerec(FILE *fp,
     fr->sc_sigma6  = pow(ir->sc_sigma,6);
     
     /* Check if we can/should do all-vs-all kernels */
-    fr->bAllvsAll       = can_use_allvsall(ir,FALSE,NULL,NULL);
+    fr->bAllvsAll       = can_use_allvsall(ir,mtop,FALSE,NULL,NULL);
     fr->AllvsAll_work   = NULL;
     fr->AllvsAll_workgb = NULL;
 
index c8c4a735431b3870b6ede8dc88233b787b69daf8..1bf8710c5c85a816f258b29e482f48c8aee049fa 100644 (file)
@@ -3635,7 +3635,7 @@ int mdrunner_membed(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
         init_parallel(fplog, cr, inputrec, mtop, state);
     }
 
-    if (can_use_allvsall(inputrec,TRUE,cr,fplog))
+    if (can_use_allvsall(inputrec,mtop,TRUE,cr,fplog))
     {
         /* All-vs-all loops do not work with domain decomposition */
         Flags |= MD_PARTDEC;