Remove inputrec from initialize_lambdas signature
authorJoe Jordan <ejjordan12@gmail.com>
Mon, 26 Apr 2021 15:33:40 +0000 (15:33 +0000)
committerJoe Jordan <ejjordan12@gmail.com>
Mon, 26 Apr 2021 15:33:40 +0000 (15:33 +0000)
src/gromacs/mdrun/md.cpp
src/gromacs/mdrun/mimic.cpp
src/gromacs/mdrun/minimize.cpp
src/gromacs/mdrun/rerun.cpp
src/gromacs/mdtypes/state.cpp
src/gromacs/mdtypes/state.h
src/gromacs/modularsimulator/freeenergyperturbationdata.cpp

index df33e9e78205f8be58675f6200d408f1c63430a8..425eda0de852359d92ff7fe84cc179def23a5793 100644 (file)
@@ -257,8 +257,15 @@ void gmx::LegacySimulator::do_md()
 
     int*                fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
     gmx::ArrayRef<real> lambda    = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
-    initialize_lambdas(
-            fplog, *ir, gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc), MASTER(cr), fep_state, lambda);
+    initialize_lambdas(fplog,
+                       ir->efep,
+                       ir->bSimTemp,
+                       *ir->fepvals,
+                       ir->simtempvals->temperatures,
+                       gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
+                       MASTER(cr),
+                       fep_state,
+                       lambda);
     Update upd(*ir, deform);
     bool   doSimulatedAnnealing = false;
     {
index 25feef03cd6163637adb110a880d630486b97d55..1ff4427896204c31bb75ee5a1cb95aebf1aa4d04 100644 (file)
@@ -229,7 +229,10 @@ void gmx::LegacySimulator::do_mimic()
     }
 
     initialize_lambdas(fplog,
-                       *ir,
+                       ir->efep,
+                       ir->bSimTemp,
+                       *ir->fepvals,
+                       ir->simtempvals->temperatures,
                        gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
                        MASTER(cr),
                        &state_global->fep_state,
index b115acdbec8a0f8740e15461362d531e2109186a..79eaaf0a57f833b720f18ffb2e4bf18ef77cf9c9 100644 (file)
@@ -394,8 +394,15 @@ static void init_em(FILE*                fplog,
     }
     int*                fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
     gmx::ArrayRef<real> lambda    = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
-    initialize_lambdas(
-            fplog, *ir, gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc), MASTER(cr), fep_state, lambda);
+    initialize_lambdas(fplog,
+                       ir->efep,
+                       ir->bSimTemp,
+                       *ir->fepvals,
+                       ir->simtempvals->temperatures,
+                       gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
+                       MASTER(cr),
+                       fep_state,
+                       lambda);
 
     if (ir->eI == IntegrationAlgorithm::NM)
     {
index 157d18e55bcf548dc19d4bcf9731587ae35c659a..74697796f8ebd1b93860b1249b784a540a22e943 100644 (file)
@@ -275,8 +275,15 @@ void gmx::LegacySimulator::do_rerun()
     }
     int*                fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
     gmx::ArrayRef<real> lambda    = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
-    initialize_lambdas(
-            fplog, *ir, gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc), MASTER(cr), fep_state, lambda);
+    initialize_lambdas(fplog,
+                       ir->efep,
+                       ir->bSimTemp,
+                       *ir->fepvals,
+                       ir->simtempvals->temperatures,
+                       gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
+                       MASTER(cr),
+                       fep_state,
+                       lambda);
     const bool        simulationsShareState = false;
     gmx_mdoutf*       outf                  = init_mdoutf(fplog,
                                    nfile,
index 9cca0564a995db5cdf1097c361476aec269b282f..fe207526272e86fd4abcc8fbb20f6717e9ab49cd 100644 (file)
@@ -344,26 +344,28 @@ void printLambdaStateToLog(FILE* fplog, gmx::ArrayRef<const real> lambda, const
     }
 }
 
-void initialize_lambdas(FILE*               fplog,
-                        const t_inputrec&   ir,
-                        gmx::ArrayRef<real> ref_t,
-                        bool                isMaster,
-                        int*                fep_state,
-                        gmx::ArrayRef<real> lambda)
+void initialize_lambdas(FILE*                            fplog,
+                        const FreeEnergyPerturbationType freeEnergyPerturbationType,
+                        const bool                       haveSimulatedTempering,
+                        const t_lambda&                  fep,
+                        gmx::ArrayRef<const real>        simulatedTemperingTemps,
+                        gmx::ArrayRef<real>              ref_t,
+                        bool                             isMaster,
+                        int*                             fep_state,
+                        gmx::ArrayRef<real>              lambda)
 {
     /* TODO: Clean up initialization of fep_state and lambda in
        t_state.  This function works, but could probably use a logic
        rewrite to keep all the different types of efep straight. */
 
-    if ((ir.efep == FreeEnergyPerturbationType::No) && (!ir.bSimTemp))
+    if ((freeEnergyPerturbationType == FreeEnergyPerturbationType::No) && (!haveSimulatedTempering))
     {
         return;
     }
 
-    const t_lambda* fep = ir.fepvals.get();
     if (isMaster)
     {
-        *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
+        *fep_state = fep.init_fep_state; /* this might overwrite the checkpoint
                                              if checkpoint is set -- a kludge is in for now
                                              to prevent this.*/
     }
@@ -372,27 +374,27 @@ void initialize_lambdas(FILE*               fplog,
     {
         double thisLambda;
         /* overwrite lambda state with init_lambda for now for backwards compatibility */
-        if (fep->init_lambda >= 0) /* if it's -1, it was never initialized */
+        if (fep.init_lambda >= 0) /* if it's -1, it was never initialized */
         {
-            thisLambda = fep->init_lambda;
+            thisLambda = fep.init_lambda;
         }
         else
         {
-            thisLambda = fep->all_lambda[i][fep->init_fep_state];
+            thisLambda = fep.all_lambda[i][fep.init_fep_state];
         }
         if (isMaster)
         {
             lambda[i] = thisLambda;
         }
     }
-    if (ir.bSimTemp)
+    if (haveSimulatedTempering)
     {
         /* need to rescale control temperatures to match current state */
         for (int i = 0; i < ref_t.ssize(); i++)
         {
             if (ref_t[i] > 0)
             {
-                ref_t[i] = ir.simtempvals->temperatures[fep->init_fep_state];
+                ref_t[i] = simulatedTemperingTemps[fep.init_fep_state];
             }
         }
     }
index 0dacb69dba5356b9987821a5d1bf87419f7061bb..4b0f8a3b9be3a844a072f67d73880ccb86e60631 100644 (file)
@@ -67,6 +67,8 @@
 #include "gromacs/utility/real.h"
 
 struct t_inputrec;
+struct t_lambda;
+enum class FreeEnergyPerturbationType;
 
 namespace gmx
 {
@@ -363,11 +365,14 @@ void printLambdaStateToLog(FILE* fplog, gmx::ArrayRef<const real> lambda, bool i
  * and lambda on master rank.
  *
  * Reports the initial lambda state to the log file. */
-void initialize_lambdas(FILE*               fplog,
-                        const t_inputrec&   ir,
-                        gmx::ArrayRef<real> ref_t,
-                        bool                isMaster,
-                        int*                fep_state,
-                        gmx::ArrayRef<real> lambda);
+void initialize_lambdas(FILE*                      fplog,
+                        FreeEnergyPerturbationType freeEnergyPerturbationType,
+                        bool                       haveSimulatedTempering,
+                        const t_lambda&            fep,
+                        gmx::ArrayRef<const real>  simulatedTemperingTemps,
+                        gmx::ArrayRef<real>        ref_t,
+                        bool                       isMaster,
+                        int*                       fep_state,
+                        gmx::ArrayRef<real>        lambda);
 
 #endif
index 6671c20d85783871dfdec2700849df9f27f2a353..de8fbdc008645ccb8d1af46c501c9aed427533eb 100644 (file)
@@ -73,7 +73,10 @@ FreeEnergyPerturbationData::FreeEnergyPerturbationData(FILE* fplog, const t_inpu
     // available on master. We have the lambda vector available everywhere, so we pass a `true`
     // for isMaster on all ranks. See #3647.
     initialize_lambdas(fplog_,
-                       inputrec_,
+                       inputrec_.efep,
+                       inputrec_.bSimTemp,
+                       *inputrec_.fepvals,
+                       inputrec_.simtempvals->temperatures,
                        gmx::arrayRefFromArray(inputrec_.opts.ref_t, inputrec_.opts.ngtc),
                        true,
                        &currentFEPState_,