Modernize t_nrnb
authorJoe Jordan <ejjordan12@gmail.com>
Tue, 23 Mar 2021 10:15:33 +0000 (10:15 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Tue, 23 Mar 2021 10:15:33 +0000 (10:15 +0000)
* Use array instead of pointer.
* Remove unused functions.
* Add doxygen to inc_nrnb.
* Use inline function instead of macro in all places where possible.
* Remove unused headers.

src/gromacs/gmxlib/nonbonded/nb_free_energy.cpp
src/gromacs/gmxlib/nrnb.cpp
src/gromacs/gmxlib/nrnb.h
src/gromacs/mdrun/runner.cpp

index f194c51d572758f0bc910f4a51093169036cb62b..795cdbce4e16b6e2bbae6cce6c6f3a209abc89a0 100644 (file)
@@ -850,8 +850,7 @@ static void nb_free_energy_kernel(const t_nblist* gmx_restrict nlist,
      * 12  flops per outer iteration
      * 150 flops per inner iteration
      */
-#pragma omp atomic
-    inc_nrnb(nrnb, eNR_NBKERNEL_FREE_ENERGY, nlist->nri * 12 + nlist->jindex[nri] * 150);
+    atomicNrnbIncrement(nrnb, eNR_NBKERNEL_FREE_ENERGY, nlist->nri * 12 + nlist->jindex[nri] * 150);
 
     if (numExcludedPairsBeyondRlist > 0)
     {
index 59ff411b15c6657ec9014d1daa4f603d25e4f2ab..5239e7506da312931ab2a16a4caa75b43a2fea71 100644 (file)
@@ -4,7 +4,7 @@
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
  * Copyright (c) 2013,2014,2015,2017,2018 by the GROMACS development team.
- * Copyright (c) 2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2019,2020,2021, 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.
 #include <cstdlib>
 #include <cstring>
 
-#include <algorithm>
-#include <vector>
-
-#include "gromacs/mdtypes/commrec.h"
-#include "gromacs/mdtypes/md_enums.h"
-#include "gromacs/utility/arraysize.h"
-
 typedef struct
 {
     const char* name;
@@ -205,8 +198,8 @@ static void pr_two(FILE* out, int c, int i)
 
 static void pr_difftime(FILE* out, double dt)
 {
-    int      ndays, nhours, nmins, nsecs;
-    gmx_bool bPrint, bPrinted;
+    int  ndays, nhours, nmins, nsecs;
+    bool bPrint, bPrinted;
 
     ndays    = static_cast<int>(dt / (24 * 3600));
     dt       = dt - 24 * 3600 * ndays;
@@ -260,29 +253,15 @@ static void pr_difftime(FILE* out, double dt)
 
 void clear_nrnb(t_nrnb* nrnb)
 {
-    int i;
-
-    for (i = 0; (i < eNRNB); i++)
+    for (int i = 0; (i < eNRNB); i++)
     {
         nrnb->n[i] = 0.0;
     }
 }
 
-void add_nrnb(t_nrnb* dest, t_nrnb* s1, t_nrnb* s2)
-{
-    int i;
-
-    for (i = 0; (i < eNRNB); i++)
-    {
-        dest->n[i] = s1->n[i] + s2->n[i];
-    }
-}
-
 void print_nrnb(FILE* out, t_nrnb* nrnb)
 {
-    int i;
-
-    for (i = 0; (i < eNRNB); i++)
+    for (int i = 0; (i < eNRNB); i++)
     {
         if (nrnb->n[i] > 0)
         {
@@ -291,35 +270,32 @@ void print_nrnb(FILE* out, t_nrnb* nrnb)
     }
 }
 
-void _inc_nrnb(t_nrnb* nrnb, int enr, int inc, char gmx_unused* file, int gmx_unused line)
-{
-    nrnb->n[enr] += inc;
-#ifdef DEBUG_NRNB
-    printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n", nbdata[enr].name, enr, inc, file, line);
-#endif
-}
-
 /* Returns in enr is the index of a full nbnxn VdW kernel */
-static gmx_bool nrnb_is_nbnxn_vdw_kernel(int enr)
+static bool nrnb_is_nbnxn_vdw_kernel(int enr)
 {
     return (enr >= eNR_NBNXN_LJ_RF && enr <= eNR_NBNXN_LJ_E);
 }
 
 /* Returns in enr is the index of an nbnxn kernel addition (LJ modification) */
-static gmx_bool nrnb_is_nbnxn_kernel_addition(int enr)
+static bool nrnb_is_nbnxn_kernel_addition(int enr)
 {
     return (enr >= eNR_NBNXN_ADD_LJ_FSW && enr <= eNR_NBNXN_ADD_LJ_EWALD_E);
 }
 
+void atomicNrnbIncrement(t_nrnb* nrnb, int index, int increment)
+{
+#pragma omp atomic
+    nrnb->n[index] += increment;
+}
+
 void print_flop(FILE* out, t_nrnb* nrnb, double* nbfs, double* mflop)
 {
-    int         i, j;
     double      mni, frac, tfrac, tflop;
     const char* myline =
             "-----------------------------------------------------------------------------";
 
     *nbfs = 0.0;
-    for (i = 0; (i < eNR_NBKERNEL_TOTAL_NR); i++)
+    for (int i = 0; (i < eNR_NBKERNEL_TOTAL_NR); i++)
     {
         if (std::strstr(nbdata[i].name, "W3-W3") != nullptr)
         {
@@ -343,7 +319,7 @@ void print_flop(FILE* out, t_nrnb* nrnb, double* nbfs, double* mflop)
         }
     }
     tflop = 0;
-    for (i = 0; (i < eNRNB); i++)
+    for (int i = 0; (i < eNRNB); i++)
     {
         tflop += 1e-6 * nrnb->n[i] * nbdata[i].flop;
     }
@@ -370,7 +346,7 @@ void print_flop(FILE* out, t_nrnb* nrnb, double* nbfs, double* mflop)
     }
     *mflop = 0.0;
     tfrac  = 0.0;
-    for (i = 0; (i < eNRNB); i++)
+    for (int i = 0; (i < eNRNB); i++)
     {
         mni = 1e-6 * nrnb->n[i];
         /* Skip empty entries and nbnxn additional flops,
@@ -384,7 +360,7 @@ void print_flop(FILE* out, t_nrnb* nrnb, double* nbfs, double* mflop)
             if (nrnb_is_nbnxn_vdw_kernel(i))
             {
                 /* Possibly add the cost of an LJ switch/Ewald function */
-                for (j = eNR_NBNXN_ADD_LJ_FSW; j <= eNR_NBNXN_ADD_LJ_EWALD; j += 2)
+                for (int j = eNR_NBNXN_ADD_LJ_FSW; j <= eNR_NBNXN_ADD_LJ_EWALD; j += 2)
                 {
                     int e_kernel_add;
 
@@ -512,131 +488,3 @@ const char* nrnb_str(int enr)
 {
     return nbdata[enr].name;
 }
-
-static const int force_index[] = {
-    eNR_BONDS,  eNR_ANGLES, eNR_PROPER, eNR_IMPROPER, eNR_RB,
-    eNR_DISRES, eNR_ORIRES, eNR_POSRES, eNR_FBPOSRES, eNR_NS,
-};
-#define NFORCE_INDEX asize(force_index)
-
-static const int constr_index[] = { eNR_SHAKE,  eNR_SHAKE_RIJ,  eNR_SETTLE,  eNR_UPDATE,
-                                    eNR_PCOUPL, eNR_CONSTR_VIR, eNR_CONSTR_V };
-#define NCONSTR_INDEX asize(constr_index)
-
-static double pr_av(FILE* log, t_commrec* cr, double fav, const double ftot[], const char* title)
-{
-    int    i, perc;
-    double dperc, unb;
-
-    unb = 0;
-    if (fav > 0)
-    {
-        fav /= cr->nnodes - cr->npmenodes;
-        fprintf(log, "\n %-26s", title);
-        for (i = 0; (i < cr->nnodes); i++)
-        {
-            dperc = (100.0 * ftot[i]) / fav;
-            unb   = std::max(unb, dperc);
-            perc  = static_cast<int>(dperc);
-            fprintf(log, "%3d ", perc);
-        }
-        if (unb > 0)
-        {
-            perc = static_cast<int>(10000.0 / unb);
-            fprintf(log, "%6d%%\n\n", perc);
-        }
-        else
-        {
-            fprintf(log, "\n\n");
-        }
-    }
-    return unb;
-}
-
-void pr_load(FILE* log, t_commrec* cr, t_nrnb nrnb[])
-{
-    t_nrnb av;
-
-    std::vector<double> ftot(cr->nnodes);
-    std::vector<double> stot(cr->nnodes);
-    for (int i = 0; (i < cr->nnodes); i++)
-    {
-        add_nrnb(&av, &av, &(nrnb[i]));
-        /* Cost due to forces */
-        for (int j = 0; (j < eNR_NBKERNEL_TOTAL_NR); j++)
-        {
-            ftot[i] += nrnb[i].n[j] * cost_nrnb(j);
-        }
-        for (int j = 0; (j < NFORCE_INDEX); j++)
-        {
-            ftot[i] += nrnb[i].n[force_index[j]] * cost_nrnb(force_index[j]);
-        }
-        /* Due to shake */
-        for (int j = 0; (j < NCONSTR_INDEX); j++)
-        {
-            stot[i] += nrnb[i].n[constr_index[j]] * cost_nrnb(constr_index[j]);
-        }
-    }
-    for (int j = 0; (j < eNRNB); j++)
-    {
-        av.n[j] = av.n[j] / static_cast<double>(cr->nnodes - cr->npmenodes);
-    }
-
-    fprintf(log, "\nDetailed load balancing info in percentage of average\n");
-
-    fprintf(log, " Type                 RANK:");
-    for (int i = 0; (i < cr->nnodes); i++)
-    {
-        fprintf(log, "%3d ", i);
-    }
-    fprintf(log, "Scaling\n");
-    fprintf(log, "---------------------------");
-    for (int i = 0; (i < cr->nnodes); i++)
-    {
-        fprintf(log, "----");
-    }
-    fprintf(log, "-------\n");
-
-    for (int j = 0; (j < eNRNB); j++)
-    {
-        double unb = 100.0;
-        if (av.n[j] > 0)
-        {
-            double dperc;
-            int    perc;
-            fprintf(log, " %-26s", nrnb_str(j));
-            for (int i = 0; (i < cr->nnodes); i++)
-            {
-                dperc = (100.0 * nrnb[i].n[j]) / av.n[j];
-                unb   = std::max(unb, dperc);
-                perc  = static_cast<int>(dperc);
-                fprintf(log, "%3d ", perc);
-            }
-            if (unb > 0)
-            {
-                perc = static_cast<int>(10000.0 / unb);
-                fprintf(log, "%6d%%\n", perc);
-            }
-            else
-            {
-                fprintf(log, "\n");
-            }
-        }
-    }
-    double fav = 0;
-    double sav = 0;
-    for (int i = 0; (i < cr->nnodes); i++)
-    {
-        fav += ftot[i];
-        sav += stot[i];
-    }
-    double uf = pr_av(log, cr, fav, ftot.data(), "Total Force");
-    double us = pr_av(log, cr, sav, stot.data(), "Total Constr.");
-
-    double unb = (uf * fav + us * sav) / (fav + sav);
-    if (unb > 0)
-    {
-        unb = 10000.0 / unb;
-        fprintf(log, "\nTotal Scaling: %.0f%% of max performance\n\n", unb);
-    }
-}
index 847c6cdb9d345ca3d168c95e7cbf42d3c7cdc0c4..a21d13ed79d7348cb79e72642dce5eaf64f99563 100644 (file)
@@ -4,7 +4,7 @@
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
  * Copyright (c) 2013,2014,2015,2017,2018 by the GROMACS development team.
- * Copyright (c) 2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2019,2020,2021, 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.
 #ifndef GMX_GMXLIB_NRNB_H
 #define GMX_GMXLIB_NRNB_H
 
-#include <stdio.h>
-
-#include "gromacs/utility/basedefinitions.h"
-
-#define eNR_NBKERNEL_NONE (-1)
+#include <array>
+#include <cstdint>
+#include <cstdio>
 
 enum
 {
@@ -175,25 +173,33 @@ enum
 
 struct t_nrnb
 {
-    double n[eNRNB] = { 0 };
+    std::array<double, eNRNB> n = { 0 };
 };
 
-struct t_commrec;
-
 void clear_nrnb(t_nrnb* nrnb);
 
-void add_nrnb(t_nrnb* dest, t_nrnb* s1, t_nrnb* s2);
-
 void print_nrnb(FILE* out, t_nrnb* nrnb);
 
-void _inc_nrnb(t_nrnb* nrnb, int enr, int inc, char* file, int line);
-
-#ifdef DEBUG_NRNB
-#    define inc_nrnb(nrnb, enr, inc) _inc_nrnb(nrnb, enr, inc, __FILE__, __LINE__)
-#else
-#    define inc_nrnb(nrnb, enr, inc) (nrnb)->n[enr] += inc
-#endif
+/*! \brief Increment the nonbonded kernel flop counters
+ *
+ * \param nrnb The nonbonded kernel flop counters.
+ * \param index Which counter to incrememnt.
+ * \param increment How much to increment the counter by.
+ */
+static inline void inc_nrnb(t_nrnb* nrnb, int index, int increment)
+{
+    nrnb->n[index] += increment;
+}
 
+/*! \brief Atomic increment of nonbonded kernel flop counters
+ *
+ * \param nrnb The nonbonded kernel flop counters.
+ * \param index Which counter to incrememnt.
+ * \param increment How much to increment the counter by.
+ *
+ * Same as inc_nrnb but includes omp atomic pragma
+ */
+void atomicNrnbIncrement(t_nrnb* nrnb, int index, int increment);
 
 void print_flop(FILE* out, t_nrnb* nrnb, double* nbfs, double* mflop);
 /* Calculates the non-bonded forces and flop count.
@@ -203,9 +209,6 @@ void print_flop(FILE* out, t_nrnb* nrnb, double* nbfs, double* mflop);
 void print_perf(FILE* out, double nodetime, double realtime, int64_t nsteps, double delta_t, double nbfs, double mflop);
 /* Prints the performance, nbfs and mflop come from print_flop */
 
-void pr_load(FILE* log, struct t_commrec* cr, t_nrnb nrnb[]);
-/* Print detailed load balancing info */
-
 int cost_nrnb(int enr);
 /* Cost in i860 cycles of this component of MD */
 
index 681acd27b272af89f9cd7b427051232dd99b069e..488951c568a3aecc545cf149caee8e6831ff7af0 100644 (file)
@@ -639,7 +639,7 @@ static void finish_run(FILE*                     fplog,
         nrnbTotalStorage = std::make_unique<t_nrnb>();
         nrnb_tot         = nrnbTotalStorage.get();
 #if GMX_MPI
-        MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
+        MPI_Allreduce(nrnb->n.data(), nrnb_tot->n.data(), eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
 #endif
     }
     else