Move physics.* to math/units.*
[alexxy/gromacs.git] / src / gromacs / mdlib / perf_est.c
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, 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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42
43 #include "perf_est.h"
44 #include "gromacs/math/vec.h"
45 #include "mtop_util.h"
46 #include "types/commrec.h"
47 #include "nbnxn_search.h"
48 #include "nbnxn_consts.h"
49
50
51 /* Computational cost of bonded, non-bonded and PME calculations.
52  * This will be machine dependent.
53  * The numbers here are accurate for Intel Core2 and AMD Athlon 64
54  * in single precision. In double precision PME mesh is slightly cheaper,
55  * although not so much that the numbers need to be adjusted.
56  */
57
58 /* Cost of a pair interaction in the "group" cut-off scheme */
59 #define C_GR_FQ        1.5
60 #define C_GR_QLJ_CUT   1.5
61 #define C_GR_QLJ_TAB   2.0
62 #define C_GR_LJ_CUT    1.0
63 #define C_GR_LJ_TAB    1.75
64 /* Cost of 1 water with one Q/LJ atom */
65 #define C_GR_QLJW_CUT  2.0
66 #define C_GR_QLJW_TAB  2.25
67 /* Cost of 1 water with one Q atom or with 1/3 water (LJ negligible) */
68 #define C_GR_QW        1.75
69
70 /* Cost of a pair interaction in the "Verlet" cut-off scheme, QEXP is Ewald */
71 #define C_VT_LJ        0.30
72 #define C_VT_QRF_LJ    0.40
73 #define C_VT_QRF       0.30
74 #define C_VT_QEXP_LJ   0.55
75 #define C_VT_QEXP      0.50
76 /* Extra cost for expensive LJ interaction, e.g. pot-switch or LJ-PME */
77 #define C_VT_LJEXP_ADD 0.20
78
79 /* Cost of PME, with all components running with SSE instructions */
80 /* Cost of particle reordering and redistribution */
81 #define C_PME_REDIST  12.0
82 /* Cost of q spreading and force interpolation per charge (mainly memory) */
83 #define C_PME_SPREAD  0.30
84 /* Cost of fft's, will be multiplied with N log(N) */
85 #define C_PME_FFT     0.20
86 /* Cost of pme_solve, will be multiplied with N */
87 #define C_PME_SOLVE   0.50
88
89 /* Cost of a bonded interaction divided by the number of (pbc_)dx nrequired */
90 #define C_BOND        5.0
91
92 int n_bonded_dx(gmx_mtop_t *mtop, gmx_bool bExcl)
93 {
94     int            mb, nmol, ftype, ndxb, ndx_excl;
95     int            ndx;
96     gmx_moltype_t *molt;
97
98     /* Count the number of pbc_rvec_sub calls required for bonded interactions.
99      * This number is also roughly proportional to the computational cost.
100      */
101     ndx      = 0;
102     ndx_excl = 0;
103     for (mb = 0; mb < mtop->nmolblock; mb++)
104     {
105         molt = &mtop->moltype[mtop->molblock[mb].type];
106         nmol = mtop->molblock[mb].nmol;
107         for (ftype = 0; ftype < F_NRE; ftype++)
108         {
109             if (interaction_function[ftype].flags & IF_BOND)
110             {
111                 switch (ftype)
112                 {
113                     case F_POSRES:
114                     case F_FBPOSRES:  ndxb = 1; break;
115                     case F_CONNBONDS: ndxb = 0; break;
116                     default:     ndxb      = NRAL(ftype) - 1; break;
117                 }
118                 ndx += nmol*ndxb*molt->ilist[ftype].nr/(1 + NRAL(ftype));
119             }
120         }
121         if (bExcl)
122         {
123             ndx_excl += nmol*(molt->excls.nra - molt->atoms.nr)/2;
124         }
125         else
126         {
127             ndx_excl = 0;
128         }
129     }
130
131     if (debug)
132     {
133         fprintf(debug, "ndx bonded %d exclusions %d\n", ndx, ndx_excl);
134     }
135
136     ndx += ndx_excl;
137
138     return ndx;
139 }
140
141 static void pp_group_load(gmx_mtop_t *mtop, t_inputrec *ir, matrix box,
142                           int *nq_tot, int *nlj_tot,
143                           double *cost_pp,
144                           gmx_bool *bChargePerturbed, gmx_bool *bTypePerturbed)
145 {
146     t_atom        *atom;
147     int            mb, nmol, atnr, cg, a, a0, ncqlj, ncq, nclj;
148     gmx_bool       bBHAM, bLJcut, bWater, bQ, bLJ;
149     int            nw, nqlj, nq, nlj;
150     float          fq, fqlj, flj, fljtab, fqljw, fqw;
151     t_iparams     *iparams;
152     gmx_moltype_t *molt;
153
154     bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
155
156     bLJcut = ((ir->vdwtype == evdwCUT) && !bBHAM);
157
158     /* Computational cost of bonded, non-bonded and PME calculations.
159      * This will be machine dependent.
160      * The numbers here are accurate for Intel Core2 and AMD Athlon 64
161      * in single precision. In double precision PME mesh is slightly cheaper,
162      * although not so much that the numbers need to be adjusted.
163      */
164     fq    = C_GR_FQ;
165     fqlj  = (bLJcut ? C_GR_QLJ_CUT : C_GR_QLJ_TAB);
166     flj   = (bLJcut ? C_GR_LJ_CUT  : C_GR_LJ_TAB);
167     /* Cost of 1 water with one Q/LJ atom */
168     fqljw = (bLJcut ? C_GR_QLJW_CUT : C_GR_QLJW_TAB);
169     /* Cost of 1 water with one Q atom or with 1/3 water (LJ negligible) */
170     fqw   = C_GR_QW;
171
172     iparams           = mtop->ffparams.iparams;
173     atnr              = mtop->ffparams.atnr;
174     nw                = 0;
175     nqlj              = 0;
176     nq                = 0;
177     nlj               = 0;
178     *bChargePerturbed = FALSE;
179     for (mb = 0; mb < mtop->nmolblock; mb++)
180     {
181         molt = &mtop->moltype[mtop->molblock[mb].type];
182         atom = molt->atoms.atom;
183         nmol = mtop->molblock[mb].nmol;
184         a    = 0;
185         for (cg = 0; cg < molt->cgs.nr; cg++)
186         {
187             bWater = !bBHAM;
188             ncqlj  = 0;
189             ncq    = 0;
190             nclj   = 0;
191             a0     = a;
192             while (a < molt->cgs.index[cg+1])
193             {
194                 bQ  = (atom[a].q != 0 || atom[a].qB != 0);
195                 bLJ = (iparams[(atnr+1)*atom[a].type].lj.c6  != 0 ||
196                        iparams[(atnr+1)*atom[a].type].lj.c12 != 0);
197                 if (atom[a].q != atom[a].qB)
198                 {
199                     *bChargePerturbed = TRUE;
200                 }
201                 if (atom[a].type != atom[a].typeB)
202                 {
203                     *bTypePerturbed = TRUE;
204                 }
205                 /* This if this atom fits into water optimization */
206                 if (!((a == a0   &&  bQ &&  bLJ) ||
207                       (a == a0+1 &&  bQ && !bLJ) ||
208                       (a == a0+2 &&  bQ && !bLJ && atom[a].q == atom[a-1].q) ||
209                       (a == a0+3 && !bQ &&  bLJ)))
210                 {
211                     bWater = FALSE;
212                 }
213                 if (bQ && bLJ)
214                 {
215                     ncqlj++;
216                 }
217                 else
218                 {
219                     if (bQ)
220                     {
221                         ncq++;
222                     }
223                     if (bLJ)
224                     {
225                         nclj++;
226                     }
227                 }
228                 a++;
229             }
230             if (bWater)
231             {
232                 nw   += nmol;
233             }
234             else
235             {
236                 nqlj += nmol*ncqlj;
237                 nq   += nmol*ncq;
238                 nlj  += nmol*nclj;
239             }
240         }
241     }
242
243     *nq_tot  = nq  + nqlj + nw*3;
244     *nlj_tot = nlj + nqlj + nw;
245
246     if (debug)
247     {
248         fprintf(debug, "nw %d nqlj %d nq %d nlj %d\n", nw, nqlj, nq, nlj);
249     }
250
251     /* For the PP non-bonded cost it is (unrealistically) assumed
252      * that all atoms are distributed homogeneously in space.
253      * Factor 3 is used because a water molecule has 3 atoms
254      * (and TIP4P effectively has 3 interactions with (water) atoms)).
255      */
256     *cost_pp = 0.5*(fqljw*nw*nqlj +
257                     fqw  *nw*(3*nw + nq) +
258                     fqlj *nqlj*nqlj +
259                     fq   *nq*(3*nw + nqlj + nq) +
260                     flj  *nlj*(nw + nqlj + nlj))
261         *4/3*M_PI*ir->rlist*ir->rlist*ir->rlist/det(box);
262 }
263
264 static void pp_verlet_load(gmx_mtop_t *mtop, t_inputrec *ir, matrix box,
265                            int *nq_tot, int *nlj_tot,
266                            double *cost_pp,
267                            gmx_bool *bChargePerturbed, gmx_bool *bTypePerturbed)
268 {
269     t_atom        *atom;
270     int            mb, nmol, atnr, cg, a, a0, nqlj, nq, nlj;
271     gmx_bool       bQRF;
272     t_iparams     *iparams;
273     gmx_moltype_t *molt;
274     real           r_eff;
275     double         c_qlj, c_q, c_lj;
276     double         nat;
277     /* Conversion factor for reference vs SIMD kernel performance.
278      * The factor is about right for SSE2/4, but should be 2 higher for AVX256.
279      */
280 #ifdef GMX_DOUBLE
281     const real     nbnxn_refkernel_fac = 4.0;
282 #else
283     const real     nbnxn_refkernel_fac = 8.0;
284 #endif
285
286     bQRF = (EEL_RF(ir->coulombtype) || ir->coulombtype == eelCUT);
287
288     iparams           = mtop->ffparams.iparams;
289     atnr              = mtop->ffparams.atnr;
290     nqlj              = 0;
291     nq                = 0;
292     *bChargePerturbed = FALSE;
293     for (mb = 0; mb < mtop->nmolblock; mb++)
294     {
295         molt = &mtop->moltype[mtop->molblock[mb].type];
296         atom = molt->atoms.atom;
297         nmol = mtop->molblock[mb].nmol;
298         a    = 0;
299         for (a = 0; a < molt->atoms.nr; a++)
300         {
301             if (atom[a].q != 0 || atom[a].qB != 0)
302             {
303                 if (iparams[(atnr+1)*atom[a].type].lj.c6  != 0 ||
304                     iparams[(atnr+1)*atom[a].type].lj.c12 != 0)
305                 {
306                     nqlj += nmol;
307                 }
308                 else
309                 {
310                     nq += nmol;
311                 }
312             }
313             if (atom[a].q != atom[a].qB)
314             {
315                 *bChargePerturbed = TRUE;
316             }
317             if (atom[a].type != atom[a].typeB)
318             {
319                 *bTypePerturbed = TRUE;
320             }
321         }
322     }
323
324     nlj = mtop->natoms - nqlj - nq;
325
326     *nq_tot  = nqlj + nq;
327     *nlj_tot = nqlj + nlj;
328
329     /* Effective cut-off for cluster pair list of 4x4 atoms */
330     r_eff = ir->rlist + nbnxn_get_rlist_effective_inc(NBNXN_CPU_CLUSTER_I_SIZE, mtop->natoms/det(box));
331
332     if (debug)
333     {
334         fprintf(debug, "nqlj %d nq %d nlj %d rlist %.3f r_eff %.3f\n",
335                 nqlj, nq, nlj, ir->rlist, r_eff);
336     }
337
338     /* Determine the cost per pair interaction */
339     c_qlj = (bQRF ? C_VT_QRF_LJ : C_VT_QEXP_LJ);
340     c_q   = (bQRF ? C_VT_QRF    : C_VT_QEXP);
341     c_lj  = C_VT_LJ;
342     if (ir->vdw_modifier == eintmodPOTSWITCH || EVDW_PME(ir->vdwtype))
343     {
344         c_qlj += C_VT_LJEXP_ADD;
345         c_lj  += C_VT_LJEXP_ADD;
346     }
347     if (EVDW_PME(ir->vdwtype) && ir->ljpme_combination_rule == eljpmeLB)
348     {
349         /* We don't have LJ-PME LB comb. rule kernels, we use slow kernels */
350         c_qlj *= nbnxn_refkernel_fac;
351         c_q   *= nbnxn_refkernel_fac;
352         c_lj  *= nbnxn_refkernel_fac;
353     }
354
355     /* For the PP non-bonded cost it is (unrealistically) assumed
356      * that all atoms are distributed homogeneously in space.
357      */
358     /* Convert mtop->natoms to double to avoid int overflow */
359     nat      = mtop->natoms;
360     *cost_pp = 0.5*nat*(nqlj*c_qlj + nq*c_q + nlj*c_lj)
361         *4/3*M_PI*r_eff*r_eff*r_eff/det(box);
362 }
363
364 float pme_load_estimate(gmx_mtop_t *mtop, t_inputrec *ir, matrix box)
365 {
366     t_atom        *atom;
367     int            mb, nmol, atnr, cg, a, a0, nq_tot, nlj_tot, f;
368     gmx_bool       bBHAM, bLJcut, bChargePerturbed, bTypePerturbed;
369     gmx_bool       bWater, bQ, bLJ;
370     double         cost_bond, cost_pp, cost_redist, cost_spread, cost_fft, cost_solve, cost_pme;
371     float          ratio;
372     t_iparams     *iparams;
373     gmx_moltype_t *molt;
374
375     /* Computational cost of bonded, non-bonded and PME calculations.
376      * This will be machine dependent.
377      * The numbers here are accurate for Intel Core2 and AMD Athlon 64
378      * in single precision. In double precision PME mesh is slightly cheaper,
379      * although not so much that the numbers need to be adjusted.
380      */
381
382     iparams = mtop->ffparams.iparams;
383     atnr    = mtop->ffparams.atnr;
384
385     cost_bond = C_BOND*n_bonded_dx(mtop, TRUE);
386
387     if (ir->cutoff_scheme == ecutsGROUP)
388     {
389         pp_group_load(mtop, ir, box,
390                       &nq_tot, &nlj_tot, &cost_pp,
391                       &bChargePerturbed, &bTypePerturbed);
392     }
393     else
394     {
395         pp_verlet_load(mtop, ir, box,
396                        &nq_tot, &nlj_tot, &cost_pp,
397                        &bChargePerturbed, &bTypePerturbed);
398     }
399
400     cost_redist = 0;
401     cost_spread = 0;
402     cost_fft    = 0;
403     cost_solve  = 0;
404
405     if (EEL_PME(ir->coulombtype))
406     {
407         f            = ((ir->efep != efepNO && bChargePerturbed) ? 2 : 1);
408         cost_redist +=   C_PME_REDIST*nq_tot;
409         cost_spread += f*C_PME_SPREAD*nq_tot*pow(ir->pme_order, 3);
410         cost_fft    += f*C_PME_FFT*ir->nkx*ir->nky*ir->nkz*log(ir->nkx*ir->nky*ir->nkz);
411         cost_solve  += f*C_PME_SOLVE*ir->nkx*ir->nky*ir->nkz;
412     }
413
414     if (EVDW_PME(ir->vdwtype))
415     {
416         f            = ((ir->efep != efepNO && bTypePerturbed) ? 2 : 1);
417         if (ir->ljpme_combination_rule == eljpmeLB)
418         {
419             /* LB combination rule: we have 7 mesh terms */
420             f       *= 7;
421         }
422         cost_redist +=   C_PME_REDIST*nlj_tot;
423         cost_spread += f*C_PME_SPREAD*nlj_tot*pow(ir->pme_order, 3);
424         cost_fft    += f*C_PME_FFT*ir->nkx*ir->nky*ir->nkz*log(ir->nkx*ir->nky*ir->nkz);
425         cost_solve  += f*C_PME_SOLVE*ir->nkx*ir->nky*ir->nkz;
426     }
427
428     cost_pme = cost_redist + cost_spread + cost_fft + cost_solve;
429
430     ratio = cost_pme/(cost_bond + cost_pp + cost_pme);
431
432     if (debug)
433     {
434         fprintf(debug,
435                 "cost_bond   %f\n"
436                 "cost_pp     %f\n"
437                 "cost_redist %f\n"
438                 "cost_spread %f\n"
439                 "cost_fft    %f\n"
440                 "cost_solve  %f\n",
441                 cost_bond, cost_pp, cost_redist, cost_spread, cost_fft, cost_solve);
442
443         fprintf(debug, "Estimate for relative PME load: %.3f\n", ratio);
444     }
445
446     return ratio;
447 }