renamed GMX_THREADS to GMX_THREAD_MPI
[alexxy/gromacs.git] / src / mdlib / gmx_wallcycle.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
36
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <string.h>
42 #include "gmx_wallcycle.h"
43 #include "gmx_cyclecounter.h"
44 #include "smalloc.h"
45 #include "gmx_fatal.h"
46
47 #ifdef GMX_LIB_MPI
48 #include <mpi.h>
49 #endif
50 #ifdef GMX_THREAD_MPI
51 #include "tmpi.h"
52 #endif
53
54 typedef struct
55 {
56     int          n;
57     gmx_cycles_t c;
58     gmx_cycles_t start;
59     gmx_cycles_t last;
60 } wallcc_t;
61
62 typedef struct gmx_wallcycle
63 {
64     wallcc_t     *wcc;
65     /* variables for testing/debugging */
66     gmx_bool         wc_barrier;
67     wallcc_t     *wcc_all;
68     int          wc_depth;
69     int          ewc_prev;
70     gmx_cycles_t cycle_prev;
71     gmx_large_int_t   reset_counters;
72 #ifdef GMX_MPI
73     MPI_Comm     mpi_comm_mygroup;
74 #endif
75     int          omp_nthreads;
76 } gmx_wallcycle_t_t;
77
78 /* Each name should not exceed 19 characters */
79 static const char *wcn[ewcNR] =
80 { "Run", "Step", "PP during PME", "Domain decomp.", "DD comm. load", "DD comm. bounds", "Vsite constr.", "Send X to PME", "Comm. coord.", "Neighbor search", "Born radii", "Force", "Wait + Comm. F", "PME mesh", "PME redist. X/F", "PME spread/gather", "PME 3D-FFT", "PME 3D-FFT Comm.", "PME solve", "Wait + Comm. X/F", "Wait + Recv. PME F", "Vsite spread", "Write traj.", "Update", "Constraints", "Comm. energies", "Enforced rotation", "Add rot. forces", "Test" };
81
82 gmx_bool wallcycle_have_counter(void)
83 {
84   return gmx_cycles_have_counter();
85 }
86
87 gmx_wallcycle_t wallcycle_init(FILE *fplog,int resetstep,t_commrec *cr, int omp_nthreads)
88 {
89     gmx_wallcycle_t wc;
90     
91     
92     if (!wallcycle_have_counter())
93     {
94         return NULL;
95     }
96
97     snew(wc,1);
98
99     wc->wc_barrier = FALSE;
100     wc->wcc_all    = NULL;
101     wc->wc_depth   = 0;
102     wc->ewc_prev   = -1;
103     wc->reset_counters = resetstep;
104     wc->omp_nthreads = omp_nthreads;
105
106 #ifdef GMX_MPI
107     if (PAR(cr) && getenv("GMX_CYCLE_BARRIER") != NULL)
108     {
109         if (fplog) 
110         {
111             fprintf(fplog,"\nWill call MPI_Barrier before each cycle start/stop call\n\n");
112         }
113         wc->wc_barrier = TRUE;
114         wc->mpi_comm_mygroup = cr->mpi_comm_mygroup;
115     }
116 #endif
117
118     snew(wc->wcc,ewcNR);
119     if (getenv("GMX_CYCLE_ALL") != NULL)
120     {
121 /*#ifndef GMX_THREAD_MPI*/
122         if (fplog) 
123         {
124             fprintf(fplog,"\nWill time all the code during the run\n\n");
125         }
126         snew(wc->wcc_all,ewcNR*ewcNR);
127 /*#else*/
128         gmx_fatal(FARGS, "GMX_CYCLE_ALL is incompatible with threaded code");
129 /*#endif*/
130     }
131     
132     return wc;
133 }
134
135 void wallcycle_destroy(gmx_wallcycle_t wc)
136 {
137     if (wc == NULL)
138     {
139         return;
140     }
141     
142     if (wc->wcc != NULL)
143     {
144         sfree(wc->wcc);
145     }
146     if (wc->wcc_all != NULL)
147     {
148         sfree(wc->wcc_all);
149     }
150     sfree(wc);
151 }
152
153 static void wallcycle_all_start(gmx_wallcycle_t wc,int ewc,gmx_cycles_t cycle)
154 {
155     wc->ewc_prev = ewc;
156     wc->cycle_prev = cycle;
157 }
158
159 static void wallcycle_all_stop(gmx_wallcycle_t wc,int ewc,gmx_cycles_t cycle)
160 {
161     wc->wcc_all[wc->ewc_prev*ewcNR+ewc].n += 1;
162     wc->wcc_all[wc->ewc_prev*ewcNR+ewc].c += cycle - wc->cycle_prev;
163 }
164
165 void wallcycle_start(gmx_wallcycle_t wc, int ewc)
166 {
167     gmx_cycles_t cycle;
168
169     if (wc == NULL)
170     {
171         return;
172     }
173
174 #ifdef GMX_MPI
175     if (wc->wc_barrier)
176     {
177         MPI_Barrier(wc->mpi_comm_mygroup);
178     }
179 #endif
180
181     cycle = gmx_cycles_read();
182     wc->wcc[ewc].start = cycle;
183     if (wc->wcc_all != NULL)
184     {
185         wc->wc_depth++;
186         if (ewc == ewcRUN)
187         {
188             wallcycle_all_start(wc,ewc,cycle);
189         }
190         else if (wc->wc_depth == 3)
191         {
192             wallcycle_all_stop(wc,ewc,cycle);
193         }
194     }
195 }
196
197 double wallcycle_stop(gmx_wallcycle_t wc, int ewc)
198 {
199     gmx_cycles_t cycle,last;
200     
201     if (wc == NULL)
202     {
203         return 0;
204     }
205     
206 #ifdef GMX_MPI
207     if (wc->wc_barrier)
208     {
209         MPI_Barrier(wc->mpi_comm_mygroup);
210     }
211 #endif
212     
213     cycle = gmx_cycles_read();
214     last = cycle - wc->wcc[ewc].start;
215     wc->wcc[ewc].c += last;
216     wc->wcc[ewc].n++;
217     if (wc->wcc_all)
218     {
219         wc->wc_depth--;
220         if (ewc == ewcRUN)
221         {
222             wallcycle_all_stop(wc,ewc,cycle);
223         }
224         else if (wc->wc_depth == 2)
225         {
226             wallcycle_all_start(wc,ewc,cycle);
227         }
228     }
229
230     return last;
231 }
232
233 void wallcycle_reset_all(gmx_wallcycle_t wc)
234 {
235     int i;
236
237     if (wc == NULL)
238     {
239         return;
240     }
241
242     for(i=0; i<ewcNR; i++)
243     {
244         wc->wcc[i].n = 0;
245         wc->wcc[i].c = 0;
246         wc->wcc[i].start = 0;
247         wc->wcc[i].last = 0;
248     }
249 }
250
251 static gmx_bool pme_subdivision(int ewc)
252 {
253     return (ewc >= ewcPME_REDISTXF && ewc <= ewcPME_SOLVE);
254 }
255
256 void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc,double cycles[])
257 {
258     wallcc_t *wcc;
259     double cycles_n[ewcNR],buf[ewcNR],*cyc_all,*buf_all;
260     int    i;
261
262     if (wc == NULL)
263     {
264         return;
265     }
266
267     wcc = wc->wcc;
268
269     if (wc->omp_nthreads>1)
270     {
271         for(i=0; i<ewcNR; i++)
272         {
273             if (pme_subdivision(i) || i==ewcPMEMESH || (i==ewcRUN && cr->duty == DUTY_PME))
274             {
275                 wcc[i].c *= wc->omp_nthreads;
276             }
277         }
278     }
279
280     if (wcc[ewcDDCOMMLOAD].n > 0)
281     {
282         wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMLOAD].c;
283     }
284     if (wcc[ewcDDCOMMBOUND].n > 0)
285     {
286         wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMBOUND].c;
287     }
288     if (wcc[ewcPME_FFTCOMM].n > 0)
289     {
290         wcc[ewcPME_FFT].c -= wcc[ewcPME_FFTCOMM].c;
291     }
292
293     if (cr->npmenodes == 0)
294     {
295         /* All nodes do PME (or no PME at all) */
296         if (wcc[ewcPMEMESH].n > 0)
297         {
298             wcc[ewcFORCE].c -= wcc[ewcPMEMESH].c;
299         }
300     }
301     else
302     {
303         /* The are PME-only nodes */
304         if (wcc[ewcPMEMESH].n > 0)
305         {
306             /* This must be a PME only node, calculate the Wait + Comm. time */
307             wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
308         }
309     }
310     
311     /* Store the cycles in a double buffer for summing */
312     for(i=0; i<ewcNR; i++)
313     {
314         cycles_n[i] = (double)wcc[i].n;
315         cycles[i]   = (double)wcc[i].c;
316     }
317     
318 #ifdef GMX_MPI
319     if (cr->nnodes > 1)
320     {
321         MPI_Allreduce(cycles_n,buf,ewcNR,MPI_DOUBLE,MPI_MAX,
322                       cr->mpi_comm_mysim);
323         for(i=0; i<ewcNR; i++)
324         {
325             wcc[i].n = (int)(buf[i] + 0.5);
326         }
327         MPI_Allreduce(cycles,buf,ewcNR,MPI_DOUBLE,MPI_SUM,
328                       cr->mpi_comm_mysim);
329         for(i=0; i<ewcNR; i++)
330         {
331             cycles[i] = buf[i];
332         }
333
334         if (wc->wcc_all != NULL)
335         {
336             snew(cyc_all,ewcNR*ewcNR);
337             snew(buf_all,ewcNR*ewcNR);
338             for(i=0; i<ewcNR*ewcNR; i++)
339             {
340                 cyc_all[i] = wc->wcc_all[i].c;
341             }
342             MPI_Allreduce(cyc_all,buf_all,ewcNR*ewcNR,MPI_DOUBLE,MPI_SUM,
343                           cr->mpi_comm_mysim);
344             for(i=0; i<ewcNR*ewcNR; i++)
345             {
346                 wc->wcc_all[i].c = buf_all[i];
347             }
348             sfree(buf_all);
349             sfree(cyc_all);
350         }
351     }
352 #endif
353 }
354
355 static void print_cycles(FILE *fplog, double c2t, const char *name, int nnodes,
356                          int n, double c, double tot)
357 {
358     char num[11];
359   
360     if (c > 0)
361     {
362         if (n > 0)
363         {
364             sprintf(num,"%10d",n);
365         }
366         else
367         {
368             sprintf(num,"          ");
369         }
370         fprintf(fplog," %-19s %4d %10s %12.3f %10.1f   %5.1f\n",
371                 name,nnodes,num,c*1e-9,c*c2t,100*c/tot);
372     }
373 }
374
375
376 void wallcycle_print(FILE *fplog, int nnodes, int npme, double realtime,
377                      gmx_wallcycle_t wc, double cycles[])
378 {
379     double c2t,tot,sum;
380     int    i,j,npp;
381     char   buf[STRLEN];
382     const char *myline = "-----------------------------------------------------------------------";
383     
384     if (wc == NULL)
385     {
386         return;
387     }
388
389     if (npme > 0)
390     {
391         npp = nnodes - npme;
392     }
393     else
394     {
395         npp  = nnodes;
396         npme = nnodes;
397     }
398     tot = cycles[ewcRUN];
399     /* PME part has to be multiplied with number of threads */
400     if (npme == 0)
401     {
402         tot += cycles[ewcPMEMESH]*(wc->omp_nthreads-1);
403     }
404     /* Conversion factor from cycles to seconds */
405     if (tot > 0)
406     {
407       c2t = (npp+npme*wc->omp_nthreads)*realtime/tot;
408     }
409     else
410     {
411       c2t = 0;
412     }
413
414     fprintf(fplog,"\n     R E A L   C Y C L E   A N D   T I M E   A C C O U N T I N G\n\n");
415
416     fprintf(fplog," Computing:         Nodes     Number     G-Cycles    Seconds     %c\n",'%');
417     fprintf(fplog,"%s\n",myline);
418     sum = 0;
419     for(i=ewcPPDURINGPME+1; i<ewcNR; i++)
420     {
421         if (!pme_subdivision(i))
422         {
423             print_cycles(fplog,c2t,wcn[i],
424                          (i==ewcPMEMESH || i==ewcPMEWAITCOMM) ? npme : npp,
425                          wc->wcc[i].n,cycles[i],tot);
426             sum += cycles[i];
427         }
428     }
429     if (wc->wcc_all != NULL)
430     {
431         for(i=0; i<ewcNR; i++)
432         {
433             for(j=0; j<ewcNR; j++)
434             {
435                 sprintf(buf,"%-9s",wcn[i]);
436                 buf[9] = ' ';
437                 sprintf(buf+10,"%-9s",wcn[j]);
438                 buf[19] = '\0';
439                 print_cycles(fplog,c2t,buf,
440                              (i==ewcPMEMESH || i==ewcPMEWAITCOMM) ? npme : npp,
441                              wc->wcc_all[i*ewcNR+j].n,
442                              wc->wcc_all[i*ewcNR+j].c,
443                              tot);
444             }
445         }
446     }
447     print_cycles(fplog,c2t,"Rest",npp,0,tot-sum,tot);
448     fprintf(fplog,"%s\n",myline);
449     print_cycles(fplog,c2t,"Total",nnodes,0,tot,tot);
450     fprintf(fplog,"%s\n",myline);
451     
452     if (wc->wcc[ewcPMEMESH].n > 0)
453     {
454         fprintf(fplog,"%s\n",myline);
455         for(i=ewcPPDURINGPME+1; i<ewcNR; i++)
456         {
457             if (pme_subdivision(i))
458             {
459                 print_cycles(fplog,c2t,wcn[i],
460                              (i>=ewcPMEMESH || i<=ewcPME_SOLVE) ? npme : npp,
461                              wc->wcc[i].n,cycles[i],tot);
462             }
463         }
464         fprintf(fplog,"%s\n",myline);
465     }
466
467     if (cycles[ewcMoveE] > tot*0.05)
468     {
469         sprintf(buf,
470                 "NOTE: %d %% of the run time was spent communicating energies,\n"
471                 "      you might want to use the -gcom option of mdrun\n",
472                 (int)(100*cycles[ewcMoveE]/tot+0.5));
473         if (fplog)
474         {
475             fprintf(fplog,"\n%s\n",buf);
476         }
477         /* Only the sim master calls this function, so always print to stderr */
478         fprintf(stderr,"\n%s\n",buf);
479     }
480 }
481
482 extern gmx_large_int_t wcycle_get_reset_counters(gmx_wallcycle_t wc)
483 {
484     if (wc == NULL)
485     {
486         return -1;
487     }
488     
489     return wc->reset_counters;
490 }
491
492 extern void wcycle_set_reset_counters(gmx_wallcycle_t wc, gmx_large_int_t reset_counters)
493 {
494     if (wc == NULL)
495         return;
496
497     wc->reset_counters = reset_counters;
498 }