Move orientation restraints checks inside init function
authorMark Abraham <mark.j.abraham@gmail.com>
Tue, 4 Feb 2014 15:29:15 +0000 (16:29 +0100)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Sun, 9 Feb 2014 16:17:35 +0000 (17:17 +0100)
This removes clutter from runner(), much like init_disre() already
does. Removed duplicate check for the number of orientation restraints
being zero. Removed an unused function parameter.

Change-Id: Icc6a08a2ed4af010e65bfe1f5e49b6c229454df1

src/gromacs/gmxlib/orires.c
src/gromacs/legacyheaders/orires.h
src/programs/mdrun/runner.c

index b354453c94cc8ee9e5a9577a224648393b940fff..4575603f995978fabc646bfec676679fe1016ac0 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
 void init_orires(FILE *fplog, const gmx_mtop_t *mtop,
                  rvec xref[],
                  const t_inputrec *ir,
-                 const gmx_multisim_t *ms, t_oriresdata *od,
-                 t_state *state)
+                 const t_commrec *cr, t_oriresdata *od,
+                 t_state *state, gmx_bool bIsParticleDecomposition)
 {
-    int                     i, j, d, ex, nmol, nr, *nr_ex;
+    int                     i, j, d, ex, nmol, *nr_ex;
     double                  mtot;
     rvec                    com;
     gmx_mtop_ilistloop_t    iloop;
     t_ilist                *il;
     gmx_mtop_atomloop_all_t aloop;
     t_atom                 *atom;
+    const gmx_multisim_t   *ms;
+
+    od->nr = gmx_mtop_ftype_count(mtop, F_ORIRES);
+    if (0 == od->nr)
+    {
+        /* Not doing orientation restraints */
+        return;
+    }
+
+    if (PAR(cr) && !bIsParticleDecomposition)
+    {
+        gmx_fatal(FARGS, "Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
+    }
+    /* Orientation restraints */
+    if (!MASTER(cr))
+    {
+        /* Nothing to do */
+        return;
+    }
+    ms = cr->ms;
 
     od->fc  = ir->orires_fc;
     od->nex = 0;
     od->S   = NULL;
-
     od->M   = NULL;
     od->eig = NULL;
     od->v   = NULL;
 
-    od->nr = gmx_mtop_ftype_count(mtop, F_ORIRES);
-    if (od->nr == 0)
-    {
-        return;
-    }
-
     nr_ex = NULL;
 
     iloop = gmx_mtop_ilistloop_init(mtop);
index 58fc412c411b9fd53955027f3c2ceddd507ba8e5..3d6c2431292802255bfac8d6fd21deba6cfe856b 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2010, by the GROMACS development team, led by
+ * Copyright (c) 2010,2014, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -48,9 +48,12 @@ extern "C" {
 void init_orires(FILE *fplog, const gmx_mtop_t *mtop,
                  rvec x[],
                  const t_inputrec *ir,
-                 const gmx_multisim_t *ms, t_oriresdata *od,
-                 t_state *state);
-/* Initializes all the orientation restraint stuff in *od */
+                 const t_commrec *cr, t_oriresdata *od,
+                 t_state *state,
+                 gmx_bool bIsParticleDecomposition);
+/* Decides whether orientation restraints can work, and initializes
+   all the orientation restraint stuff in *od (and assumes *od is
+   already allocated. */
 
 real calc_orires_dev(const gmx_multisim_t *ms,
                      int nfa, const t_iatom fa[], const t_iparams ip[],
index 2a29f9cda305e551a52ad289fbe14678ea13285e..bdebcd9b139c62a83b2bfd3435b2b6ef970e0dd2 100644 (file)
@@ -1358,19 +1358,8 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
     /* This needs to be called before read_checkpoint to extend the state */
     init_disres(fplog, mtop, inputrec, cr, Flags & MD_PARTDEC, fcd, state, repl_ex_nst > 0);
 
-    if (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0)
-    {
-        if (PAR(cr) && !(Flags & MD_PARTDEC))
-        {
-            gmx_fatal(FARGS, "Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
-        }
-        /* Orientation restraints */
-        if (MASTER(cr))
-        {
-            init_orires(fplog, mtop, state->x, inputrec, cr->ms, &(fcd->orires),
-                        state);
-        }
-    }
+    init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
+                state, Flags & MD_PARTDEC);
 
     if (DEFORM(*inputrec))
     {