Print lambda state to log file if using AWH for FEP.
authorMagnus Lundborg <magnus.lundborg@scilifelab.se>
Tue, 1 Sep 2020 11:24:20 +0000 (13:24 +0200)
committerPaul Bauer <paul.bauer.q@gmail.com>
Thu, 3 Sep 2020 11:35:06 +0000 (11:35 +0000)
Change-Id: I08c7e6d51a299eb53556d5bde63c29f7c7640f6e

src/gromacs/awh/awh.h
src/gromacs/mdrun/md.cpp
src/gromacs/mdtypes/state.cpp
src/gromacs/mdtypes/state.h

index d8bdfb71ed1e0dd864ae98be6ca5471a4d07c7c3..efd7c80a43c38ced4a439d3cdbf8644013412a30 100644 (file)
@@ -264,6 +264,10 @@ public:
      */
     bool needForeignEnergyDifferences(int64_t step) const;
 
+    /*! \brief Returns true if AWH has a bias with a free energy lambda state dimension at all.
+     */
+    bool hasFepLambdaDimension() const;
+
 private:
     /*! \brief Returns whether we need to write output at the current step.
      *
@@ -271,10 +275,6 @@ private:
      */
     bool isOutputStep(int64_t step) const;
 
-    /*! \brief Returns true if AWH has a bias with a free energy lambda state dimension at all.
-     */
-    bool hasFepLambdaDimension() const;
-
 private:
     std::vector<BiasCoupledToSystem> biasCoupledToSystem_; /**< AWH biases and definitions of their coupling to the system. */
     const int64_t    seed_;   /**< Random seed for MC jumping with umbrella type bias potential. */
index 4b907ec397e919af19fe3a3b242d414e9dbebc39..3d1c95a411a8b2e951db9efa242174ce97623f8c 100644 (file)
@@ -1532,6 +1532,11 @@ void gmx::LegacySimulator::do_md()
                                                    do_log ? fplog : nullptr, step, t,
                                                    fr->fcdata.get(), awh.get());
             }
+            if (do_log && ir->bDoAwh && awh->hasFepLambdaDimension())
+            {
+                const bool isInitialOutput = false;
+                printLambdaStateToLog(fplog, state->lambda, isInitialOutput);
+            }
 
             if (ir->bPull)
             {
index 61c280f997638ae99c8587c0f599acba707fd351..5e80a742083019c6035677a379cf248bd59cc7d8 100644 (file)
@@ -280,6 +280,19 @@ void preserve_box_shape(const t_inputrec* ir, matrix box_rel, matrix box)
     }
 }
 
+void printLambdaStateToLog(FILE* fplog, const gmx::ArrayRef<real> lambda, const bool isInitialOutput)
+{
+    if (fplog != nullptr)
+    {
+        fprintf(fplog, "%s vector of lambda components:[ ", isInitialOutput ? "Initial" : "Current");
+        for (const auto& l : lambda)
+        {
+            fprintf(fplog, "%10.4f ", l);
+        }
+        fprintf(fplog, "]\n%s", isInitialOutput ? "" : "\n");
+    }
+}
+
 void initialize_lambdas(FILE*               fplog,
                         const t_inputrec&   ir,
                         bool                isMaster,
@@ -338,13 +351,6 @@ void initialize_lambdas(FILE*               fplog,
     }
 
     /* Send to the log the information on the current lambdas */
-    if (fplog != nullptr)
-    {
-        fprintf(fplog, "Initial vector of lambda components:[ ");
-        for (const auto& l : lambda)
-        {
-            fprintf(fplog, "%10.4f ", l);
-        }
-        fprintf(fplog, "]\n");
-    }
+    const bool isInitialOutput = true;
+    printLambdaStateToLog(fplog, lambda, isInitialOutput);
 }
index 016711df7e3137df3fd11bd18d7be0e83e997353..7be1407fbf8198967cee46bce723ab0667283a13 100644 (file)
@@ -341,6 +341,15 @@ static inline gmx::ArrayRef<const gmx::RVec> positionsFromStatePointer(const t_s
     }
 };
 
+/*! \brief Prints the current lambda state to the log file.
+ *
+ * \param[in] fplog  The log file. If fplog == nullptr there will be no output.
+ * \param[in] lambda The array of lambda values.
+ * \param[in] isInitialOutput Whether this output is the initial lambda state or not.
+ */
+void printLambdaStateToLog(FILE* fplog, gmx::ArrayRef<real> lambda, bool isInitialOutput);
+
+
 /*! \brief Fills fep_state, lambda, and lam0 if needed
  *
  * If FEP or simulated tempering is in use: