3f30f596fe014513c413239f61154a5ac8d9685b
[alexxy/gromacs.git] / src / mdlib / perf_est.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  *                This source code is part of
5  * 
6  *                 G   R   O   M   A   C   S
7  * 
8  *          GROningen MAchine for Chemical Simulations
9  * 
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2008, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14  
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  * 
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  * 
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <math.h>
40
41 #include "perf_est.h"
42 #include "physics.h"
43 #include "vec.h"
44 #include "mtop_util.h"
45 #include "types/commrec.h"
46 #include "nbnxn_search.h"
47 #include "nbnxn_consts.h"
48
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" */
70 #define C_VT_LJ       0.30
71 #define C_VT_QLJ_RF   0.40
72 #define C_VT_Q_RF     0.30
73 #define C_VT_QLJ_TAB  0.55
74 #define C_VT_Q_TAB    0.50
75
76 /* Cost of PME, with all components running with SSE instructions */
77 /* Cost of particle reordering and redistribution */
78 #define C_PME_REDIST  12.0
79 /* Cost of q spreading and force interpolation per charge (mainly memory) */
80 #define C_PME_SPREAD  0.30
81 /* Cost of fft's, will be multiplied with N log(N) */
82 #define C_PME_FFT     0.20
83 /* Cost of pme_solve, will be multiplied with N */
84 #define C_PME_SOLVE   0.50
85
86 /* Cost of a bonded interaction divided by the number of (pbc_)dx nrequired */
87 #define C_BOND        5.0
88
89 int n_bonded_dx(gmx_mtop_t *mtop,gmx_bool bExcl)
90 {
91   int mb,nmol,ftype,ndxb,ndx_excl;
92   int ndx;
93   gmx_moltype_t *molt;
94
95   /* Count the number of pbc_rvec_sub calls required for bonded interactions.
96    * This number is also roughly proportional to the computational cost.
97    */
98   ndx = 0;
99   ndx_excl = 0;
100   for(mb=0; mb<mtop->nmolblock; mb++) {
101     molt = &mtop->moltype[mtop->molblock[mb].type];
102     nmol = mtop->molblock[mb].nmol;
103     for(ftype=0; ftype<F_NRE; ftype++) {
104       if (interaction_function[ftype].flags & IF_BOND) {
105         switch (ftype) {
106                 case F_POSRES:    ndxb = 1; break;
107                 case F_CONNBONDS: ndxb = 0; break;
108                 default:     ndxb = NRAL(ftype) - 1; break;
109         }
110         ndx += nmol*ndxb*molt->ilist[ftype].nr/(1 + NRAL(ftype));
111       }
112     }
113     if (bExcl) {
114       ndx_excl += nmol*(molt->excls.nra - molt->atoms.nr)/2;
115     } else {
116       ndx_excl = 0;
117     }
118   }
119
120   if (debug)
121     fprintf(debug,"ndx bonded %d exclusions %d\n",ndx,ndx_excl);
122
123   ndx += ndx_excl;
124
125   return ndx;
126 }
127
128 static void pp_group_load(gmx_mtop_t *mtop,t_inputrec *ir,matrix box,
129                           int *nq_tot,
130                           double *cost_pp,
131                           gmx_bool *bChargePerturbed)
132 {
133     t_atom *atom;
134     int  mb,nmol,atnr,cg,a,a0,ncqlj,ncq,nclj;
135     gmx_bool bBHAM,bLJcut,bWater,bQ,bLJ;
136     int nw,nqlj,nq,nlj;
137     float fq,fqlj,flj,fljtab,fqljw,fqw;
138     t_iparams *iparams;
139     gmx_moltype_t *molt;
140
141     bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
142
143     bLJcut = ((ir->vdwtype == evdwCUT) && !bBHAM);
144
145     /* Computational cost of bonded, non-bonded and PME calculations.
146      * This will be machine dependent.
147      * The numbers here are accurate for Intel Core2 and AMD Athlon 64
148      * in single precision. In double precision PME mesh is slightly cheaper,
149      * although not so much that the numbers need to be adjusted.
150      */
151     fq    = C_GR_FQ;
152     fqlj  = (bLJcut ? C_GR_QLJ_CUT : C_GR_QLJ_TAB);
153     flj   = (bLJcut ? C_GR_LJ_CUT  : C_GR_LJ_TAB);
154     /* Cost of 1 water with one Q/LJ atom */
155     fqljw = (bLJcut ? C_GR_QLJW_CUT : C_GR_QLJW_TAB);
156     /* Cost of 1 water with one Q atom or with 1/3 water (LJ negligible) */
157     fqw   = C_GR_QW;
158
159     iparams = mtop->ffparams.iparams;
160     atnr = mtop->ffparams.atnr;
161     nw   = 0;
162     nqlj = 0;
163     nq   = 0;
164     nlj  = 0;
165     *bChargePerturbed = FALSE;
166     for(mb=0; mb<mtop->nmolblock; mb++)
167         {
168         molt = &mtop->moltype[mtop->molblock[mb].type];
169         atom = molt->atoms.atom;
170         nmol = mtop->molblock[mb].nmol;
171         a = 0;
172         for(cg=0; cg<molt->cgs.nr; cg++)
173         {
174             bWater = !bBHAM;
175             ncqlj = 0;
176             ncq   = 0;
177             nclj  = 0;
178             a0    = a;
179             while (a < molt->cgs.index[cg+1])
180             {
181                 bQ  = (atom[a].q != 0 || atom[a].qB != 0);
182                 bLJ = (iparams[(atnr+1)*atom[a].type].lj.c6  != 0 ||
183                        iparams[(atnr+1)*atom[a].type].lj.c12 != 0);
184                 if (atom[a].q != atom[a].qB)
185                 {
186                     *bChargePerturbed = TRUE;
187                 }
188                 /* This if this atom fits into water optimization */
189                 if (!((a == a0   &&  bQ &&  bLJ) ||
190                       (a == a0+1 &&  bQ && !bLJ) ||
191                       (a == a0+2 &&  bQ && !bLJ && atom[a].q == atom[a-1].q) ||
192                       (a == a0+3 && !bQ &&  bLJ)))
193                     bWater = FALSE;
194                 if (bQ && bLJ)
195                 {
196                     ncqlj++;
197                 }
198                 else
199                 {
200                     if (bQ)
201                     {
202                         ncq++;
203                     }
204                     if (bLJ)
205                     {
206                         nclj++;
207                     }
208                 }
209                 a++;
210             }
211             if (bWater)
212             {
213                 nw   += nmol;
214             }
215             else
216             {
217                 nqlj += nmol*ncqlj;
218                 nq   += nmol*ncq;
219                 nlj  += nmol*nclj;
220             }
221         }
222     }
223
224     *nq_tot = nq + nqlj + nw*3;
225
226     if (debug)
227     {
228       fprintf(debug,"nw %d nqlj %d nq %d nlj %d\n",nw,nqlj,nq,nlj);
229     }
230
231     /* For the PP non-bonded cost it is (unrealistically) assumed
232      * that all atoms are distributed homogeneously in space.
233      * Factor 3 is used because a water molecule has 3 atoms
234      * (and TIP4P effectively has 3 interactions with (water) atoms)).
235      */
236     *cost_pp = 0.5*(fqljw*nw*nqlj +
237                     fqw  *nw*(3*nw + nq) +
238                     fqlj *nqlj*nqlj +
239                     fq   *nq*(3*nw + nqlj + nq) +
240                     flj  *nlj*(nw + nqlj + nlj))
241         *4/3*M_PI*ir->rlist*ir->rlist*ir->rlist/det(box);
242 }
243
244 static void pp_verlet_load(gmx_mtop_t *mtop,t_inputrec *ir,matrix box,
245                            int *nq_tot,
246                            double *cost_pp,
247                            gmx_bool *bChargePerturbed)
248 {
249     t_atom *atom;
250     int  mb,nmol,atnr,cg,a,a0,nqlj,nq,nlj;
251     gmx_bool bQRF;
252     t_iparams *iparams;
253     gmx_moltype_t *molt;
254     float r_eff;
255     double nat;
256
257     bQRF = (EEL_RF(ir->coulombtype) || ir->coulombtype == eelCUT);
258
259     iparams = mtop->ffparams.iparams;
260     atnr = mtop->ffparams.atnr;
261     nqlj = 0;
262     nq   = 0;
263     *bChargePerturbed = FALSE;
264     for(mb=0; mb<mtop->nmolblock; mb++)
265         {
266         molt = &mtop->moltype[mtop->molblock[mb].type];
267         atom = molt->atoms.atom;
268         nmol = mtop->molblock[mb].nmol;
269         a = 0;
270         for(a=0; a<molt->atoms.nr; a++)
271         {
272             if (atom[a].q != 0 || atom[a].qB != 0)
273             {
274                 if (iparams[(atnr+1)*atom[a].type].lj.c6  != 0 ||
275                     iparams[(atnr+1)*atom[a].type].lj.c12 != 0)
276                 {
277                     nqlj += nmol;
278                 }
279                 else
280                 {
281                     nq += nmol;
282                 }
283             }
284             if (atom[a].q != atom[a].qB)
285             {
286                 *bChargePerturbed = TRUE;
287             }
288         }
289     }
290
291     nlj = mtop->natoms - nqlj - nq;
292
293     *nq_tot = nqlj + nq;
294
295     /* Effective cut-off for cluster pair list of 4x4 atoms */
296     r_eff = ir->rlist + nbnxn_get_rlist_effective_inc(NBNXN_CPU_CLUSTER_I_SIZE,mtop->natoms/det(box));
297
298     if (debug)
299     {
300         fprintf(debug,"nqlj %d nq %d nlj %d rlist %.3f r_eff %.3f\n",
301                 nqlj,nq,nlj,ir->rlist,r_eff);
302     }
303
304     /* For the PP non-bonded cost it is (unrealistically) assumed
305      * that all atoms are distributed homogeneously in space.
306      */
307     /* Convert mtop->natoms to double to avoid int overflow */
308     nat = mtop->natoms;
309     *cost_pp = 0.5*(nqlj*nat*(bQRF ? C_VT_QLJ_RF : C_VT_QLJ_TAB) +
310                     nq*nat*(bQRF ? C_VT_Q_RF : C_VT_Q_TAB) +
311                     nlj*nat*C_VT_LJ)
312         *4/3*M_PI*r_eff*r_eff*r_eff/det(box);
313 }
314
315 float pme_load_estimate(gmx_mtop_t *mtop,t_inputrec *ir,matrix box)
316 {
317   t_atom *atom;
318   int  mb,nmol,atnr,cg,a,a0,nq_tot;
319   gmx_bool bBHAM,bLJcut,bChargePerturbed,bWater,bQ,bLJ;
320   double cost_bond,cost_pp,cost_redist,cost_spread,cost_fft,cost_solve,cost_pme;
321   float ratio;
322   t_iparams *iparams;
323   gmx_moltype_t *molt;
324
325   /* Computational cost of bonded, non-bonded and PME calculations.
326    * This will be machine dependent.
327    * The numbers here are accurate for Intel Core2 and AMD Athlon 64
328    * in single precision. In double precision PME mesh is slightly cheaper,
329    * although not so much that the numbers need to be adjusted.
330    */
331
332   iparams = mtop->ffparams.iparams;
333   atnr = mtop->ffparams.atnr;
334
335   cost_bond = C_BOND*n_bonded_dx(mtop,TRUE);
336
337   if (ir->cutoff_scheme == ecutsGROUP)
338   {
339       pp_group_load(mtop,ir,box,&nq_tot,&cost_pp,&bChargePerturbed);
340   }
341   else
342   {
343       pp_verlet_load(mtop,ir,box,&nq_tot,&cost_pp,&bChargePerturbed);
344   }
345   
346   cost_redist = C_PME_REDIST*nq_tot;
347   cost_spread = C_PME_SPREAD*nq_tot*pow(ir->pme_order,3);
348   cost_fft    = C_PME_FFT*ir->nkx*ir->nky*ir->nkz*log(ir->nkx*ir->nky*ir->nkz);
349   cost_solve  = C_PME_SOLVE*ir->nkx*ir->nky*ir->nkz;
350
351   if (ir->efep != efepNO && bChargePerturbed) {
352     /* All PME work, except redist & spline coefficient calculation, doubles */
353     cost_spread *= 2;
354     cost_fft    *= 2;
355     cost_solve  *= 2;
356   }
357
358   cost_pme = cost_redist + cost_spread + cost_fft + cost_solve;
359
360   ratio = cost_pme/(cost_bond + cost_pp + cost_pme);
361
362   if (debug) {
363     fprintf(debug,
364             "cost_bond   %f\n"
365             "cost_pp     %f\n"
366             "cost_redist %f\n"
367             "cost_spread %f\n"
368             "cost_fft    %f\n"
369             "cost_solve  %f\n",
370             cost_bond,cost_pp,cost_redist,cost_spread,cost_fft,cost_solve);
371
372     fprintf(debug,"Estimate for relative PME load: %.3f\n",ratio);
373   }
374
375   return ratio;
376 }