Add _real user defined literal
authorRoland Schulz <roland.schulz@intel.com>
Mon, 6 Aug 2018 06:34:40 +0000 (23:34 -0700)
committerMark Abraham <mark.j.abraham@gmail.com>
Sun, 12 Aug 2018 13:15:47 +0000 (15:15 +0200)
Change-Id: Ia656d894db2f663a3b8cfc1f6897f20d97759edc

docs/dev-manual/language-features.rst
src/gromacs/gmxana/gmx_disre.cpp
src/gromacs/gmxana/gmx_genpr.cpp
src/gromacs/gmxana/gmx_rmsdist.cpp
src/gromacs/gmxana/gmx_xpm2ps.cpp
src/gromacs/math/invertmatrix.cpp
src/gromacs/mdlib/nbnxn_search.cpp
src/gromacs/trajectoryanalysis/modules/sasa.cpp
src/gromacs/utility/real.h

index e74699aeb1b6c7c81b37b7a0388ae4fd35c42b14..22a2d164d56b0aae3da7b7e63ced621d0fef24f3 100644 (file)
@@ -62,7 +62,8 @@ a release.
 * Don't use C-style casts; use ``const_cast``, ``static_cast`` or
   ``reinterpret_cast as appropriate``. See the point on RTTI for
   ``dynamic_cast``. For emphasizing type (e.g. intentional integer division)
-  use constructor syntax.
+  use constructor syntax. For creating real constants use the user-defined literal
+  _real (e.g. 2.5_real instead of static_cast<real>(2.5)).
 * Avoid overloading functions unless all variants really do the same
   thing, just with different types. Instead, consider making the
   function names more descriptive.
index 3f07b144f03a0b99a355f85c2a82d19f99ab2799..e12dd7e566db7e065384c68f450db8c05a4f754f 100644 (file)
@@ -618,7 +618,7 @@ static void dump_disre_matrix(const char *fn, t_dr_result *dr, int ndr,
             {
                 fprintf(debug, "DR %d, atoms %d, %d, distance %g\n", i, ai, aj, rav);
             }
-            rviol           = std::max(static_cast<real>(0.0), rav-idef->iparams[tp].disres.up1);
+            rviol           = std::max(0.0_real, rav-idef->iparams[tp].disres.up1);
             matrix[ri][rj] += w_dr[i]*rviol;
             matrix[rj][ri] += w_dr[i]*rviol;
             hi              = std::max(hi, matrix[ri][rj]);
index ff988d408c6aac3407bac842b5fc2baba81719d3..e22380046e84d8064489dd7f532b21bdcabaef5b 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,2014,2015,2017, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2017,2018, 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.
@@ -231,7 +231,7 @@ int gmx_genpr(int argc, char *argv[])
                         {
                             dd = disre_dist;
                         }
-                        lo = std::max(static_cast<real>(0.0), d-dd);
+                        lo = std::max(0.0_real, d-dd);
                         hi = d+dd;
                         fprintf(out, "%5d %5d %1d %5d %10d %10g %10g %10g %10g\n",
                                 ind_grp[i]+1, ind_grp[j]+1, 1, k, 1,
index 7027c241ff49cde99e3b05c2728c0361e05712b2..6cf2eae1edda8939f375bee485bd5b74b714590c 100644 (file)
@@ -590,7 +590,7 @@ static void calc_rms(int nind, int nframes,
         {
             mean  = dtot[i][j]/nframes;
             mean2 = dtot2[i][j]/nframes;
-            rms   = std::sqrt(std::max(static_cast<real>(0.0), mean2-mean*mean));
+            rms   = std::sqrt(std::max(0.0_real, mean2-mean*mean));
             rmsc  = rms/mean;
             if (mean > *meanmax)
             {
index ef76c0abd87311849766434f7aab935c10311037..5c04390a80b32b232d429183d88ec1ef14fac010 100644 (file)
@@ -756,7 +756,7 @@ static void tick_spacing(int n, real axis[], real offset, char axisnm,
     {
         for (f = 0; f < NFACT && bTryAgain; f++)
         {
-            space = std::pow(static_cast<real>(10.0), static_cast<real>(t)) * major_fact[f];
+            space = std::pow(10.0_real, static_cast<real>(t)) * major_fact[f];
             /* count how many ticks we would get: */
             i = 0;
             for (j = 0; j < n; j++)
index d4b3cefa1dc5401703422108ff8f0f59013b5511..5c181f672e48e26477821f18f68985622d6df84e 100644 (file)
@@ -73,12 +73,12 @@ void invertBoxMatrix(const matrix src, matrix dest)
 
 void invertMatrix(const matrix src, matrix dest)
 {
-    const real smallreal = static_cast<real>(1.0e-24);
-    const real largereal = static_cast<real>(1.0e24);
+    const real smallreal = 1.0e-24_real;
+    const real largereal = 1.0e24_real;
 
     real       determinant = det(src);
-    real       c           = static_cast<real>(1.0)/determinant;
-    real       fc          = static_cast<real>(std::fabs(c));
+    real       c           = 1.0_real/determinant;
+    real       fc          = std::fabs(c);
 
     if ((fc <= smallreal) || (fc >= largereal))
     {
index 81979e6431e86d3dc603cd98d0d87eaa441e5fd7..aa45047efee9742d2b6099ca33567035c4984c5d 100644 (file)
@@ -2756,7 +2756,7 @@ static void get_nsubpair_target(const nbnxn_search   *nbs,
          * groups of atoms we'll anyhow be limited by nsubpair_target_min,
          * so this overestimation will not matter.
          */
-        nsp_est = std::max(nsp_est, grid->nsubc_tot*static_cast<real>(14));
+        nsp_est = std::max(nsp_est, grid->nsubc_tot*14._real);
 
         if (debug)
         {
index 4e16edf92875f3f7a8bf3ecb6d7ad97da65d17bb..41ce97d2f2f6a801bc1ddf15961215da7319c2a8 100644 (file)
@@ -830,8 +830,7 @@ void computeAreas(const Selection &surfaceSel, const Selection &sel,
 
     if (bResAt)
     {
-        std::fill(resAreaWork->begin(), resAreaWork->end(),
-                  static_cast<real>(0.0));
+        std::fill(resAreaWork->begin(), resAreaWork->end(), 0.0_real);
     }
     for (int i = 0; i < sel.posCount(); ++i)
     {
@@ -938,7 +937,7 @@ Sasa::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
         if (surfaceSel.isDynamic())
         {
             std::fill(frameData.atomAreas_.begin(), frameData.atomAreas_.end(),
-                      static_cast<real>(0.0));
+                      0.0_real);
             for (size_t i = 0; i < frameData.index_.size(); ++i)
             {
                 frameData.atomAreas_[frameData.index_[i]] = area[i];
index 682514f26dba086018c76b9772105054d3fa6964..34ba0b070505a8e68dbb7aba15204a39db05e719 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,2014,2015,2017, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2017,2018, 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.
@@ -140,4 +140,19 @@ typedef float           real;
 
 #endif /* GMX_DOUBLE */
 
+#ifdef __cplusplus
+/*! \brief User defined literal for real numbers.
+ *
+ * Examples: 2._real, 2.5_real, .5_real. The number is always of type real.
+ *
+ * It is best to use a real constant whenever it is used only with operands which are real.
+ * If a constant is double than the compiler is forced to do operations directly involving the constant
+ * in double even if all variables are real. A constant shouldn't be real when used with double operands,
+ * because then the constant is less accurate with GMX_DOUBLE=no.
+ *
+ * See https://en.cppreference.com/w/cpp/language/user_literal for details on this lanuage feature.
+ */
+constexpr real operator"" _real(long double x) { return real(x); }
+#endif
+
 #endif