Added message to warn that group kernels are deprecated.
authorErik Lindahl <erik@kth.se>
Mon, 24 Feb 2014 23:06:24 +0000 (00:06 +0100)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Wed, 26 Feb 2014 15:30:40 +0000 (16:30 +0100)
Nothing will change for now, but both grompp and mdrun now
warn that we will remove the support for group kernels in
a future version when all interaction forms are supported by
the verlet cutoff scheme (tentatively a year from now).

Change-Id: Id00a67e44e60818cb0e6f1897d6939f53ea01540

src/gromacs/gmxpreprocess/readir.c
src/gromacs/mdlib/forcerec.c

index 63b5e8f4b4190c3a082a65cf831ed742b25d5b10..b63b5fc843c8a3d3603eee0514a3d94d56690811 100644 (file)
@@ -272,6 +272,11 @@ void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts,
 
     if (ir->cutoff_scheme == ecutsGROUP)
     {
+        warning_note(wi,
+                     "The group cutoff scheme is deprecated in Gromacs 5.0 and will be removed in a future "
+                     "release when all interaction forms are supported for the verlet scheme. The verlet "
+                     "scheme already scales better, and it is compatible with GPUs and other accelerators.");
+
         /* BASIC CUT-OFF STUFF */
         if (ir->rlist == 0 ||
             !((ir_coulomb_might_be_zero_at_cutoff(ir) && ir->rcoulomb > ir->rlist) ||
index e24126491b8cc0b5ac566049b155ead80ef40334..0c7db6ed40fbdadb4eb34f350c6366292927c389 100644 (file)
@@ -2472,6 +2472,21 @@ void init_forcerec(FILE              *fp,
     fr->bGrid         = (ir->ns_type == ensGRID);
     fr->ePBC          = ir->ePBC;
 
+    if (fr->cutoff_scheme == ecutsGROUP)
+    {
+        const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
+            "removed in a future release when 'verlet' supports all interaction forms.\n";
+
+        if (MASTER(cr))
+        {
+            fprintf(stderr, "\n%s\n", note);
+        }
+        if (fp != NULL)
+        {
+            fprintf(fp, "\n%s\n", note);
+        }
+    }
+
     /* Determine if we will do PBC for distances in bonded interactions */
     if (fr->ePBC == epbcNONE)
     {