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