Remove unnecessary config.h includes
[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 "gmxpre.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 "gromacs/legacyheaders/perf_est.h"
46 #include "gromacs/legacyheaders/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 #if __ICC == 1400 || __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, fljtab, 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     for (mb = 0; mb < mtop->nmolblock; mb++)
182     {
183         molt = &mtop->moltype[mtop->molblock[mb].type];
184         atom = molt->atoms.atom;
185         nmol = mtop->molblock[mb].nmol;
186         a    = 0;
187         for (cg = 0; cg < molt->cgs.nr; cg++)
188         {
189             bWater = !bBHAM;
190             ncqlj  = 0;
191             ncq    = 0;
192             nclj   = 0;
193             a0     = a;
194             while (a < molt->cgs.index[cg+1])
195             {
196                 bQ  = (atom[a].q != 0 || atom[a].qB != 0);
197                 bLJ = (iparams[(atnr+1)*atom[a].type].lj.c6  != 0 ||
198                        iparams[(atnr+1)*atom[a].type].lj.c12 != 0);
199                 if (atom[a].q != atom[a].qB)
200                 {
201                     *bChargePerturbed = TRUE;
202                 }
203                 if (atom[a].type != atom[a].typeB)
204                 {
205                     *bTypePerturbed = TRUE;
206                 }
207                 /* This if this atom fits into water optimization */
208                 if (!((a == a0   &&  bQ &&  bLJ) ||
209                       (a == a0+1 &&  bQ && !bLJ) ||
210                       (a == a0+2 &&  bQ && !bLJ && atom[a].q == atom[a-1].q) ||
211                       (a == a0+3 && !bQ &&  bLJ)))
212                 {
213                     bWater = FALSE;
214                 }
215                 if (bQ && bLJ)
216                 {
217                     ncqlj++;
218                 }
219                 else
220                 {
221                     if (bQ)
222                     {
223                         ncq++;
224                     }
225                     if (bLJ)
226                     {
227                         nclj++;
228                     }
229                 }
230                 a++;
231             }
232             if (bWater)
233             {
234                 nw   += nmol;
235             }
236             else
237             {
238                 nqlj += nmol*ncqlj;
239                 nq   += nmol*ncq;
240                 nlj  += nmol*nclj;
241             }
242         }
243     }
244
245     *nq_tot  = nq  + nqlj + nw*3;
246     *nlj_tot = nlj + nqlj + nw;
247
248     if (debug)
249     {
250         fprintf(debug, "nw %d nqlj %d nq %d nlj %d\n", nw, nqlj, nq, nlj);
251     }
252
253     /* For the PP non-bonded cost it is (unrealistically) assumed
254      * that all atoms are distributed homogeneously in space.
255      * Factor 3 is used because a water molecule has 3 atoms
256      * (and TIP4P effectively has 3 interactions with (water) atoms)).
257      */
258     *cost_pp = 0.5*(fqljw*nw*nqlj +
259                     fqw  *nw*(3*nw + nq) +
260                     fqlj *nqlj*nqlj +
261                     fq   *nq*(3*nw + nqlj + nq) +
262                     flj  *nlj*(nw + nqlj + nlj))
263         *4/3*M_PI*ir->rlist*ir->rlist*ir->rlist/det(box);
264 }
265
266 static void pp_verlet_load(gmx_mtop_t *mtop, t_inputrec *ir, matrix box,
267                            int *nq_tot, int *nlj_tot,
268                            double *cost_pp,
269                            gmx_bool *bChargePerturbed, gmx_bool *bTypePerturbed)
270 {
271     t_atom        *atom;
272     int            mb, nmol, atnr, cg, a, a0, nqlj, nq, nlj;
273     gmx_bool       bQRF;
274     t_iparams     *iparams;
275     gmx_moltype_t *molt;
276     real           r_eff;
277     double         c_qlj, c_q, c_lj;
278     double         nat;
279     /* Conversion factor for reference vs SIMD kernel performance.
280      * The factor is about right for SSE2/4, but should be 2 higher for AVX256.
281      */
282 #ifdef GMX_DOUBLE
283     const real     nbnxn_refkernel_fac = 4.0;
284 #else
285     const real     nbnxn_refkernel_fac = 8.0;
286 #endif
287
288     bQRF = (EEL_RF(ir->coulombtype) || ir->coulombtype == eelCUT);
289
290     iparams           = mtop->ffparams.iparams;
291     atnr              = mtop->ffparams.atnr;
292     nqlj              = 0;
293     nq                = 0;
294     *bChargePerturbed = FALSE;
295     for (mb = 0; mb < mtop->nmolblock; mb++)
296     {
297         molt = &mtop->moltype[mtop->molblock[mb].type];
298         atom = molt->atoms.atom;
299         nmol = mtop->molblock[mb].nmol;
300         a    = 0;
301         for (a = 0; a < molt->atoms.nr; a++)
302         {
303             if (atom[a].q != 0 || atom[a].qB != 0)
304             {
305                 if (iparams[(atnr+1)*atom[a].type].lj.c6  != 0 ||
306                     iparams[(atnr+1)*atom[a].type].lj.c12 != 0)
307                 {
308                     nqlj += nmol;
309                 }
310                 else
311                 {
312                     nq += nmol;
313                 }
314             }
315             if (atom[a].q != atom[a].qB)
316             {
317                 *bChargePerturbed = TRUE;
318             }
319             if (atom[a].type != atom[a].typeB)
320             {
321                 *bTypePerturbed = TRUE;
322             }
323         }
324     }
325
326     nlj = mtop->natoms - nqlj - nq;
327
328     *nq_tot  = nqlj + nq;
329     *nlj_tot = nqlj + nlj;
330
331     /* Effective cut-off for cluster pair list of 4x4 atoms */
332     r_eff = ir->rlist + nbnxn_get_rlist_effective_inc(NBNXN_CPU_CLUSTER_I_SIZE, mtop->natoms/det(box));
333
334     if (debug)
335     {
336         fprintf(debug, "nqlj %d nq %d nlj %d rlist %.3f r_eff %.3f\n",
337                 nqlj, nq, nlj, ir->rlist, r_eff);
338     }
339
340     /* Determine the cost per pair interaction */
341     c_qlj = (bQRF ? C_VT_QRF_LJ : C_VT_QEXP_LJ);
342     c_q   = (bQRF ? C_VT_QRF    : C_VT_QEXP);
343     c_lj  = C_VT_LJ;
344     if (ir->vdw_modifier == eintmodPOTSWITCH || EVDW_PME(ir->vdwtype))
345     {
346         c_qlj += C_VT_LJEXP_ADD;
347         c_lj  += C_VT_LJEXP_ADD;
348     }
349     if (EVDW_PME(ir->vdwtype) && ir->ljpme_combination_rule == eljpmeLB)
350     {
351         /* We don't have LJ-PME LB comb. rule kernels, we use slow kernels */
352         c_qlj *= nbnxn_refkernel_fac;
353         c_q   *= nbnxn_refkernel_fac;
354         c_lj  *= nbnxn_refkernel_fac;
355     }
356
357     /* For the PP non-bonded cost it is (unrealistically) assumed
358      * that all atoms are distributed homogeneously in space.
359      */
360     /* Convert mtop->natoms to double to avoid int overflow */
361     nat      = mtop->natoms;
362     *cost_pp = 0.5*nat*(nqlj*c_qlj + nq*c_q + nlj*c_lj)
363         *4/3*M_PI*r_eff*r_eff*r_eff/det(box);
364 }
365
366 float pme_load_estimate(gmx_mtop_t *mtop, t_inputrec *ir, matrix box)
367 {
368     t_atom        *atom;
369     int            mb, nmol, atnr, cg, a, a0, nq_tot, nlj_tot, f;
370     gmx_bool       bBHAM, bLJcut, bChargePerturbed, bTypePerturbed;
371     gmx_bool       bWater, bQ, bLJ;
372     double         cost_bond, cost_pp, cost_redist, cost_spread, cost_fft, cost_solve, cost_pme;
373     float          ratio;
374     t_iparams     *iparams;
375     gmx_moltype_t *molt;
376
377     /* Computational cost of bonded, non-bonded and PME calculations.
378      * This will be machine dependent.
379      * The numbers here are accurate for Intel Core2 and AMD Athlon 64
380      * in single precision. In double precision PME mesh is slightly cheaper,
381      * although not so much that the numbers need to be adjusted.
382      */
383
384     iparams = mtop->ffparams.iparams;
385     atnr    = mtop->ffparams.atnr;
386
387     cost_bond = C_BOND*n_bonded_dx(mtop, TRUE);
388
389     if (ir->cutoff_scheme == ecutsGROUP)
390     {
391         pp_group_load(mtop, ir, box,
392                       &nq_tot, &nlj_tot, &cost_pp,
393                       &bChargePerturbed, &bTypePerturbed);
394     }
395     else
396     {
397         pp_verlet_load(mtop, ir, box,
398                        &nq_tot, &nlj_tot, &cost_pp,
399                        &bChargePerturbed, &bTypePerturbed);
400     }
401
402     cost_redist = 0;
403     cost_spread = 0;
404     cost_fft    = 0;
405     cost_solve  = 0;
406
407     if (EEL_PME(ir->coulombtype))
408     {
409         f            = ((ir->efep != efepNO && bChargePerturbed) ? 2 : 1);
410         cost_redist +=   C_PME_REDIST*nq_tot;
411         cost_spread += f*C_PME_SPREAD*nq_tot*pow(ir->pme_order, 3);
412         cost_fft    += f*C_PME_FFT*ir->nkx*ir->nky*ir->nkz*log(ir->nkx*ir->nky*ir->nkz);
413         cost_solve  += f*C_PME_SOLVE*ir->nkx*ir->nky*ir->nkz;
414     }
415
416     if (EVDW_PME(ir->vdwtype))
417     {
418         f            = ((ir->efep != efepNO && bTypePerturbed) ? 2 : 1);
419         if (ir->ljpme_combination_rule == eljpmeLB)
420         {
421             /* LB combination rule: we have 7 mesh terms */
422             f       *= 7;
423         }
424         cost_redist +=   C_PME_REDIST*nlj_tot;
425         cost_spread += f*C_PME_SPREAD*nlj_tot*pow(ir->pme_order, 3);
426         cost_fft    += f*C_PME_FFT*ir->nkx*ir->nky*ir->nkz*log(ir->nkx*ir->nky*ir->nkz);
427         cost_solve  += f*C_PME_SOLVE*ir->nkx*ir->nky*ir->nkz;
428     }
429
430     cost_pme = cost_redist + cost_spread + cost_fft + cost_solve;
431
432     ratio = cost_pme/(cost_bond + cost_pp + cost_pme);
433
434     if (debug)
435     {
436         fprintf(debug,
437                 "cost_bond   %f\n"
438                 "cost_pp     %f\n"
439                 "cost_redist %f\n"
440                 "cost_spread %f\n"
441                 "cost_fft    %f\n"
442                 "cost_solve  %f\n",
443                 cost_bond, cost_pp, cost_redist, cost_spread, cost_fft, cost_solve);
444
445         fprintf(debug, "Estimate for relative PME load: %.3f\n", ratio);
446     }
447
448     return ratio;
449 }