Re-organize BlueGene toolchain files
[gromacs.git] / src / gmxlib / nrnb.c
blob1fd315c14a2e6601fcfd965790188ced5086bcf8
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
42 #include <string.h>
43 #include "types/commrec.h"
44 #include "sysstuff.h"
45 #include "gmx_fatal.h"
46 #include "names.h"
47 #include "macros.h"
48 #include "nrnb.h"
49 #include "main.h"
50 #include "smalloc.h"
51 #include "copyrite.h"
57 typedef struct {
58 const char *name;
59 int flop;
60 } t_nrnb_data;
63 static const t_nrnb_data nbdata[eNRNB] = {
64 /* These are re-used for different NB kernels, since there are so many.
65 * The actual number of flops is set dynamically.
67 { "NB VdW [V&F]", 1 },
68 { "NB VdW [F]", 1 },
69 { "NB Elec. [V&F]", 1 },
70 { "NB Elec. [F]", 1 },
71 { "NB Elec. [W3,V&F]", 1 },
72 { "NB Elec. [W3,F]", 1 },
73 { "NB Elec. [W3-W3,V&F]", 1 },
74 { "NB Elec. [W3-W3,F]", 1 },
75 { "NB Elec. [W4,V&F]", 1 },
76 { "NB Elec. [W4,F]", 1 },
77 { "NB Elec. [W4-W4,V&F]", 1 },
78 { "NB Elec. [W4-W4,F]", 1 },
79 { "NB VdW & Elec. [V&F]", 1 },
80 { "NB VdW & Elec. [F]", 1 },
81 { "NB VdW & Elec. [W3,V&F]", 1 },
82 { "NB VdW & Elec. [W3,F]", 1 },
83 { "NB VdW & Elec. [W3-W3,V&F]", 1 },
84 { "NB VdW & Elec. [W3-W3,F]", 1 },
85 { "NB VdW & Elec. [W4,V&F]", 1 },
86 { "NB VdW & Elec. [W4,F]", 1 },
87 { "NB VdW & Elec. [W4-W4,V&F]", 1 },
88 { "NB VdW & Elec. [W4-W4,F]", 1 },
90 { "NB Generic kernel", 1 },
91 { "NB Free energy kernel", 1 },
92 { "NB All-vs-all", 1 },
93 { "NB All-vs-all, GB", 1 },
95 { "Pair Search distance check", 9 }, /* nbnxn pair dist. check */
96 /* nbnxn kernel flops are based on inner-loops without exclusion checks.
97 * Plain Coulomb runs through the RF kernels, except with CUDA.
98 * invsqrt is counted as 6 flops: 1 for _mm_rsqt_ps + 5 for iteration.
99 * The flops are equal for plain-C, x86 SIMD and CUDA, except for:
100 * - plain-C kernel uses one flop more for Coulomb-only (F) than listed
101 * - x86 SIMD LJ geom-comb.rule kernels (fastest) use 2 more flops
102 * - x86 SIMD LJ LB-comb.rule kernels (fast) use 3 (8 for F+E) more flops
103 * - GPU always does exclusions, which requires 2-4 flops, but as invsqrt
104 * is always counted as 6 flops, this roughly compensates.
106 { "NxN RF Elec. + VdW [F]", 38 }, /* nbnxn kernel LJ+RF, no ener */
107 { "NxN RF Elec. + VdW [V&F]", 54 },
108 { "NxN QSTab Elec. + VdW [F]", 41 }, /* nbnxn kernel LJ+tab, no en */
109 { "NxN QSTab Elec. + VdW [V&F]", 59 },
110 { "NxN Ewald Elec. + VdW [F]", 66 }, /* nbnxn kernel LJ+Ewald, no en */
111 { "NxN Ewald Elec. + VdW [V&F]", 107 },
112 { "NxN VdW [F]", 33 }, /* nbnxn kernel LJ, no ener */
113 { "NxN VdW [V&F]", 43 },
114 { "NxN RF Electrostatics [F]", 31 }, /* nbnxn kernel RF, no ener */
115 { "NxN RF Electrostatics [V&F]", 36 },
116 { "NxN QSTab Elec. [F]", 34 }, /* nbnxn kernel tab, no ener */
117 { "NxN QSTab Elec. [V&F]", 41 },
118 { "NxN Ewald Elec. [F]", 61 }, /* nbnxn kernel Ewald, no ener */
119 { "NxN Ewald Elec. [V&F]", 84 },
120 { "1,4 nonbonded interactions", 90 },
121 { "Born radii (Still)", 47 },
122 { "Born radii (HCT/OBC)", 183 },
123 { "Born force chain rule", 15 },
124 { "All-vs-All Still radii", 1 },
125 { "All-vs-All HCT/OBC radii", 1 },
126 { "All-vs-All Born chain rule", 1 },
127 { "Calc Weights", 36 },
128 { "Spread Q", 6 },
129 { "Spread Q Bspline", 2 },
130 { "Gather F", 23 },
131 { "Gather F Bspline", 6 },
132 { "3D-FFT", 8 },
133 { "Convolution", 4 },
134 { "Solve PME", 64 },
135 { "NS-Pairs", 21 },
136 { "Reset In Box", 3 },
137 { "Shift-X", 6 },
138 { "CG-CoM", 3 },
139 { "Sum Forces", 1 },
140 { "Bonds", 59 },
141 { "G96Bonds", 44 },
142 { "FENE Bonds", 58 },
143 { "Tab. Bonds", 62 },
144 { "Restraint Potential", 86 },
145 { "Linear Angles", 57 },
146 { "Angles", 168 },
147 { "G96Angles", 150 },
148 { "Quartic Angles", 160 },
149 { "Tab. Angles", 169 },
150 { "Propers", 229 },
151 { "Impropers", 208 },
152 { "RB-Dihedrals", 247 },
153 { "Four. Dihedrals", 247 },
154 { "Tab. Dihedrals", 227 },
155 { "Dist. Restr.", 200 },
156 { "Orient. Restr.", 200 },
157 { "Dihedral Restr.", 200 },
158 { "Pos. Restr.", 50 },
159 { "Angle Restr.", 191 },
160 { "Angle Restr. Z", 164 },
161 { "Morse Potent.", 83 },
162 { "Cubic Bonds", 54 },
163 { "Walls", 31 },
164 { "Polarization", 59 },
165 { "Anharmonic Polarization", 72 },
166 { "Water Pol.", 62 },
167 { "Thole Pol.", 296 },
168 { "Virial", 18 },
169 { "Update", 31 },
170 { "Ext.ens. Update", 54 },
171 { "Stop-CM", 10 },
172 { "P-Coupling", 6 },
173 { "Calc-Ekin", 27 },
174 { "Lincs", 60 },
175 { "Lincs-Mat", 4 },
176 { "Shake", 30 },
177 { "Constraint-V", 8 },
178 { "Shake-Init", 10 },
179 { "Constraint-Vir", 24 },
180 { "Settle", 323 },
181 { "Virtual Site 2", 23 },
182 { "Virtual Site 3", 37 },
183 { "Virtual Site 3fd", 95 },
184 { "Virtual Site 3fad", 176 },
185 { "Virtual Site 3out", 87 },
186 { "Virtual Site 4fd", 110 },
187 { "Virtual Site 4fdn", 254 },
188 { "Virtual Site N", 15 },
189 { "Mixed Generalized Born stuff", 10 }
193 void init_nrnb(t_nrnb *nrnb)
195 int i;
197 for(i=0; (i<eNRNB); i++)
198 nrnb->n[i]=0.0;
201 void cp_nrnb(t_nrnb *dest, t_nrnb *src)
203 int i;
205 for(i=0; (i<eNRNB); i++)
206 dest->n[i]=src->n[i];
209 void add_nrnb(t_nrnb *dest, t_nrnb *s1, t_nrnb *s2)
211 int i;
213 for(i=0; (i<eNRNB); i++)
214 dest->n[i]=s1->n[i]+s2->n[i];
217 void print_nrnb(FILE *out, t_nrnb *nrnb)
219 int i;
221 for(i=0; (i<eNRNB); i++)
222 if (nrnb->n[i] > 0)
223 fprintf(out," %-26s %10.0f.\n",nbdata[i].name,nrnb->n[i]);
226 void _inc_nrnb(t_nrnb *nrnb,int enr,int inc,char *file,int line)
228 nrnb->n[enr]+=inc;
229 #ifdef DEBUG_NRNB
230 printf("nrnb %15s(%2d) incremented with %8d from file %s line %d\n",
231 nbdata[enr].name,enr,inc,file,line);
232 #endif
235 void print_flop(FILE *out,t_nrnb *nrnb,double *nbfs,double *mflop)
237 int i;
238 double mni,frac,tfrac,tflop;
239 const char *myline = "-----------------------------------------------------------------------------";
241 *nbfs = 0.0;
242 for(i=0; (i<eNR_NBKERNEL_ALLVSALLGB); i++) {
243 if (strstr(nbdata[i].name,"W3-W3") != NULL)
244 *nbfs += 9e-6*nrnb->n[i];
245 else if (strstr(nbdata[i].name,"W3") != NULL)
246 *nbfs += 3e-6*nrnb->n[i];
247 else if (strstr(nbdata[i].name,"W4-W4") != NULL)
248 *nbfs += 10e-6*nrnb->n[i];
249 else if (strstr(nbdata[i].name,"W4") != NULL)
250 *nbfs += 4e-6*nrnb->n[i];
251 else
252 *nbfs += 1e-6*nrnb->n[i];
254 tflop=0;
255 for(i=0; (i<eNRNB); i++)
256 tflop+=1e-6*nrnb->n[i]*nbdata[i].flop;
258 if (tflop == 0) {
259 fprintf(out,"No MEGA Flopsen this time\n");
260 return;
262 if (out) {
263 fprintf(out,"\n\tM E G A - F L O P S A C C O U N T I N G\n\n");
266 if (out)
268 fprintf(out," NB=Group-cutoff nonbonded kernels NxN=N-by-N cluster Verlet kernels\n");
269 fprintf(out," RF=Reaction-Field VdW=Van der Waals QSTab=quadratic-spline table\n");
270 fprintf(out," W3=SPC/TIP3p W4=TIP4p (single or pairs)\n");
271 fprintf(out," V&F=Potential and force V=Potential only F=Force only\n\n");
273 fprintf(out," %-32s %16s %15s %7s\n",
274 "Computing:","M-Number","M-Flops","% Flops");
275 fprintf(out,"%s\n",myline);
277 *mflop=0.0;
278 tfrac=0.0;
279 for(i=0; (i<eNRNB); i++) {
280 mni = 1e-6*nrnb->n[i];
281 *mflop += mni*nbdata[i].flop;
282 frac = 100.0*mni*nbdata[i].flop/tflop;
283 tfrac += frac;
284 if (out && mni != 0)
285 fprintf(out," %-32s %16.6f %15.3f %6.1f\n",
286 nbdata[i].name,mni,mni*nbdata[i].flop,frac);
288 if (out) {
289 fprintf(out,"%s\n",myline);
290 fprintf(out," %-32s %16s %15.3f %6.1f\n",
291 "Total","",*mflop,tfrac);
292 fprintf(out,"%s\n\n",myline);
296 void print_perf(FILE *out,double nodetime,double realtime,int nprocs,
297 gmx_large_int_t nsteps,real delta_t,
298 double nbfs,double mflop,
299 int omp_nth_pp)
301 real runtime;
303 fprintf(out,"\n");
305 if (realtime > 0)
307 fprintf(out,"%12s %12s %12s %10s\n","","Core t (s)","Wall t (s)","(%)");
308 fprintf(out,"%12s %12.3f %12.3f %10.1f\n","Time:",
309 nodetime, realtime, 100.0*nodetime/realtime);
310 /* only print day-hour-sec format if realtime is more than 30 min */
311 if (realtime > 30*60)
313 fprintf(out,"%12s %12s","","");
314 pr_difftime(out,realtime);
316 if (delta_t > 0)
318 mflop = mflop/realtime;
319 runtime = nsteps*delta_t;
321 if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
323 fprintf(out,"%12s %12s %12s\n",
324 "","(ns/day)","(hour/ns)");
325 fprintf(out,"%12s %12.3f %12.3f\n","Performance:",
326 runtime*24*3.6/realtime,1000*realtime/(3600*runtime));
328 else
330 fprintf(out,"%12s %12s %12s %12s %12s\n",
331 "","(Mnbf/s)",(mflop > 1000) ? "(GFlops)" : "(MFlops)",
332 "(ns/day)","(hour/ns)");
333 fprintf(out,"%12s %12.3f %12.3f %12.3f %12.3f\n","Performance:",
334 nbfs/realtime,(mflop > 1000) ? (mflop/1000) : mflop,
335 runtime*24*3.6/realtime,1000*realtime/(3600*runtime));
338 else
340 if (getenv("GMX_DETAILED_PERF_STATS") == NULL)
342 fprintf(out,"%12s %14s\n",
343 "","(steps/hour)");
344 fprintf(out,"%12s %14.1f\n","Performance:",
345 nsteps*3600.0/realtime);
347 else
349 fprintf(out,"%12s %12s %12s %14s\n",
350 "","(Mnbf/s)",(mflop > 1000) ? "(GFlops)" : "(MFlops)",
351 "(steps/hour)");
352 fprintf(out,"%12s %12.3f %12.3f %14.1f\n","Performance:",
353 nbfs/realtime,(mflop > 1000) ? (mflop/1000) : mflop,
354 nsteps*3600.0/realtime);
360 int cost_nrnb(int enr)
362 return nbdata[enr].flop;
365 const char *nrnb_str(int enr)
367 return nbdata[enr].name;
370 static const int force_index[]={
371 eNR_BONDS, eNR_ANGLES, eNR_PROPER, eNR_IMPROPER,
372 eNR_RB, eNR_DISRES, eNR_ORIRES, eNR_POSRES,
373 eNR_NS,
375 #define NFORCE_INDEX asize(force_index)
377 static const int constr_index[]={
378 eNR_SHAKE, eNR_SHAKE_RIJ, eNR_SETTLE, eNR_UPDATE, eNR_PCOUPL,
379 eNR_CONSTR_VIR,eNR_CONSTR_V
381 #define NCONSTR_INDEX asize(constr_index)
383 static double pr_av(FILE *log,t_commrec *cr,
384 double fav,double ftot[],const char *title)
386 int i,perc;
387 double dperc,unb;
389 unb=0;
390 if (fav > 0) {
391 fav /= cr->nnodes - cr->npmenodes;
392 fprintf(log,"\n %-26s",title);
393 for(i=0; (i<cr->nnodes); i++) {
394 dperc=(100.0*ftot[i])/fav;
395 unb=max(unb,dperc);
396 perc=dperc;
397 fprintf(log,"%3d ",perc);
399 if (unb > 0) {
400 perc=10000.0/unb;
401 fprintf(log,"%6d%%\n\n",perc);
403 else
404 fprintf(log,"\n\n");
406 return unb;
409 void pr_load(FILE *log,t_commrec *cr,t_nrnb nrnb[])
411 int i,j,perc;
412 double dperc,unb,uf,us;
413 double *ftot,fav;
414 double *stot,sav;
415 t_nrnb *av;
417 snew(av,1);
418 snew(ftot,cr->nnodes);
419 snew(stot,cr->nnodes);
420 init_nrnb(av);
421 for(i=0; (i<cr->nnodes); i++) {
422 add_nrnb(av,av,&(nrnb[i]));
423 /* Cost due to forces */
424 for(j=0; (j<eNR_NBKERNEL_ALLVSALLGB); j++)
425 ftot[i]+=nrnb[i].n[j]*cost_nrnb(j);
426 for(j=0; (j<NFORCE_INDEX); j++)
427 ftot[i]+=nrnb[i].n[force_index[j]]*cost_nrnb(force_index[j]);
428 /* Due to shake */
429 for(j=0; (j<NCONSTR_INDEX); j++) {
430 stot[i]+=nrnb[i].n[constr_index[j]]*cost_nrnb(constr_index[j]);
433 for(j=0; (j<eNRNB); j++)
434 av->n[j]=av->n[j]/(double)(cr->nnodes - cr->npmenodes);
436 fprintf(log,"\nDetailed load balancing info in percentage of average\n");
438 fprintf(log," Type NODE:");
439 for(i=0; (i<cr->nnodes); i++)
440 fprintf(log,"%3d ",i);
441 fprintf(log,"Scaling\n");
442 fprintf(log,"---------------------------");
443 for(i=0; (i<cr->nnodes); i++)
444 fprintf(log,"----");
445 fprintf(log,"-------\n");
447 for(j=0; (j<eNRNB); j++) {
448 unb=100.0;
449 if (av->n[j] > 0) {
450 fprintf(log," %-26s",nrnb_str(j));
451 for(i=0; (i<cr->nnodes); i++) {
452 dperc=(100.0*nrnb[i].n[j])/av->n[j];
453 unb=max(unb,dperc);
454 perc=dperc;
455 fprintf(log,"%3d ",perc);
457 if (unb > 0) {
458 perc=10000.0/unb;
459 fprintf(log,"%6d%%\n",perc);
461 else
462 fprintf(log,"\n");
465 fav=sav=0;
466 for(i=0; (i<cr->nnodes); i++) {
467 fav+=ftot[i];
468 sav+=stot[i];
470 uf=pr_av(log,cr,fav,ftot,"Total Force");
471 us=pr_av(log,cr,sav,stot,"Total Constr.");
473 unb=(uf*fav+us*sav)/(fav+sav);
474 if (unb > 0) {
475 unb=10000.0/unb;
476 fprintf(log,"\nTotal Scaling: %.0f%% of max performance\n\n",unb);