More sc-gapsys templates plus fix unused warnigs.
authorSebastian Kehl <sebastian.kehl@mpcdf.mpg.de>
Fri, 1 Oct 2021 17:41:55 +0000 (19:41 +0200)
committerMagnus Lundborg <magnus.lundborg@scilifelab.se>
Tue, 5 Oct 2021 11:38:46 +0000 (11:38 +0000)
src/gromacs/gmxlib/nonbonded/nb_softcore.h

index dabbaa28534a0026ab965c312dbc03690d7a31b7..501ddb1eb9dfe5e1288b27f0817438bda87b18ba 100644 (file)
 #include "gromacs/simd/simd_math.h"
 
 /* linearized electrostatics */
-template<class RealType, class BoolType>
-static inline void quadraticApproximationCoulomb(const RealType qq,
-                                                 const RealType rInvQ,
-                                                 const RealType r,
-                                                 const real     lambdaFac,
-                                                 const real     dLambdaFac,
-                                                 RealType*      force,
-                                                 RealType*      potential,
-                                                 RealType*      dvdl,
-                                                 BoolType       dvdlMask)
+template<bool computeForces, class RealType, class BoolType>
+static inline void quadraticApproximationCoulomb(const RealType       qq,
+                                                 const RealType       rInvQ,
+                                                 const RealType       r,
+                                                 const real           lambdaFac,
+                                                 const real           dLambdaFac,
+                                                 RealType gmx_unused* force,
+                                                 RealType*            potential,
+                                                 RealType*            dvdl,
+                                                 BoolType             dvdlMask)
 {
     RealType constFac = qq * rInvQ;
     RealType linFac   = constFac * r * rInvQ;
     RealType quadrFac = linFac * r * rInvQ;
 
     /* Computing Coulomb force and potential energy */
-    *force = -2 * quadrFac + 3 * linFac;
+    if constexpr (computeForces)
+    {
+        *force = -2 * quadrFac + 3 * linFac;
+    }
 
     *potential = quadrFac - 3 * (linFac - constFac);
 
@@ -69,19 +72,19 @@ static inline void quadraticApproximationCoulomb(const RealType qq,
 
 /* reaction-field linearized electrostatics */
 template<bool computeForces, class RealType, class BoolType>
-static inline void reactionFieldQuadraticPotential(const RealType qq,
-                                                   const real     facel,
-                                                   const RealType r,
-                                                   const real     rCutoff,
-                                                   const real     lambdaFac,
-                                                   const real     dLambdaFac,
-                                                   const RealType alphaEff,
-                                                   const real     krf,
-                                                   const real     potentialShift,
-                                                   RealType*      force,
-                                                   RealType*      potential,
-                                                   RealType*      dvdl,
-                                                   BoolType       mask)
+static inline void reactionFieldQuadraticPotential(const RealType       qq,
+                                                   const real           facel,
+                                                   const RealType       r,
+                                                   const real           rCutoff,
+                                                   const real           lambdaFac,
+                                                   const real           dLambdaFac,
+                                                   const RealType       alphaEff,
+                                                   const real           krf,
+                                                   const real           potentialShift,
+                                                   RealType gmx_unused* force,
+                                                   RealType*            potential,
+                                                   RealType*            dvdl,
+                                                   BoolType             mask)
 {
     /* check if we have to use the hardcore values */
     BoolType computeValues = mask && (lambdaFac < 1 && 0 < alphaEff && facel != 0);
@@ -110,7 +113,7 @@ static inline void reactionFieldQuadraticPotential(const RealType qq,
             RealType forceQuad(0);
             RealType potentialQuad(0);
             RealType dvdlQuad(0);
-            quadraticApproximationCoulomb(
+            quadraticApproximationCoulomb<computeForces>(
                     qq, rInvQ, r, lambdaFac, dLambdaFac, &forceQuad, &potentialQuad, &dvdlQuad, computeValues);
 
             // rf modification
@@ -120,7 +123,7 @@ static inline void reactionFieldQuadraticPotential(const RealType qq,
             // update
             if constexpr (computeForces)
             {
-                *force     = gmx::blend(*force, forceQuad, computeValues);
+                *force = gmx::blend(*force, forceQuad, computeValues);
             }
             *potential = gmx::blend(*potential, potentialQuad, computeValues);
             *dvdl      = *dvdl + gmx::selectByMask(dvdlQuad, computeValues && withinCutoff);
@@ -130,18 +133,18 @@ static inline void reactionFieldQuadraticPotential(const RealType qq,
 
 /* ewald linearized electrostatics */
 template<bool computeForces, class RealType, class BoolType>
-static inline void ewaldQuadraticPotential(const RealType qq,
-                                           const real     facel,
-                                           const RealType r,
-                                           const real     rCutoff,
-                                           const real     lambdaFac,
-                                           const real     dLambdaFac,
-                                           const RealType alphaEff,
-                                           const real     potentialShift,
-                                           RealType*      force,
-                                           RealType*      potential,
-                                           RealType*      dvdl,
-                                           BoolType       mask)
+static inline void ewaldQuadraticPotential(const RealType       qq,
+                                           const real           facel,
+                                           const RealType       r,
+                                           const real           rCutoff,
+                                           const real           lambdaFac,
+                                           const real           dLambdaFac,
+                                           const RealType       alphaEff,
+                                           const real           potentialShift,
+                                           RealType gmx_unused* force,
+                                           RealType*            potential,
+                                           RealType*            dvdl,
+                                           BoolType             mask)
 {
     /* check if we have to use the hardcore values */
     BoolType computeValues = mask && (lambdaFac < 1 && 0 < alphaEff && facel != 0);
@@ -170,7 +173,7 @@ static inline void ewaldQuadraticPotential(const RealType qq,
             RealType forceQuad(0);
             RealType potentialQuad(0);
             RealType dvdlQuad(0);
-            quadraticApproximationCoulomb(
+            quadraticApproximationCoulomb<computeForces>(
                     qq, rInvQ, r, lambdaFac, dLambdaFac, &forceQuad, &potentialQuad, &dvdlQuad, computeValues);
 
             // ewald modification
@@ -179,7 +182,7 @@ static inline void ewaldQuadraticPotential(const RealType qq,
             // update
             if constexpr (computeForces)
             {
-                *force     = gmx::blend(*force, forceQuad, computeValues);
+                *force = gmx::blend(*force, forceQuad, computeValues);
             }
             *potential = gmx::blend(*potential, potentialQuad, computeValues);
             *dvdl      = *dvdl + gmx::selectByMask(dvdlQuad, computeValues && withinCutoff);
@@ -189,20 +192,20 @@ static inline void ewaldQuadraticPotential(const RealType qq,
 
 /* cutoff LJ with quadratic appximation of lj-potential */
 template<bool computeForces, class RealType, class BoolType>
-static inline void lennardJonesQuadraticPotential(const RealType c6,
-                                                  const RealType c12,
-                                                  const RealType r,
-                                                  const RealType rsq,
-                                                  const real     lambdaFac,
-                                                  const real     dLambdaFac,
-                                                  const RealType sigma6,
-                                                  const RealType alphaEff,
-                                                  const real     repulsionShift,
-                                                  const real     dispersionShift,
-                                                  RealType*      force,
-                                                  RealType*      potential,
-                                                  RealType*      dvdl,
-                                                  BoolType       mask)
+static inline void lennardJonesQuadraticPotential(const RealType       c6,
+                                                  const RealType       c12,
+                                                  const RealType       r,
+                                                  const RealType       rsq,
+                                                  const real           lambdaFac,
+                                                  const real           dLambdaFac,
+                                                  const RealType       sigma6,
+                                                  const RealType       alphaEff,
+                                                  const real           repulsionShift,
+                                                  const real           dispersionShift,
+                                                  RealType gmx_unused* force,
+                                                  RealType*            potential,
+                                                  RealType*            dvdl,
+                                                  BoolType             mask)
 {
     constexpr real c_twentySixSeventh = 26.0_real / 7.0_real;
     constexpr real c_oneSixth         = 1.0_real / 6.0_real;
@@ -249,9 +252,9 @@ static inline void lennardJonesQuadraticPotential(const RealType c6,
             constFac  = 91 * rInv12C - 28 * rInv6C;
 
             /* Computing LJ force and potential energy */
-            RealType forceQuad     = -quadrFac + linearFac;
-            RealType potentialQuad = c_half * quadrFac - linearFac + constFac;
-            RealType dvdlQuad      = dLambdaFac * 28 * (lambdaFac * lambdaFacRevInv)
+            RealType gmx_unused forceQuad     = -quadrFac + linearFac;
+            RealType            potentialQuad = c_half * quadrFac - linearFac + constFac;
+            RealType            dvdlQuad      = dLambdaFac * 28 * (lambdaFac * lambdaFacRevInv)
                                 * ((6.5_real * rInv14C - rInv8C) - (13 * rInv13C - 2 * rInv7C)
                                    + (6.5_real * rInv12C - rInv6C));
 
@@ -260,7 +263,7 @@ static inline void lennardJonesQuadraticPotential(const RealType c6,
                                                 computeValues);
             if constexpr (computeForces)
             {
-                *force     = gmx::blend(*force, forceQuad, computeValues);
+                *force = gmx::blend(*force, forceQuad, computeValues);
             }
             *potential = gmx::blend(*potential, potentialQuad, computeValues);
             *dvdl      = *dvdl + gmx::selectByMask(dvdlQuad, computeValues);