Introduce TimeStep and Offset strong types
authorPascal Merz <pascal.merz@me.com>
Sat, 28 Nov 2020 19:13:39 +0000 (12:13 -0700)
committerAndrey Alekseenko <al42and@gmail.com>
Wed, 21 Apr 2021 11:23:53 +0000 (11:23 +0000)
When building integration algorithms in modular simulator, the time
step (for propagations) and the scheduling offset (for algorithms called
periodically) are the most frequently used type of arguments. They are
also currently the only arguments which are passed by simple types
(int and double, respectively).

The current change introduces strong types, increasing both safety and
readability. Note that these strong types can easily be adapted to
inheriting from a general base class if one is introduced (refs #4005).

src/gromacs/modularsimulator/modularsimulator.cpp
src/gromacs/modularsimulator/modularsimulatorinterfaces.h
src/gromacs/modularsimulator/parrinellorahmanbarostat.cpp
src/gromacs/modularsimulator/parrinellorahmanbarostat.h
src/gromacs/modularsimulator/propagator.cpp
src/gromacs/modularsimulator/propagator.h
src/gromacs/modularsimulator/velocityscalingtemperaturecoupling.cpp
src/gromacs/modularsimulator/velocityscalingtemperaturecoupling.h

index d9587f7f8bc3ea1d44c4a0e00946db745d5543a7..503fd87547d9a92dbe99904cf3c27866980f0784 100644 (file)
@@ -115,13 +115,13 @@ void ModularSimulator::addIntegrationElements(ModularSimulatorAlgorithmBuilder*
             || legacySimulatorData_->inputrec->etc == TemperatureCoupling::Berendsen
             || legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
         {
-            builder->add<VelocityScalingTemperatureCoupling>(-1,
+            builder->add<VelocityScalingTemperatureCoupling>(Offset(-1),
                                                              UseFullStepKE::No,
                                                              ReportPreviousStepConservedEnergy::No,
                                                              PropagatorTag("LeapFrogPropagator"));
         }
         builder->add<Propagator<IntegrationStage::LeapFrog>>(
-                PropagatorTag("LeapFrogPropagator"), legacySimulatorData_->inputrec->delta_t);
+                PropagatorTag("LeapFrogPropagator"), TimeStep(legacySimulatorData_->inputrec->delta_t));
         if (legacySimulatorData_->constr)
         {
             builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
@@ -129,7 +129,7 @@ void ModularSimulator::addIntegrationElements(ModularSimulatorAlgorithmBuilder*
         builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>>();
         if (legacySimulatorData_->inputrec->epc == PressureCoupling::ParrinelloRahman)
         {
-            builder->add<ParrinelloRahmanBarostat>(-1, PropagatorTag("LeapFrogPropagator"));
+            builder->add<ParrinelloRahmanBarostat>(Offset(-1), PropagatorTag("LeapFrogPropagator"));
         }
     }
     else if (legacySimulatorData_->inputrec->eI == IntegrationAlgorithm::VV)
@@ -137,7 +137,7 @@ void ModularSimulator::addIntegrationElements(ModularSimulatorAlgorithmBuilder*
         // The velocity verlet integration algorithm
         builder->add<ForceElement>();
         builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
-                PropagatorTag("VelocityHalfStep"), 0.5 * legacySimulatorData_->inputrec->delta_t);
+                PropagatorTag("VelocityHalfStep"), TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
         if (legacySimulatorData_->constr)
         {
             builder->add<ConstraintsElement<ConstraintVariable::Velocities>>();
@@ -148,13 +148,14 @@ void ModularSimulator::addIntegrationElements(ModularSimulatorAlgorithmBuilder*
             || legacySimulatorData_->inputrec->etc == TemperatureCoupling::Berendsen)
         {
             builder->add<VelocityScalingTemperatureCoupling>(
-                    0,
+                    Offset(0),
                     UseFullStepKE::Yes,
                     ReportPreviousStepConservedEnergy::Yes,
                     PropagatorTag("VelocityHalfAndPositionFullStep"));
         }
         builder->add<Propagator<IntegrationStage::VelocityVerletPositionsAndVelocities>>(
-                PropagatorTag("VelocityHalfAndPositionFullStep"), legacySimulatorData_->inputrec->delta_t);
+                PropagatorTag("VelocityHalfAndPositionFullStep"),
+                TimeStep(legacySimulatorData_->inputrec->delta_t));
         if (legacySimulatorData_->constr)
         {
             builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
@@ -162,7 +163,7 @@ void ModularSimulator::addIntegrationElements(ModularSimulatorAlgorithmBuilder*
         builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
         if (legacySimulatorData_->inputrec->epc == PressureCoupling::ParrinelloRahman)
         {
-            builder->add<ParrinelloRahmanBarostat>(-1, PropagatorTag("VelocityHalfStep"));
+            builder->add<ParrinelloRahmanBarostat>(Offset(-1), PropagatorTag("VelocityHalfStep"));
         }
     }
     else
index fc66d253745c77371a00dc65d849a289e280969e..44fac47aaa4e257ed02d3f331e32ac9dd47f81e1 100644 (file)
@@ -474,6 +474,36 @@ private:
     const std::string name_;
 };
 
+/*! \internal
+ * \brief Strong type used to denote propagation time steps
+ */
+struct TimeStep
+{
+    //! Explicit constructor
+    explicit TimeStep(real timeStep) : timeStep_(timeStep) {}
+    //! Can be used as underlying type
+    operator const real&() const { return timeStep_; }
+
+private:
+    //! The time step
+    real timeStep_;
+};
+
+/*! \internal
+ * \brief Strong type used to denote scheduling offsets
+ */
+struct Offset
+{
+    //! Explicit constructor
+    explicit Offset(int offset) : offset_(offset) {}
+    //! Can be used as underlying type
+    operator const int&() const { return offset_; }
+
+private:
+    //! The offset
+    int offset_;
+};
+
 /*! \internal
  * \brief Information needed to connect a propagator to a thermostat
  */
index f35d1846ae27e7110e5e3c91b662093424fcfba7..72f2d3bf8b3295824b8c540d5a59d13f7d83f498 100644 (file)
@@ -329,7 +329,7 @@ ISimulatorElement* ParrinelloRahmanBarostat::getElementPointerImpl(
         EnergyData*                             energyData,
         FreeEnergyPerturbationData gmx_unused* freeEnergyPerturbationData,
         GlobalCommunicationHelper gmx_unused* globalCommunicationHelper,
-        int                                   offset,
+        Offset                                offset,
         const PropagatorTag&                  propagatorTag)
 {
     auto* element  = builderHelper->storeElement(std::make_unique<ParrinelloRahmanBarostat>(
index 91bc1a1223f5036e169acd0acc27f339c7ea9deb..211ac669dd5b0677bcf5f0f032d10f9e4b0316d1 100644 (file)
@@ -131,7 +131,7 @@ public:
                           EnergyData*                             energyData,
                           FreeEnergyPerturbationData gmx_unused* freeEnergyPerturbationData,
                           GlobalCommunicationHelper gmx_unused* globalCommunicationHelper,
-                          int                                   offset,
+                          Offset                                offset,
                           const PropagatorTag&                  propagatorTag);
 
 private:
index 476d74c3a3ec1161a17ebd0f8b7c8e03aa3f888f..2b2df893482483c8709076ad48052442993a803a 100644 (file)
@@ -907,7 +907,7 @@ ISimulatorElement* Propagator<integrationStage>::getElementPointerImpl(
         FreeEnergyPerturbationData gmx_unused* freeEnergyPerturbationData,
         GlobalCommunicationHelper gmx_unused* globalCommunicationHelper,
         const PropagatorTag&                  propagatorTag,
-        double                                timestep)
+        TimeStep                              timestep)
 {
     GMX_RELEASE_ASSERT(!(integrationStage == IntegrationStage::ScaleVelocities
                          || integrationStage == IntegrationStage::ScalePositions)
@@ -952,7 +952,7 @@ ISimulatorElement* Propagator<integrationStage>::getElementPointerImpl(
                                  freeEnergyPerturbationData,
                                  globalCommunicationHelper,
                                  propagatorTag,
-                                 0.0);
+                                 TimeStep(0.0));
 }
 //! \endcond
 
index 19cd106a22dadcbc748f2f49027508386df360bf..631a98f133782fface5d4f07495e118585b1a6b1 100644 (file)
@@ -191,7 +191,7 @@ public:
                                                     FreeEnergyPerturbationData* freeEnergyPerturbationData,
                                                     GlobalCommunicationHelper* globalCommunicationHelper,
                                                     const PropagatorTag&       propagatorTag,
-                                                    double                     timestep);
+                                                    TimeStep                   timestep);
 
     /*! \brief Factory method implementation
      *
index 5e77ae63b912211aa213c22f6525227824bc5540..9eef0d2707655c7f7e30c34573e2cfd84fd2d106 100644 (file)
@@ -658,7 +658,7 @@ ISimulatorElement* VelocityScalingTemperatureCoupling::getElementPointerImpl(
         EnergyData*                     energyData,
         FreeEnergyPerturbationData gmx_unused* freeEnergyPerturbationData,
         GlobalCommunicationHelper gmx_unused* globalCommunicationHelper,
-        int                                   offset,
+        Offset                                offset,
         UseFullStepKE                         useFullStepKE,
         ReportPreviousStepConservedEnergy     reportPreviousStepConservedEnergy,
         const PropagatorTag&                  propagatorTag)
index 6b70509a91f15bed981552bbcf79a2b12902c6ed..8f068fa3315b68e2d031c0df066736dd4f1fdde3 100644 (file)
@@ -152,7 +152,7 @@ public:
                           EnergyData*                             energyData,
                           FreeEnergyPerturbationData*             freeEnergyPerturbationData,
                           GlobalCommunicationHelper*              globalCommunicationHelper,
-                          int                                     offset,
+                          Offset                                  offset,
                           UseFullStepKE                           useFullStepKE,
                           ReportPreviousStepConservedEnergy       reportPreviousStepConservedEnergy,
                           const PropagatorTag&                    propagatorTag);