Merge branch release-5-1
[alexxy/gromacs.git] / src / gromacs / mdlib / perf_est.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2008, The GROMACS development team.
6  * Copyright (c) 2012,2014,2015, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "gromacs/legacyheaders/perf_est.h"
40
41 #include <cmath>
42
43 #include "gromacs/legacyheaders/types/commrec.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/mdlib/nbnxn_consts.h"
46 #include "gromacs/mdlib/nbnxn_search.h"
47 #include "gromacs/topology/topology.h"
48 #include "gromacs/utility/fatalerror.h"
49
50 /* Computational cost of bonded, non-bonded and PME calculations.
51  * This will be machine dependent.
52  * The numbers here are accurate for Intel Core2 and AMD Athlon 64
53  * in single precision. In double precision PME mesh is slightly cheaper,
54  * although not so much that the numbers need to be adjusted.
55  */
56
57 /* Cost of a pair interaction in the "group" cut-off scheme */
58 #define C_GR_FQ        1.5
59 #define C_GR_QLJ_CUT   1.5
60 #define C_GR_QLJ_TAB   2.0
61 #define C_GR_LJ_CUT    1.0
62 #define C_GR_LJ_TAB    1.75
63 /* Cost of 1 water with one Q/LJ atom */
64 #define C_GR_QLJW_CUT  2.0
65 #define C_GR_QLJW_TAB  2.25
66 /* Cost of 1 water with one Q atom or with 1/3 water (LJ negligible) */
67 #define C_GR_QW        1.75
68
69 /* Cost of a pair interaction in the "Verlet" cut-off scheme, QEXP is Ewald */
70 #define C_VT_LJ        0.30
71 #define C_VT_QRF_LJ    0.40
72 #define C_VT_QRF       0.30
73 #define C_VT_QEXP_LJ   0.55
74 #define C_VT_QEXP      0.50
75 /* Extra cost for expensive LJ interaction, e.g. pot-switch or LJ-PME */
76 #define C_VT_LJEXP_ADD 0.20
77
78 /* Cost of PME, with all components running with SSE instructions */
79 /* Cost of particle reordering and redistribution */
80 #define C_PME_REDIST  12.0
81 /* Cost of q spreading and force interpolation per charge (mainly memory) */
82 #define C_PME_SPREAD  0.30
83 /* Cost of fft's, will be multiplied with N log(N) */
84 #define C_PME_FFT     0.20
85 /* Cost of pme_solve, will be multiplied with N */
86 #define C_PME_SOLVE   0.50
87
88 /* Cost of a bonded interaction divided by the number of (pbc_)dx nrequired */
89 #define C_BOND        5.0
90
91 int n_bonded_dx(gmx_mtop_t *mtop, gmx_bool bExcl)
92 {
93     int            mb, nmol, ftype, ndxb, ndx_excl;
94     int            ndx;
95     gmx_moltype_t *molt;
96
97     /* Count the number of pbc_rvec_sub calls required for bonded interactions.
98      * This number is also roughly proportional to the computational cost.
99      */
100     ndx      = 0;
101     ndx_excl = 0;
102 #if defined _ICC && __ICC == 1400 || defined __ICL && __ICL == 1400
103 #pragma novector /* Work-around for incorrect vectorization */
104 #endif
105     for (mb = 0; mb < mtop->nmolblock; mb++)
106     {
107         molt = &mtop->moltype[mtop->molblock[mb].type];
108         nmol = mtop->molblock[mb].nmol;
109         for (ftype = 0; ftype < F_NRE; ftype++)
110         {
111             if (interaction_function[ftype].flags & IF_BOND)
112             {
113                 switch (ftype)
114                 {
115                     case F_POSRES:
116                     case F_FBPOSRES:  ndxb = 1; break;
117                     case F_CONNBONDS: ndxb = 0; break;
118                     default:     ndxb      = NRAL(ftype) - 1; break;
119                 }
120                 ndx += nmol*ndxb*molt->ilist[ftype].nr/(1 + NRAL(ftype));
121             }
122         }
123         if (bExcl)
124         {
125             ndx_excl += nmol*(molt->excls.nra - molt->atoms.nr)/2;
126         }
127         else
128         {
129             ndx_excl = 0;
130         }
131     }
132
133     if (debug)
134     {
135         fprintf(debug, "ndx bonded %d exclusions %d\n", ndx, ndx_excl);
136     }
137
138     ndx += ndx_excl;
139
140     return ndx;
141 }
142
143 static void pp_group_load(gmx_mtop_t *mtop, t_inputrec *ir, matrix box,
144                           int *nq_tot, int *nlj_tot,
145                           double *cost_pp,
146                           gmx_bool *bChargePerturbed, gmx_bool *bTypePerturbed)
147 {
148     t_atom        *atom;
149     int            mb, nmol, atnr, cg, a, a0, ncqlj, ncq, nclj;
150     gmx_bool       bBHAM, bLJcut, bWater, bQ, bLJ;
151     int            nw, nqlj, nq, nlj;
152     float          fq, fqlj, flj, fqljw, fqw;
153     t_iparams     *iparams;
154     gmx_moltype_t *molt;
155
156     bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
157
158     bLJcut = ((ir->vdwtype == evdwCUT) && !bBHAM);
159
160     /* Computational cost of bonded, non-bonded and PME calculations.
161      * This will be machine dependent.
162      * The numbers here are accurate for Intel Core2 and AMD Athlon 64
163      * in single precision. In double precision PME mesh is slightly cheaper,
164      * although not so much that the numbers need to be adjusted.
165      */
166     fq    = C_GR_FQ;
167     fqlj  = (bLJcut ? C_GR_QLJ_CUT : C_GR_QLJ_TAB);
168     flj   = (bLJcut ? C_GR_LJ_CUT  : C_GR_LJ_TAB);
169     /* Cost of 1 water with one Q/LJ atom */
170     fqljw = (bLJcut ? C_GR_QLJW_CUT : C_GR_QLJW_TAB);
171     /* Cost of 1 water with one Q atom or with 1/3 water (LJ negligible) */
172     fqw   = C_GR_QW;
173
174     iparams           = mtop->ffparams.iparams;
175     atnr              = mtop->ffparams.atnr;
176     nw                = 0;
177     nqlj              = 0;
178     nq                = 0;
179     nlj               = 0;
180     *bChargePerturbed = FALSE;
181     *bTypePerturbed   = FALSE;
182     for (mb = 0; mb < mtop->nmolblock; mb++)
183     {
184         molt = &mtop->moltype[mtop->molblock[mb].type];
185         atom = molt->atoms.atom;
186         nmol = mtop->molblock[mb].nmol;
187         a    = 0;
188         for (cg = 0; cg < molt->cgs.nr; cg++)
189         {
190             bWater = !bBHAM;
191             ncqlj  = 0;
192             ncq    = 0;
193             nclj   = 0;
194             a0     = a;
195             while (a < molt->cgs.index[cg+1])
196             {
197                 bQ  = (atom[a].q != 0 || atom[a].qB != 0);
198                 bLJ = (iparams[(atnr+1)*atom[a].type].lj.c6  != 0 ||
199                        iparams[(atnr+1)*atom[a].type].lj.c12 != 0);
200                 if (atom[a].q != atom[a].qB)
201                 {
202                     *bChargePerturbed = TRUE;
203                 }
204                 if (atom[a].type != atom[a].typeB)
205                 {
206                     *bTypePerturbed = TRUE;
207                 }
208                 /* This if this atom fits into water optimization */
209                 if (!((a == a0   &&  bQ &&  bLJ) ||
210                       (a == a0+1 &&  bQ && !bLJ) ||
211                       (a == a0+2 &&  bQ && !bLJ && atom[a].q == atom[a-1].q) ||
212                       (a == a0+3 && !bQ &&  bLJ)))
213                 {
214                     bWater = FALSE;
215                 }
216                 if (bQ && bLJ)
217                 {
218                     ncqlj++;
219                 }
220                 else
221                 {
222                     if (bQ)
223                     {
224                         ncq++;
225                     }
226                     if (bLJ)
227                     {
228                         nclj++;
229                     }
230                 }
231                 a++;
232             }
233             if (bWater)
234             {
235                 nw   += nmol;
236             }
237             else
238             {
239                 nqlj += nmol*ncqlj;
240                 nq   += nmol*ncq;
241                 nlj  += nmol*nclj;
242             }
243         }
244     }
245
246     *nq_tot  = nq  + nqlj + nw*3;
247     *nlj_tot = nlj + nqlj + nw;
248
249     if (debug)
250     {
251         fprintf(debug, "nw %d nqlj %d nq %d nlj %d\n", nw, nqlj, nq, nlj);
252     }
253
254     /* For the PP non-bonded cost it is (unrealistically) assumed
255      * that all atoms are distributed homogeneously in space.
256      * Factor 3 is used because a water molecule has 3 atoms
257      * (and TIP4P effectively has 3 interactions with (water) atoms)).
258      */
259     *cost_pp = 0.5*(fqljw*nw*nqlj +
260                     fqw  *nw*(3*nw + nq) +
261                     fqlj *nqlj*nqlj +
262                     fq   *nq*(3*nw + nqlj + nq) +
263                     flj  *nlj*(nw + nqlj + nlj))
264         *4/3*M_PI*ir->rlist*ir->rlist*ir->rlist/det(box);
265 }
266
267 static void pp_verlet_load(gmx_mtop_t *mtop, t_inputrec *ir, matrix box,
268                            int *nq_tot, int *nlj_tot,
269                            double *cost_pp,
270                            gmx_bool *bChargePerturbed, gmx_bool *bTypePerturbed)
271 {
272     t_atom        *atom;
273     int            mb, nmol, atnr, a, nqlj, nq, nlj;
274     gmx_bool       bQRF;
275     t_iparams     *iparams;
276     gmx_moltype_t *molt;
277     real           r_eff;
278     double         c_qlj, c_q, c_lj;
279     double         nat;
280     /* Conversion factor for reference vs SIMD kernel performance.
281      * The factor is about right for SSE2/4, but should be 2 higher for AVX256.
282      */
283 #ifdef GMX_DOUBLE
284     const real     nbnxn_refkernel_fac = 4.0;
285 #else
286     const real     nbnxn_refkernel_fac = 8.0;
287 #endif
288
289     bQRF = (EEL_RF(ir->coulombtype) || ir->coulombtype == eelCUT);
290
291     iparams           = mtop->ffparams.iparams;
292     atnr              = mtop->ffparams.atnr;
293     nqlj              = 0;
294     nq                = 0;
295     *bChargePerturbed = FALSE;
296     *bTypePerturbed   = FALSE;
297     for (mb = 0; mb < mtop->nmolblock; mb++)
298     {
299         molt = &mtop->moltype[mtop->molblock[mb].type];
300         atom = molt->atoms.atom;
301         nmol = mtop->molblock[mb].nmol;
302         for (a = 0; a < molt->atoms.nr; a++)
303         {
304             if (atom[a].q != 0 || atom[a].qB != 0)
305             {
306                 if (iparams[(atnr+1)*atom[a].type].lj.c6  != 0 ||
307                     iparams[(atnr+1)*atom[a].type].lj.c12 != 0)
308                 {
309                     nqlj += nmol;
310                 }
311                 else
312                 {
313                     nq += nmol;
314                 }
315             }
316             if (atom[a].q != atom[a].qB)
317             {
318                 *bChargePerturbed = TRUE;
319             }
320             if (atom[a].type != atom[a].typeB)
321             {
322                 *bTypePerturbed = TRUE;
323             }
324         }
325     }
326
327     nlj = mtop->natoms - nqlj - nq;
328
329     *nq_tot  = nqlj + nq;
330     *nlj_tot = nqlj + nlj;
331
332     /* Effective cut-off for cluster pair list of 4x4 atoms */
333     r_eff = ir->rlist + nbnxn_get_rlist_effective_inc(NBNXN_CPU_CLUSTER_I_SIZE, mtop->natoms/det(box));
334
335     if (debug)
336     {
337         fprintf(debug, "nqlj %d nq %d nlj %d rlist %.3f r_eff %.3f\n",
338                 nqlj, nq, nlj, ir->rlist, r_eff);
339     }
340
341     /* Determine the cost per pair interaction */
342     c_qlj = (bQRF ? C_VT_QRF_LJ : C_VT_QEXP_LJ);
343     c_q   = (bQRF ? C_VT_QRF    : C_VT_QEXP);
344     c_lj  = C_VT_LJ;
345     if (ir->vdw_modifier == eintmodPOTSWITCH || EVDW_PME(ir->vdwtype))
346     {
347         c_qlj += C_VT_LJEXP_ADD;
348         c_lj  += C_VT_LJEXP_ADD;
349     }
350     if (EVDW_PME(ir->vdwtype) && ir->ljpme_combination_rule == eljpmeLB)
351     {
352         /* We don't have LJ-PME LB comb. rule kernels, we use slow kernels */
353         c_qlj *= nbnxn_refkernel_fac;
354         c_q   *= nbnxn_refkernel_fac;
355         c_lj  *= nbnxn_refkernel_fac;
356     }
357
358     /* For the PP non-bonded cost it is (unrealistically) assumed
359      * that all atoms are distributed homogeneously in space.
360      */
361     /* Convert mtop->natoms to double to avoid int overflow */
362     nat      = mtop->natoms;
363     *cost_pp = 0.5*nat*(nqlj*c_qlj + nq*c_q + nlj*c_lj)
364         *4/3*M_PI*r_eff*r_eff*r_eff/det(box);
365 }
366
367 float pme_load_estimate(gmx_mtop_t *mtop, t_inputrec *ir, matrix box)
368 {
369     int            nq_tot, nlj_tot, f;
370     gmx_bool       bChargePerturbed, bTypePerturbed;
371     double         cost_bond, cost_pp, cost_redist, cost_spread, cost_fft, cost_solve, cost_pme;
372     float          ratio;
373
374     /* Computational cost of bonded, non-bonded and PME calculations.
375      * This will be machine dependent.
376      * The numbers here are accurate for Intel Core2 and AMD Athlon 64
377      * in single precision. In double precision PME mesh is slightly cheaper,
378      * although not so much that the numbers need to be adjusted.
379      */
380
381     cost_bond = C_BOND*n_bonded_dx(mtop, TRUE);
382
383     if (ir->cutoff_scheme == ecutsGROUP)
384     {
385         pp_group_load(mtop, ir, box,
386                       &nq_tot, &nlj_tot, &cost_pp,
387                       &bChargePerturbed, &bTypePerturbed);
388     }
389     else
390     {
391         pp_verlet_load(mtop, ir, box,
392                        &nq_tot, &nlj_tot, &cost_pp,
393                        &bChargePerturbed, &bTypePerturbed);
394     }
395
396     cost_redist = 0;
397     cost_spread = 0;
398     cost_fft    = 0;
399     cost_solve  = 0;
400
401     if (EEL_PME(ir->coulombtype))
402     {
403         f            = ((ir->efep != efepNO && bChargePerturbed) ? 2 : 1);
404         cost_redist +=   C_PME_REDIST*nq_tot;
405         cost_spread += f*C_PME_SPREAD*nq_tot*std::pow(static_cast<real>(ir->pme_order), static_cast<real>(3.0));
406         cost_fft    += f*C_PME_FFT*ir->nkx*ir->nky*ir->nkz*std::log(static_cast<real>(ir->nkx*ir->nky*ir->nkz));
407         cost_solve  += f*C_PME_SOLVE*ir->nkx*ir->nky*ir->nkz;
408     }
409
410     if (EVDW_PME(ir->vdwtype))
411     {
412         f            = ((ir->efep != efepNO && bTypePerturbed) ? 2 : 1);
413         if (ir->ljpme_combination_rule == eljpmeLB)
414         {
415             /* LB combination rule: we have 7 mesh terms */
416             f       *= 7;
417         }
418         cost_redist +=   C_PME_REDIST*nlj_tot;
419         cost_spread += f*C_PME_SPREAD*nlj_tot*std::pow(static_cast<real>(ir->pme_order), static_cast<real>(3.0));
420         cost_fft    += f*C_PME_FFT*ir->nkx*ir->nky*ir->nkz*std::log(static_cast<real>(ir->nkx*ir->nky*ir->nkz));
421         cost_solve  += f*C_PME_SOLVE*ir->nkx*ir->nky*ir->nkz;
422     }
423
424     cost_pme = cost_redist + cost_spread + cost_fft + cost_solve;
425
426     ratio = cost_pme/(cost_bond + cost_pp + cost_pme);
427
428     if (debug)
429     {
430         fprintf(debug,
431                 "cost_bond   %f\n"
432                 "cost_pp     %f\n"
433                 "cost_redist %f\n"
434                 "cost_spread %f\n"
435                 "cost_fft    %f\n"
436                 "cost_solve  %f\n",
437                 cost_bond, cost_pp, cost_redist, cost_spread, cost_fft, cost_solve);
438
439         fprintf(debug, "Estimate for relative PME load: %.3f\n", ratio);
440     }
441
442     return ratio;
443 }