Pre-emptively fix some C++ warnings.
authorTeemu Murtola <teemu.murtola@gmail.com>
Tue, 2 Jul 2013 16:51:38 +0000 (19:51 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Fri, 19 Jul 2013 20:07:36 +0000 (22:07 +0200)
Moved sepdvdlformat from force.h into a source file to avoid warnings
about an unused variable.  These appeared as soon as a C++ source file
included force.h.  As an additional benefit, makes printing this
information more typesafe; this is how #1294 was found.

Change-Id: Id856fbe193d9b4a74ecfb49cba9eab03bad0a7ba

src/gromacs/legacyheaders/force.h
src/gromacs/mdlib/force.c
src/gromacs/mdlib/minimize.c
src/gromacs/mdlib/sim_util.c
src/programs/mdrun/md.c

index c31e93f04735b67985f4b5a5a08bea32c757a2cf..f456fa0fd65ef5b7f78e6d05e77961dbc606af69 100644 (file)
@@ -50,7 +50,7 @@
 extern "C" {
 #endif
 
-static const char *sepdvdlformat = "  %-30s V %12.5e  dVdl %12.5e\n";
+void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda);
 
 void calc_vir(FILE *fplog, int nxf, rvec x[], rvec f[], tensor vir,
               gmx_bool bScrewPBC, matrix box);
index 4f69c437b3404c8407f150ad0e90b537cdfe04da..6b5589e719eae0ad75666142a78f51f031107d64 100644 (file)
@@ -139,6 +139,11 @@ static void reduce_thread_forces(int n, rvec *f,
     }
 }
 
+void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda)
+{
+    fprintf(fplog, "  %-30s V %12.5e  dVdl %12.5e\n", s, v, dvdlambda);
+}
+
 void do_force_lowlevel(FILE       *fplog,   gmx_large_int_t step,
                        t_forcerec *fr,      t_inputrec *ir,
                        t_idef     *idef,    t_commrec  *cr,
@@ -182,8 +187,7 @@ void do_force_lowlevel(FILE       *fplog,   gmx_large_int_t step,
     double  t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
 #endif
 
-#define PRINT_SEPDVDL(s, v, dvdlambda) if (bSepDVDL) {fprintf(fplog, sepdvdlformat, s, v, dvdlambda); }
-
+#define PRINT_SEPDVDL(s, v, dvdlambda) if (bSepDVDL) { gmx_print_sepdvdl(fplog, s, v, dvdlambda); }
 
     set_pbc(&pbc, fr->ePBC, box);
 
index 10e4ae43fc38261baeb3068805bc2b008dbc2e42..82cd24dd1a0c21407518acd65bb938fa1efffa15 100644 (file)
@@ -796,7 +796,7 @@ static void evaluate_energy(FILE *fplog, gmx_bool bVerbose, t_commrec *cr,
                   NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
         if (fr->bSepDVDL && fplog)
         {
-            fprintf(fplog, sepdvdlformat, "Constraints", t, dvdl_constr);
+            gmx_print_sepdvdl(fplog, "Constraints", t, dvdl_constr);
         }
         enerd->term[F_DVDL_CONSTR] += dvdl_constr;
         m_add(force_vir, shake_vir, vir);
index c0fd4cd0078aeb570273dd2269eb679a6b61ac6c..5723a848db68dc970ef8068ac8afa0285abbd7d8 100644 (file)
@@ -426,8 +426,7 @@ static void posres_wrapper(FILE *fplog,
                   fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
     if (bSepDVDL)
     {
-        fprintf(fplog, sepdvdlformat,
-                interaction_function[F_POSRES].longname, v, dvdl);
+        gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
     }
     enerd->term[F_POSRES] += v;
     /* If just the force constant changes, the FEP term is linear,
@@ -480,7 +479,7 @@ static void pull_potential_wrapper(FILE *fplog,
                        cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
     if (bSepDVDL)
     {
-        fprintf(fplog, sepdvdlformat, "Com pull", enerd->term[F_COM_PULL], dvdl);
+        gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
     }
     enerd->dvdl_lin[efptRESTRAINT] += dvdl;
 }
@@ -507,7 +506,7 @@ static void pme_receive_force_ener(FILE           *fplog,
                       &cycles_seppme);
     if (bSepDVDL)
     {
-        fprintf(fplog, sepdvdlformat, "PME mesh", e, dvdl);
+        gmx_print_sepdvdl(fplog, "PME mesh", e, dvdl);
     }
     enerd->term[F_COUL_RECIP] += e;
     enerd->dvdl_lin[efptCOUL] += dvdl;
@@ -1694,12 +1693,6 @@ void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
            groups, &(inputrec->opts), top, mdatoms,
            cr, nrnb, lambda, dvdlambda, &enerd->grpp, bFillGrid,
            bDoLongRangeNS);
-        if (bSepDVDL)
-        {
-            fprintf(fplog, sepdvdlformat, "LR non-bonded", 0.0, dvdlambda);
-        }
-        enerd->dvdl_lin[efptVDW]  += dvdlambda[efptVDW];
-        enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
 
         wallcycle_stop(wcycle, ewcNS);
     }
@@ -2384,8 +2377,7 @@ void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
 
         if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
         {
-            fprintf(fplog, sepdvdlformat, "Dispersion correction",
-                    *enercorr, dvdlambda);
+            gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
         }
         if (fr->efep != efepNO)
         {
index c75b17e76cca036936ed938ec2823d7e9ea69897..503375ec066b121fd124fe43d45cd56078f6d044 100644 (file)
@@ -1780,7 +1780,7 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
 
                 if (fr->bSepDVDL && fplog && do_log)
                 {
-                    fprintf(fplog, sepdvdlformat, "Constraint dV/dl", 0.0, dvdl_constr);
+                    gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
                 }
                 if (bVV)
                 {