1 /* ----------------------------------------------------------------------
2    SPARTA - Stochastic PArallel Rarefied-gas Time-accurate Analyzer
3    http://sparta.sandia.gov
4    Steve Plimpton, sjplimp@sandia.gov, Michael Gallis, magalli@sandia.gov
5    Sandia National Laboratories
6 
7    Copyright (2014) Sandia Corporation.  Under the terms of Contract
8    DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
9    certain rights in this software.  This software is distributed under
10    the GNU General Public License.
11 
12    See the README file in the top-level SPARTA directory.
13 ------------------------------------------------------------------------- */
14 
15 #include "spatype.h"
16 #include "mpi.h"
17 #include "math.h"
18 #include "stdlib.h"
19 #include "string.h"
20 #include "stats.h"
21 #include "update.h"
22 #include "particle.h"
23 #include "collide.h"
24 #include "domain.h"
25 #include "grid.h"
26 #include "surf.h"
27 #include "surf_collide.h"
28 #include "surf_react.h"
29 #include "modify.h"
30 #include "fix.h"
31 #include "compute.h"
32 #include "input.h"
33 #include "variable.h"
34 #include "output.h"
35 #include "timer.h"
36 #include "memory.h"
37 #include "error.h"
38 
39 using namespace SPARTA_NS;
40 
41 // customize a new keyword by adding to this list:
42 
43 // step,elapsed,elaplong,dt,cpu,tpcpu,spcpu,wall,
44 // np,ntouch,ncomm,nbound,nexit,nscoll,nscheck,ncoll,nattempt,nreact,nsreact,
45 // npave,ntouchave,ncommave,nboundave,nexitave,nscollave,nscheckave,
46 // ncollave,nattemptave,nreactave,nsreactave,
47 // ngrid,nsplit,maxlevel,
48 // vol,lx,ly,lz,xlo,xhi,ylo,yhi,zlo,zhi
49 
50 enum{INT,FLOAT,BIGINT};
51 enum{SCALAR,VECTOR,ARRAY};
52 
53 #define INVOKED_SCALAR 1
54 #define INVOKED_VECTOR 2
55 #define INVOKED_ARRAY 4
56 
57 #define MAXLINE 8192               // make this 4x longer than Input::MAXLINE
58 #define DELTA 8
59 
60 /* ---------------------------------------------------------------------- */
61 
Stats(SPARTA * sparta)62 Stats::Stats(SPARTA *sparta) : Pointers(sparta)
63 {
64   MPI_Comm_rank(world,&me);
65 
66   MPI_Barrier(world);
67   wall0 = MPI_Wtime();
68 
69   line = new char[MAXLINE];
70 
71   keyword = NULL;
72   vfunc = NULL;
73   vtype = NULL;
74 
75   field2index = NULL;
76   argindex1 = NULL;
77   argindex2 = NULL;
78 
79   // default args
80 
81   char **arg = new char*[3];
82   arg[0] = (char *) "step";
83   arg[1] = (char *) "cpu";
84   arg[2] = (char *) "np";
85 
86   nfield = 3;
87   allocate();
88   set_fields(3,arg);
89 
90   delete [] arg;
91 
92   // stats_modify defaults
93 
94   flushflag = 0;
95 
96   // format strings
97 
98   char *bigint_format = (char *) BIGINT_FORMAT;
99 
100   format_float_def = (char *) "%12.8g";
101   format_int_def = (char *) "%8d";
102   sprintf(format_bigint_def,"%%8%s",&bigint_format[1]);
103 
104   format_line_user = NULL;
105   format_float_user = NULL;
106   format_int_user = NULL;
107   format_bigint_user = NULL;
108 }
109 
110 /* ---------------------------------------------------------------------- */
111 
~Stats()112 Stats::~Stats()
113 {
114   delete [] line;
115   deallocate();
116 
117   // format strings
118 
119   delete [] format_line_user;
120   delete [] format_float_user;
121   delete [] format_int_user;
122   delete [] format_bigint_user;
123 }
124 
125 /* ---------------------------------------------------------------------- */
126 
init()127 void Stats::init()
128 {
129   int i,n;
130 
131   // set format string for each field
132   // add trailing '/n' to last value
133 
134   char *format_line = NULL;
135   if (format_line_user) {
136     int n = strlen(format_line_user) + 1;
137     format_line = new char[n];
138     strcpy(format_line,format_line_user);
139   }
140 
141   char *ptr,*format_line_ptr;
142   for (i = 0; i < nfield; i++) {
143     format[i][0] = '\0';
144 
145     if (format_line) {
146       if (i == 0) format_line_ptr = strtok(format_line," \0");
147       else format_line_ptr = strtok(NULL," \0");
148     }
149 
150     if (format_column_user[i]) ptr = format_column_user[i];
151     else if (vtype[i] == FLOAT) {
152       if (format_float_user) ptr = format_float_user;
153       else if (format_line_user) ptr = format_line_ptr;
154       else ptr = format_float_def;
155     } else if (vtype[i] == INT) {
156       if (format_int_user) ptr = format_int_user;
157       else if (format_line_user) ptr = format_line_ptr;
158       else ptr = format_int_def;
159     } else if (vtype[i] == BIGINT) {
160       if (format_bigint_user) ptr = format_bigint_user;
161       else if (format_line_user) ptr = format_line_ptr;
162       else ptr = format_bigint_def;
163     }
164 
165     n = strlen(format[i]);
166     sprintf(&format[i][n],"%s ",ptr);
167   }
168   strcat(format[nfield-1],"\n");
169 
170   delete [] format_line;
171 
172   // find current ptr for each SurfCollide and SurfReact ID
173 
174   int m;
175 
176   for (int i = 0; i < nsurfcollide; i++) {
177     m = surf->find_collide(id_surf_collide[i]);
178     if (m < 0) error->all(FLERR,"Could not find stats surf collide ID");
179     sc[i] = surf->sc[m];
180   }
181   for (int i = 0; i < nsurfreact; i++) {
182     m = surf->find_react(id_surf_react[i]);
183     if (m < 0) error->all(FLERR,"Could not find stats surf collide ID");
184     sr[i] = surf->sr[m];
185   }
186 
187   // find current ptr for each Compute ID
188 
189   for (int i = 0; i < ncompute; i++) {
190     m = modify->find_compute(id_compute[i]);
191     if (m < 0) error->all(FLERR,"Could not find stats compute ID");
192     computes[i] = modify->compute[m];
193   }
194 
195   // find current ptr for each Fix ID
196   // check that fix frequency is acceptable with stats output frequency
197 
198   for (int i = 0; i < nfix; i++) {
199     m = modify->find_fix(id_fix[i]);
200     if (m < 0) error->all(FLERR,"Could not find stats fix ID");
201     fixes[i] = modify->fix[m];
202     if (output->stats_every % fixes[i]->global_freq)
203       error->all(FLERR,"Stats and fix not computed at compatible times");
204   }
205 
206   // find current ptr for each Variable ID
207 
208   for (int i = 0; i < nvariable; i++) {
209     m = input->variable->find(id_variable[i]);
210     if (m < 0) error->all(FLERR,"Could not find stats variable name");
211     variables[i] = m;
212   }
213 }
214 
215 /* ---------------------------------------------------------------------- */
216 
header()217 void Stats::header()
218 {
219   int loc = 0;
220   for (int i = 0; i < nfield; i++)
221     loc += sprintf(&line[loc],"%s ",keyword[i]);
222   sprintf(&line[loc],"\n");
223 
224   if (me == 0) {
225     if (screen) fprintf(screen,"%s",line);
226     if (logfile) fprintf(logfile,"%s",line);
227   }
228 }
229 
230 /* ---------------------------------------------------------------------- */
231 
compute(int flag)232 void Stats::compute(int flag)
233 {
234   int i;
235 
236   firststep = flag;
237 
238   // invoke Compute methods needed for stats keywords
239 
240   for (i = 0; i < ncompute; i++)
241     if (compute_which[i] == SCALAR) {
242       if (!(computes[i]->invoked_flag & INVOKED_SCALAR)) {
243 	computes[i]->compute_scalar();
244 	computes[i]->invoked_flag |= INVOKED_SCALAR;
245       }
246     } else if (compute_which[i] == VECTOR) {
247       if (!(computes[i]->invoked_flag & INVOKED_VECTOR)) {
248 	computes[i]->compute_vector();
249 	computes[i]->invoked_flag |= INVOKED_VECTOR;
250       }
251     } else if (compute_which[i] == ARRAY) {
252       if (!(computes[i]->invoked_flag & INVOKED_ARRAY)) {
253 	computes[i]->compute_array();
254 	computes[i]->invoked_flag |= INVOKED_ARRAY;
255       }
256     }
257 
258   // add each stat value to line with its specific format
259 
260   int loc = 0;
261   for (ifield = 0; ifield < nfield; ifield++) {
262     (this->*vfunc[ifield])();
263     if (vtype[ifield] == FLOAT)
264       loc += sprintf(&line[loc],format[ifield],dvalue);
265     else if (vtype[ifield] == INT)
266       loc += sprintf(&line[loc],format[ifield],ivalue);
267     else if (vtype[ifield] == BIGINT) {
268       loc += sprintf(&line[loc],format[ifield],bivalue);
269     }
270   }
271 
272   // print line to screen and logfile
273 
274   if (me == 0) {
275     if (screen) fprintf(screen,"%s",line);
276     if (logfile) {
277       fprintf(logfile,"%s",line);
278       if (flushflag) fflush(logfile);
279     }
280   }
281 }
282 
283 /* ----------------------------------------------------------------------
284    modify stats parameters
285 ------------------------------------------------------------------------- */
286 
modify_params(int narg,char ** arg)287 void Stats::modify_params(int narg, char **arg)
288 {
289   if (narg == 0) error->all(FLERR,"Illegal stats_modify command");
290 
291   int iarg = 0;
292   while (iarg < narg) {
293     if (strcmp(arg[iarg],"every") == 0) {
294       if (iarg+2 > narg) error->all(FLERR,"Illegal stats_modify command");
295       if (strstr(arg[iarg+1],"v_") == arg[iarg+1]) {
296 	delete [] output->var_stats;
297 	int n = strlen(&arg[iarg+1][2]) + 1;
298 	output->var_stats = new char[n];
299 	strcpy(output->var_stats,&arg[iarg+1][2]);
300       } else error->all(FLERR,"Illegal stats_modify command");
301       output->stats_every = 0;
302       iarg += 2;
303 
304     } else if (strcmp(arg[iarg],"flush") == 0) {
305       if (iarg+2 > narg) error->all(FLERR,"Illegal stats_modify command");
306       if (strcmp(arg[iarg+1],"no") == 0) flushflag = 0;
307       else if (strcmp(arg[iarg+1],"yes") == 0) flushflag = 1;
308       else error->all(FLERR,"Illegal stats_modify command");
309       iarg += 2;
310 
311     } else if (strcmp(arg[iarg],"format") == 0) {
312       if (iarg+2 > narg) error->all(FLERR,"Illegal stats_modify command");
313 
314       if (strcmp(arg[iarg+1],"none") == 0) {
315         delete [] format_line_user;
316         delete [] format_int_user;
317         delete [] format_bigint_user;
318         delete [] format_float_user;
319         format_line_user = NULL;
320         format_int_user = NULL;
321         format_bigint_user = NULL;
322         format_float_user = NULL;
323         for (int i = 0; i < nfield; i++) {
324           delete [] format_column_user[i];
325           format_column_user[i] = NULL;
326         }
327         iarg += 2;
328         continue;
329       }
330 
331       if (iarg+3 > narg) error->all(FLERR,"Illegal stats_modify command");
332 
333       if (strcmp(arg[iarg+1],"line") == 0) {
334         delete [] format_line_user;
335         int n = strlen(arg[iarg+2]) + 1;
336         format_line_user = new char[n];
337         strcpy(format_line_user,arg[iarg+2]);
338       } else if (strcmp(arg[iarg+1],"int") == 0) {
339         if (format_int_user) delete [] format_int_user;
340         int n = strlen(arg[iarg+2]) + 1;
341         format_int_user = new char[n];
342         strcpy(format_int_user,arg[iarg+2]);
343         if (format_bigint_user) delete [] format_bigint_user;
344         n = strlen(format_int_user) + 8;
345         format_bigint_user = new char[n];
346         // replace "d" in format_int_user with bigint format specifier
347         // use of &str[1] removes leading '%' from BIGINT_FORMAT string
348         char *ptr = strchr(format_int_user,'d');
349         if (ptr == NULL)
350           error->all(FLERR,
351                      "Stats_modify int format does not contain d character");
352         char str[8];
353         sprintf(str,"%s",BIGINT_FORMAT);
354         *ptr = '\0';
355         sprintf(format_bigint_user,"%s%s%s",format_int_user,&str[1],ptr+1);
356         *ptr = 'd';
357       } else if (strcmp(arg[iarg+1],"float") == 0) {
358         if (format_float_user) delete [] format_float_user;
359         int n = strlen(arg[iarg+2]) + 1;
360         format_float_user = new char[n];
361         strcpy(format_float_user,arg[iarg+2]);
362       } else {
363         int i = input->inumeric(FLERR,arg[iarg+1]) - 1;
364         if (i < 0 || i >= nfield)
365           error->all(FLERR,"Illegal stats_modify command");
366         if (format_column_user[i]) delete [] format_column_user[i];
367         int n = strlen(arg[iarg+2]) + 1;
368         format_column_user[i] = new char[n];
369         strcpy(format_column_user[i],arg[iarg+2]);
370       }
371       iarg += 3;
372 
373     } else error->all(FLERR,"Illegal stats_modify command");
374   }
375 }
376 
377 /* ----------------------------------------------------------------------
378    allocate all per-field memory
379 ------------------------------------------------------------------------- */
380 
allocate()381 void Stats::allocate()
382 {
383   int n = nfield;
384 
385   keyword = new char*[n];
386   for (int i = 0; i < n; i++) keyword[i] = new char[32];
387   vfunc = new FnPtr[n];
388   vtype = new int[n];
389 
390   format = new char*[n];
391   for (int i = 0; i < n; i++) format[i] = new char[32];
392   format_column_user = new char*[n];
393   for (int i = 0; i < n; i++) format_column_user[i] = NULL;
394 
395   field2index = new int[n];
396   argindex1 = new int[n];
397   argindex2 = new int[n];
398 
399   // memory for computes, fixes, variables
400 
401   ncompute = 0;
402   id_compute = new char*[n];
403   compute_which = new int[n];
404   computes = new Compute*[n];
405 
406   nfix = 0;
407   id_fix = new char*[n];
408   fixes = new Fix*[n];
409 
410   nsurfcollide = 0;
411   id_surf_collide = new char*[n];
412   sc = new SurfCollide*[n];
413 
414   nsurfreact = 0;
415   id_surf_react = new char*[n];
416   sr = new SurfReact*[n];
417 
418   nvariable = 0;
419   id_variable = new char*[n];
420   variables = new int[n];
421 }
422 
423 /* ----------------------------------------------------------------------
424    deallocate all per-field memory
425 ------------------------------------------------------------------------- */
426 
deallocate()427 void Stats::deallocate()
428 {
429   int n = nfield;
430 
431   for (int i = 0; i < n; i++) delete [] keyword[i];
432   delete [] keyword;
433   delete [] vfunc;
434   delete [] vtype;
435 
436   for (int i = 0; i < n; i++) delete [] format[i];
437   delete [] format;
438   for (int i = 0; i < n; i++) delete [] format_column_user[i];
439   delete [] format_column_user;
440 
441   delete [] field2index;
442   delete [] argindex1;
443   delete [] argindex2;
444 
445   for (int i = 0; i < ncompute; i++) delete [] id_compute[i];
446   delete [] id_compute;
447   delete [] compute_which;
448   delete [] computes;
449 
450   for (int i = 0; i < nfix; i++) delete [] id_fix[i];
451   delete [] id_fix;
452   delete [] fixes;
453 
454   for (int i = 0; i < nsurfcollide; i++) delete [] id_surf_collide[i];
455   delete [] id_surf_collide;
456   delete [] sc;
457 
458   for (int i = 0; i < nsurfreact; i++) delete [] id_surf_react[i];
459   delete [] id_surf_react;
460   delete [] sr;
461 
462   for (int i = 0; i < nvariable; i++) delete [] id_variable[i];
463   delete [] id_variable;
464   delete [] variables;
465 }
466 
467 /* ----------------------------------------------------------------------
468    set fields of stats output from args
469    called by constructor with default fields
470    called by Input with stats_style command args
471 ------------------------------------------------------------------------- */
472 
set_fields(int narg,char ** arg)473 void Stats::set_fields(int narg, char **arg)
474 {
475   deallocate();
476 
477   // expand args if any have wildcard character "*"
478 
479   int expand = 0;
480   char **earg;
481   int nargnew = input->expand_args(narg,arg,0,earg);
482 
483   if (earg != arg) expand = 1;
484   arg = earg;
485 
486   nfield = nargnew;
487   allocate();
488   nfield = 0;
489 
490   // customize a new keyword by adding to if statement
491 
492   for (int i = 0; i < nargnew; i++) {
493     if (strcmp(arg[i],"step") == 0) {
494       addfield("Step",&Stats::compute_step,BIGINT);
495     } else if (strcmp(arg[i],"elapsed") == 0) {
496       addfield("Elapsed",&Stats::compute_elapsed,BIGINT);
497     } else if (strcmp(arg[i],"elaplong") == 0) {
498       addfield("Elapsed",&Stats::compute_elaplong,BIGINT);
499     } else if (strcmp(arg[i],"dt") == 0) {
500       addfield("Dt",&Stats::compute_dt,FLOAT);
501     } else if (strcmp(arg[i],"cpu") == 0) {
502       addfield("CPU",&Stats::compute_cpu,FLOAT);
503     } else if (strcmp(arg[i],"tpcpu") == 0) {
504       addfield("T/CPU",&Stats::compute_tpcpu,FLOAT);
505     } else if (strcmp(arg[i],"spcpu") == 0) {
506       addfield("S/CPU",&Stats::compute_spcpu,FLOAT);
507     } else if (strcmp(arg[i],"wall") == 0) {
508       addfield("WALL",&Stats::compute_wall,FLOAT);
509 
510     } else if (strcmp(arg[i],"np") == 0) {
511       addfield("Np",&Stats::compute_np,BIGINT);
512     } else if (strcmp(arg[i],"ntouch") == 0) {
513       addfield("Ntouch",&Stats::compute_ntouch,BIGINT);
514     } else if (strcmp(arg[i],"ncomm") == 0) {
515       addfield("Ncomm",&Stats::compute_ncomm,BIGINT);
516     } else if (strcmp(arg[i],"nbound") == 0) {
517       addfield("Nbound",&Stats::compute_nbound,BIGINT);
518     } else if (strcmp(arg[i],"nexit") == 0) {
519       addfield("Nexit",&Stats::compute_nexit,BIGINT);
520     } else if (strcmp(arg[i],"nscoll") == 0) {
521       addfield("Nscoll",&Stats::compute_nscoll,BIGINT);
522     } else if (strcmp(arg[i],"nscheck") == 0) {
523       addfield("Nscheck",&Stats::compute_nscheck,BIGINT);
524     } else if (strcmp(arg[i],"ncoll") == 0) {
525       addfield("Ncoll",&Stats::compute_ncoll,BIGINT);
526     } else if (strcmp(arg[i],"nattempt") == 0) {
527       addfield("Natt",&Stats::compute_nattempt,BIGINT);
528     } else if (strcmp(arg[i],"nreact") == 0) {
529       addfield("Nreact",&Stats::compute_nreact,BIGINT);
530     } else if (strcmp(arg[i],"nsreact") == 0) {
531       addfield("Nsreact",&Stats::compute_nsreact,BIGINT);
532 
533     } else if (strcmp(arg[i],"npave") == 0) {
534       addfield("Npave",&Stats::compute_npave,FLOAT);
535     } else if (strcmp(arg[i],"ntouchave") == 0) {
536       addfield("Ntouchave",&Stats::compute_ntouchave,FLOAT);
537     } else if (strcmp(arg[i],"ncommave") == 0) {
538       addfield("Ncommave",&Stats::compute_ncommave,FLOAT);
539     } else if (strcmp(arg[i],"nboundave") == 0) {
540       addfield("Nboundave",&Stats::compute_nboundave,FLOAT);
541     } else if (strcmp(arg[i],"nexitave") == 0) {
542       addfield("Nexitave",&Stats::compute_nexitave,FLOAT);
543     } else if (strcmp(arg[i],"nscollave") == 0) {
544       addfield("Nscollave",&Stats::compute_nscollave,FLOAT);
545     } else if (strcmp(arg[i],"nscheckave") == 0) {
546       addfield("Nschckave",&Stats::compute_nscheckave,FLOAT);
547     } else if (strcmp(arg[i],"ncollave") == 0) {
548       addfield("Ncollave",&Stats::compute_ncollave,FLOAT);
549     } else if (strcmp(arg[i],"nattemptave") == 0) {
550       addfield("Nattave",&Stats::compute_nattemptave,FLOAT);
551     } else if (strcmp(arg[i],"nreactave") == 0) {
552       addfield("Nreactave",&Stats::compute_nreactave,FLOAT);
553     } else if (strcmp(arg[i],"nsreactave") == 0) {
554       addfield("Nsreactave",&Stats::compute_nsreactave,FLOAT);
555 
556     } else if (strcmp(arg[i],"ngrid") == 0) {
557       addfield("Ngrid",&Stats::compute_ngrid,BIGINT);
558     } else if (strcmp(arg[i],"nsplit") == 0) {
559       addfield("Nsplit",&Stats::compute_nsplit,INT);
560     } else if (strcmp(arg[i],"maxlevel") == 0) {
561       addfield("Maxlevel",&Stats::compute_maxlevel,INT);
562 
563     } else if (strcmp(arg[i],"vol") == 0) {
564       addfield("Volume",&Stats::compute_vol,FLOAT);
565     } else if (strcmp(arg[i],"lx") == 0) {
566       addfield("Lx",&Stats::compute_lx,FLOAT);
567     } else if (strcmp(arg[i],"ly") == 0) {
568       addfield("Ly",&Stats::compute_ly,FLOAT);
569     } else if (strcmp(arg[i],"lz") == 0) {
570       addfield("Lz",&Stats::compute_lz,FLOAT);
571 
572     } else if (strcmp(arg[i],"xlo") == 0) {
573       addfield("Xlo",&Stats::compute_xlo,FLOAT);
574     } else if (strcmp(arg[i],"xhi") == 0) {
575       addfield("Xhi",&Stats::compute_xhi,FLOAT);
576     } else if (strcmp(arg[i],"ylo") == 0) {
577       addfield("Ylo",&Stats::compute_ylo,FLOAT);
578     } else if (strcmp(arg[i],"yhi") == 0) {
579       addfield("Yhi",&Stats::compute_yhi,FLOAT);
580     } else if (strcmp(arg[i],"zlo") == 0) {
581       addfield("Zlo",&Stats::compute_zlo,FLOAT);
582     } else if (strcmp(arg[i],"zhi") == 0) {
583       addfield("Zhi",&Stats::compute_zhi,FLOAT);
584 
585     // surf collide value = s_ID, surf react value = r_ID
586     // count trailing [] and store int arguments
587     // copy = at most 8 chars of ID to pass to addfield
588 
589     } else if ((strncmp(arg[i],"s_",2) == 0) ||
590 	       (strncmp(arg[i],"r_",2) == 0)) {
591 
592       int n = strlen(arg[i]);
593       char *id = new char[n];
594       strcpy(id,&arg[i][2]);
595 
596       // parse zero or one or two trailing brackets from ID
597       // argindex1,argindex2 = int inside each bracket pair, 0 if no bracket
598 
599       char *ptr = strchr(id,'[');
600       if (ptr == NULL) argindex1[nfield] = 0;
601       else {
602 	*ptr = '\0';
603 	argindex1[nfield] = input->variable->int_between_brackets(ptr,0);
604 	ptr++;
605 	if (*ptr == '[') {
606 	  argindex2[nfield] = input->variable->int_between_brackets(ptr,0);
607 	  ptr++;
608 	} else argindex2[nfield] = 0;
609       }
610 
611       if (arg[i][0] == 's') {
612 	n = surf->find_collide(id);
613 	if (n < 0) error->all(FLERR,"Could not find stats surf collide ID");
614 	if (argindex1[nfield] == 0 || argindex2[nfield] > 0)
615 	  error->all(FLERR,"Stats surf collide is not indexed correctly");
616 	if (surf->sc[n]->vector_flag == 0)
617 	  error->all(FLERR,"Stats surf collide does not compute vector");
618 	if (argindex1[nfield] > surf->sc[n]->size_vector)
619 	  error->all(FLERR,"Stats surf collide vector "
620                      "is accessed out-of-range");
621 
622 	field2index[nfield] = add_surf_collide(id);
623 	addfield(arg[i],&Stats::compute_surf_collide,FLOAT);
624 
625       } else if (arg[i][0] == 'r') {
626 	n = surf->find_react(id);
627 	if (n < 0) error->all(FLERR,"Could not find stats surf react ID");
628 	if (argindex1[nfield] == 0 || argindex2[nfield] > 0)
629 	  error->all(FLERR,"Stats surf react is not indexed correctly");
630 	if (surf->sr[n]->vector_flag == 0)
631 	  error->all(FLERR,"Stats surf react does not compute vector");
632 	if (argindex1[nfield] > surf->sr[n]->size_vector)
633 	  error->all(FLERR,"Stats surf react vector is accessed out-of-range");
634 
635 	field2index[nfield] = add_surf_react(id);
636 	addfield(arg[i],&Stats::compute_surf_react,FLOAT);
637       }
638 
639       delete [] id;
640 
641     // compute value = c_ID, fix value = f_ID, variable value = v_ID
642     // count trailing [] and store int arguments
643     // copy = at most 8 chars of ID to pass to addfield
644 
645     } else if ((strncmp(arg[i],"c_",2) == 0) ||
646 	       (strncmp(arg[i],"f_",2) == 0) ||
647 	       (strncmp(arg[i],"v_",2) == 0)) {
648 
649       int n = strlen(arg[i]);
650       char *id = new char[n];
651       strcpy(id,&arg[i][2]);
652 
653       // parse zero or one or two trailing brackets from ID
654       // argindex1,argindex2 = int inside each bracket pair, 0 if no bracket
655 
656       char *ptr = strchr(id,'[');
657       if (ptr == NULL) argindex1[nfield] = 0;
658       else {
659 	*ptr = '\0';
660 	argindex1[nfield] = input->variable->int_between_brackets(ptr,0);
661 	ptr++;
662 	if (*ptr == '[') {
663 	  argindex2[nfield] = input->variable->int_between_brackets(ptr,0);
664 	  ptr++;
665 	} else argindex2[nfield] = 0;
666       }
667 
668       if (arg[i][0] == 'c') {
669 	n = modify->find_compute(id);
670 	if (n < 0) error->all(FLERR,"Could not find stats compute ID");
671 	if (argindex1[nfield] == 0 && modify->compute[n]->scalar_flag == 0)
672 	  error->all(FLERR,"Stats compute does not compute scalar");
673 	if (argindex1[nfield] > 0 && argindex2[nfield] == 0) {
674 	  if (modify->compute[n]->vector_flag == 0)
675 	    error->all(FLERR,"Stats compute does not compute vector");
676 	  if (argindex1[nfield] > modify->compute[n]->size_vector)
677 	    error->all(FLERR,"Stats compute vector is accessed out-of-range");
678 	}
679 	if (argindex1[nfield] > 0 && argindex2[nfield] > 0) {
680 	  if (modify->compute[n]->array_flag == 0)
681 	    error->all(FLERR,"Stats compute does not compute array");
682 	  if (argindex1[nfield] > modify->compute[n]->size_array_rows ||
683 	      argindex2[nfield] > modify->compute[n]->size_array_cols)
684 	    error->all(FLERR,"Stats compute array is accessed out-of-range");
685 	}
686 
687 	if (argindex1[nfield] == 0)
688 	  field2index[nfield] = add_compute(id,SCALAR);
689 	else if (argindex2[nfield] == 0)
690 	  field2index[nfield] = add_compute(id,VECTOR);
691 	else
692 	  field2index[nfield] = add_compute(id,ARRAY);
693 	addfield(arg[i],&Stats::compute_compute,FLOAT);
694 
695       } else if (arg[i][0] == 'f') {
696 	n = modify->find_fix(id);
697 	if (n < 0) error->all(FLERR,"Could not find stats fix ID");
698 	if (argindex1[nfield] == 0 && modify->fix[n]->scalar_flag == 0)
699 	  error->all(FLERR,"Stats fix does not compute scalar");
700 	if (argindex1[nfield] > 0 && argindex2[nfield] == 0) {
701 	  if (modify->fix[n]->vector_flag == 0)
702 	    error->all(FLERR,"Stats fix does not compute vector");
703 	  if (argindex1[nfield] > modify->fix[n]->size_vector)
704 	    error->all(FLERR,"Stats fix vector is accessed out-of-range");
705 	}
706 	if (argindex1[nfield] > 0 && argindex2[nfield] > 0) {
707 	  if (modify->fix[n]->array_flag == 0)
708 	    error->all(FLERR,"Stats fix does not compute array");
709 	  if (argindex1[nfield] > modify->fix[n]->size_array_rows ||
710 	      argindex2[nfield] > modify->fix[n]->size_array_cols)
711 	    error->all(FLERR,"Stats fix array is accessed out-of-range");
712 	}
713 
714 	field2index[nfield] = add_fix(id);
715 	addfield(arg[i],&Stats::compute_fix,FLOAT);
716 
717       } else if (arg[i][0] == 'v') {
718 	n = input->variable->find(id);
719 	if (n < 0) error->all(FLERR,"Could not find stats variable name");
720 	if (input->variable->equal_style(n) == 0)
721 	  error->all(FLERR,"Stats variable is not equal-style variable");
722 	if (argindex1[nfield])
723 	  error->all(FLERR,"Stats variable cannot be indexed");
724 
725 	field2index[nfield] = add_variable(id);
726 	addfield(arg[i],&Stats::compute_variable,FLOAT);
727       }
728 
729       delete [] id;
730 
731     } else error->all(FLERR,"Invalid keyword in stats_style command");
732   }
733 
734   // if wildcard expansion occurred, free earg memory from expand_args()
735 
736   if (expand) {
737     for (int i = 0; i < nargnew; i++) delete [] earg[i];
738     memory->sfree(earg);
739   }
740 }
741 
742 /* ----------------------------------------------------------------------
743    add field to list of quantities to print
744 ------------------------------------------------------------------------- */
745 
addfield(const char * key,FnPtr func,int typeflag)746 void Stats::addfield(const char *key, FnPtr func, int typeflag)
747 {
748   strcpy(keyword[nfield],key);
749   vfunc[nfield] = func;
750   vtype[nfield] = typeflag;
751   nfield++;
752 }
753 
754 /* ----------------------------------------------------------------------
755    add compute ID to list of Compute objects to call
756    return location of where this Compute is in list
757    if already in list with same which, do not add, just return index
758 ------------------------------------------------------------------------- */
759 
add_compute(const char * id,int which)760 int Stats::add_compute(const char *id, int which)
761 {
762   int icompute;
763   for (icompute = 0; icompute < ncompute; icompute++)
764     if ((strcmp(id,id_compute[icompute]) == 0) &&
765 	which == compute_which[icompute]) break;
766   if (icompute < ncompute) return icompute;
767 
768   int n = strlen(id) + 1;
769   id_compute[ncompute] = new char[n];
770   strcpy(id_compute[ncompute],id);
771   compute_which[ncompute] = which;
772   ncompute++;
773   return ncompute-1;
774 }
775 
776 /* ----------------------------------------------------------------------
777    add fix ID to list of Fix objects to call
778 ------------------------------------------------------------------------- */
779 
add_fix(const char * id)780 int Stats::add_fix(const char *id)
781 {
782   int n = strlen(id) + 1;
783   id_fix[nfix] = new char[n];
784   strcpy(id_fix[nfix],id);
785   nfix++;
786   return nfix-1;
787 }
788 
789 /* ----------------------------------------------------------------------
790    add surface collide ID to list of SurfCollide objects to call
791 ------------------------------------------------------------------------- */
792 
add_surf_collide(const char * id)793 int Stats::add_surf_collide(const char *id)
794 {
795   int n = strlen(id) + 1;
796   id_surf_collide[nsurfcollide] = new char[n];
797   strcpy(id_surf_collide[nsurfcollide],id);
798   nsurfcollide++;
799   return nsurfcollide-1;
800 }
801 
802 /* ----------------------------------------------------------------------
803    add surface react ID to list of SurfReact objects to call
804 ------------------------------------------------------------------------- */
805 
add_surf_react(const char * id)806 int Stats::add_surf_react(const char *id)
807 {
808   int n = strlen(id) + 1;
809   id_surf_react[nsurfreact] = new char[n];
810   strcpy(id_surf_react[nsurfreact],id);
811   nsurfreact++;
812   return nsurfreact-1;
813 }
814 
815 /* ----------------------------------------------------------------------
816    add variable ID to list of Variables to evaluate
817 ------------------------------------------------------------------------- */
818 
add_variable(const char * id)819 int Stats::add_variable(const char *id)
820 {
821   int n = strlen(id) + 1;
822   id_variable[nvariable] = new char[n];
823   strcpy(id_variable[nvariable],id);
824   nvariable++;
825   return nvariable-1;
826 }
827 
828 /* ----------------------------------------------------------------------
829    compute a single stats value, word is any stats keyword
830    called when a variable is evaluated by Variable class
831    return value as double in answer
832    return 0 if str is recoginzed keyword, 1 if unrecognized
833    customize a new keyword by adding to if statement
834 ------------------------------------------------------------------------- */
835 
evaluate_keyword(char * word,double * answer)836 int Stats::evaluate_keyword(char *word, double *answer)
837 {
838   // invoke a lo-level stats routine to compute the variable value
839 
840   if (strcmp(word,"step") == 0) {
841     compute_step();
842     dvalue = bivalue;
843 
844   } else if (strcmp(word,"elapsed") == 0) {
845     if (update->runflag == 0)
846       error->all(FLERR,
847 		 "Variable stats keyword cannot be used between runs");
848     compute_elapsed();
849     dvalue = bivalue;
850 
851   } else if (strcmp(word,"elaplong") == 0) {
852     if (update->runflag == 0)
853       error->all(FLERR,
854 		 "Variable stats keyword cannot be used between runs");
855     compute_elaplong();
856     dvalue = bivalue;
857 
858   } else if (strcmp(word,"dt") == 0) {
859     compute_dt();
860 
861   } else if (strcmp(word,"cpu") == 0) {
862     if (update->runflag == 0)
863       error->all(FLERR,
864 		 "Variable stats keyword cannot be used between runs");
865     compute_cpu();
866 
867   } else if (strcmp(word,"tpcpu") == 0) {
868     if (update->runflag == 0)
869       error->all(FLERR,
870 		 "Variable stats keyword cannot be used between runs");
871     compute_tpcpu();
872 
873   } else if (strcmp(word,"spcpu") == 0) {
874     if (update->runflag == 0)
875       error->all(FLERR,
876 		 "Variable stats keyword cannot be used between runs");
877     compute_spcpu();
878 
879   } else if (strcmp(word,"wall") == 0) {
880     compute_wall();
881 
882   } else if (strcmp(word,"np") == 0) {
883     compute_np();
884     dvalue = bivalue;
885   } else if (strcmp(word,"ntouch") == 0) {
886     compute_ntouch();
887     dvalue = bivalue;
888   } else if (strcmp(word,"ncomm") == 0) {
889     compute_ncomm();
890     dvalue = bivalue;
891   } else if (strcmp(word,"nbound") == 0) {
892     compute_nbound();
893     dvalue = bivalue;
894   } else if (strcmp(word,"nexit") == 0) {
895     compute_nexit();
896     dvalue = bivalue;
897   } else if (strcmp(word,"nscoll") == 0) {
898     compute_nscoll();
899     dvalue = bivalue;
900   } else if (strcmp(word,"nscheck") == 0) {
901     compute_nscheck();
902     dvalue = bivalue;
903   } else if (strcmp(word,"ncoll") == 0) {
904     compute_ncoll();
905     dvalue = bivalue;
906   } else if (strcmp(word,"nattempt") == 0) {
907     compute_nattempt();
908     dvalue = bivalue;
909   } else if (strcmp(word,"nreact") == 0) {
910     compute_nreact();
911     dvalue = bivalue;
912   } else if (strcmp(word,"nsreact") == 0) {
913     compute_nsreact();
914     dvalue = bivalue;
915   }
916 
917   else if (strcmp(word,"npave") == 0) compute_npave();
918   else if (strcmp(word,"ntouchave") == 0) compute_ntouchave();
919   else if (strcmp(word,"ncommave") == 0) compute_ncommave();
920   else if (strcmp(word,"nboundave") == 0) compute_nboundave();
921   else if (strcmp(word,"nexitave") == 0) compute_nexitave();
922   else if (strcmp(word,"nscollave") == 0) compute_nscollave();
923   else if (strcmp(word,"nscheckave") == 0) compute_nscheckave();
924   else if (strcmp(word,"ncollave") == 0) compute_ncollave();
925   else if (strcmp(word,"nattemptave") == 0) compute_nattemptave();
926   else if (strcmp(word,"nreactave") == 0) compute_nreactave();
927   else if (strcmp(word,"nsreactave") == 0) compute_nsreactave();
928 
929   else if (strcmp(word,"ngrid") == 0) {
930     compute_ngrid();
931     dvalue = bivalue;
932   }
933   else if (strcmp(word,"nsplit") == 0) compute_nsplit();
934   else if (strcmp(word,"maxlevel") == 0) compute_maxlevel();
935 
936   else if (strcmp(word,"vol") == 0) compute_vol();
937   else if (strcmp(word,"lx") == 0) compute_lx();
938   else if (strcmp(word,"ly") == 0) compute_ly();
939   else if (strcmp(word,"lz") == 0) compute_lz();
940 
941   else if (strcmp(word,"xlo") == 0) compute_xlo();
942   else if (strcmp(word,"xhi") == 0) compute_xhi();
943   else if (strcmp(word,"ylo") == 0) compute_ylo();
944   else if (strcmp(word,"yhi") == 0) compute_yhi();
945   else if (strcmp(word,"zlo") == 0) compute_zlo();
946   else if (strcmp(word,"zhi") == 0) compute_zhi();
947 
948   else return 1;
949 
950   *answer = dvalue;
951   return 0;
952 }
953 
954 /* ----------------------------------------------------------------------
955    extraction of SurfCollide, SurfReact, Compute, Fix, Variable results
956 ------------------------------------------------------------------------- */
957 
compute_surf_collide()958 void Stats::compute_surf_collide()
959 {
960   int m = field2index[ifield];
961   SurfCollide *one = sc[m];
962   dvalue = one->compute_vector(argindex1[ifield]-1);
963 }
964 
965 /* ---------------------------------------------------------------------- */
966 
compute_surf_react()967 void Stats::compute_surf_react()
968 {
969   int m = field2index[ifield];
970   SurfReact *one = sr[m];
971   dvalue = one->compute_vector(argindex1[ifield]-1);
972 }
973 
974 /* ---------------------------------------------------------------------- */
975 
compute_compute()976 void Stats::compute_compute()
977 {
978   int m = field2index[ifield];
979   Compute *compute = computes[m];
980 
981   if (compute_which[m] == SCALAR)
982     dvalue = compute->scalar;
983   else if (compute_which[m] == VECTOR)
984     dvalue = compute->vector[argindex1[ifield]-1];
985   else
986     dvalue = compute->array[argindex1[ifield]-1][argindex2[ifield]-1];
987 }
988 
989 /* ---------------------------------------------------------------------- */
990 
compute_fix()991 void Stats::compute_fix()
992 {
993   int m = field2index[ifield];
994   Fix *fix = fixes[m];
995 
996   if (argindex1[ifield] == 0)
997     dvalue = fix->compute_scalar();
998   else if (argindex2[ifield] == 0)
999     dvalue = fix->compute_vector(argindex1[ifield]-1);
1000   else
1001     dvalue = fix->compute_array(argindex1[ifield]-1,argindex2[ifield]-1);
1002 }
1003 
1004 /* ---------------------------------------------------------------------- */
1005 
compute_variable()1006 void Stats::compute_variable()
1007 {
1008   dvalue = input->variable->compute_equal(variables[field2index[ifield]]);
1009 }
1010 
1011 /* ----------------------------------------------------------------------
1012    one method for every keyword stats can output
1013    called by compute() or evaluate_keyword()
1014    compute will have already been called
1015    set ivalue/dvalue/bivalue if value is int/double/bigint
1016    customize a new keyword by adding a method
1017 ------------------------------------------------------------------------- */
1018 
compute_step()1019 void Stats::compute_step()
1020 {
1021   bivalue = update->ntimestep;
1022 }
1023 
1024 /* ---------------------------------------------------------------------- */
1025 
compute_elapsed()1026 void Stats::compute_elapsed()
1027 {
1028   bivalue = update->ntimestep - update->firststep;
1029 }
1030 
1031 /* ---------------------------------------------------------------------- */
1032 
compute_elaplong()1033 void Stats::compute_elaplong()
1034 {
1035   bivalue = update->ntimestep - update->beginstep;
1036 }
1037 
1038 /* ---------------------------------------------------------------------- */
1039 
compute_dt()1040 void Stats::compute_dt()
1041 {
1042   dvalue = update->dt;
1043 }
1044 
1045 /* ---------------------------------------------------------------------- */
1046 
compute_cpu()1047 void Stats::compute_cpu()
1048 {
1049   if (firststep == 0) dvalue = 0.0;
1050   else dvalue = timer->elapsed(TIME_LOOP);
1051 }
1052 
1053 /* ---------------------------------------------------------------------- */
1054 
compute_tpcpu()1055 void Stats::compute_tpcpu()
1056 {
1057   double new_cpu;
1058   double new_time = update->ntimestep * update->dt;
1059 
1060   if (firststep == 0) {
1061     new_cpu = 0.0;
1062     dvalue = 0.0;
1063   } else {
1064     new_cpu = timer->elapsed(TIME_LOOP);
1065     double cpu_diff = new_cpu - last_tpcpu;
1066     double time_diff = new_time - last_time;
1067     if (time_diff > 0.0 && cpu_diff > 0.0) dvalue = time_diff/cpu_diff;
1068     else dvalue = 0.0;
1069   }
1070 
1071   last_time = new_time;
1072   last_tpcpu = new_cpu;
1073 }
1074 
1075 /* ---------------------------------------------------------------------- */
1076 
compute_spcpu()1077 void Stats::compute_spcpu()
1078 {
1079   double new_cpu;
1080   int new_step = update->ntimestep;
1081 
1082   if (firststep == 0) {
1083     new_cpu = 0.0;
1084     dvalue = 0.0;
1085   } else {
1086     new_cpu = timer->elapsed(TIME_LOOP);
1087     double cpu_diff = new_cpu - last_spcpu;
1088     int step_diff = new_step - last_step;
1089     if (cpu_diff > 0.0) dvalue = step_diff/cpu_diff;
1090     else dvalue = 0.0;
1091   }
1092 
1093   last_step = new_step;
1094   last_spcpu = new_cpu;
1095 }
1096 
1097 /* ---------------------------------------------------------------------- */
1098 
compute_wall()1099 void Stats::compute_wall()
1100 {
1101   dvalue = MPI_Wtime() - wall0;
1102 }
1103 
1104 /* ---------------------------------------------------------------------- */
1105 
compute_np()1106 void Stats::compute_np()
1107 {
1108   bigint n = particle->nlocal;
1109   MPI_Allreduce(&n,&particle->nglobal,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1110   bivalue = particle->nglobal;
1111 }
1112 
1113 /* ---------------------------------------------------------------------- */
1114 
compute_ntouch()1115 void Stats::compute_ntouch()
1116 {
1117   bigint n = update->ntouch_one;
1118   MPI_Allreduce(&n,&bivalue,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1119 }
1120 
1121 /* ---------------------------------------------------------------------- */
1122 
compute_ncomm()1123 void Stats::compute_ncomm()
1124 {
1125   bigint n = update->ncomm_one;
1126   MPI_Allreduce(&n,&bivalue,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1127 }
1128 
1129 /* ---------------------------------------------------------------------- */
1130 
compute_nbound()1131 void Stats::compute_nbound()
1132 {
1133   bigint n = update->nboundary_one;
1134   MPI_Allreduce(&n,&bivalue,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1135 }
1136 
1137 /* ---------------------------------------------------------------------- */
1138 
compute_nexit()1139 void Stats::compute_nexit()
1140 {
1141   bigint n = update->nexit_one;
1142   MPI_Allreduce(&n,&bivalue,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1143 }
1144 
1145 /* ---------------------------------------------------------------------- */
1146 
compute_nscoll()1147 void Stats::compute_nscoll()
1148 {
1149   bigint n = update->nscollide_one;
1150   MPI_Allreduce(&n,&bivalue,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1151 }
1152 
1153 /* ---------------------------------------------------------------------- */
1154 
compute_nscheck()1155 void Stats::compute_nscheck()
1156 {
1157   bigint n = update->nscheck_one;
1158   MPI_Allreduce(&n,&bivalue,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1159 }
1160 
1161 /* ---------------------------------------------------------------------- */
1162 
compute_ncoll()1163 void Stats::compute_ncoll()
1164 {
1165   if (!collide) bivalue = 0;
1166   else {
1167     bigint n = collide->ncollide_one;
1168     MPI_Allreduce(&n,&bivalue,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1169   }
1170 }
1171 
1172 /* ---------------------------------------------------------------------- */
1173 
compute_nattempt()1174 void Stats::compute_nattempt()
1175 {
1176   if (!collide) bivalue = 0;
1177   else {
1178     bigint n = collide->nattempt_one;
1179     MPI_Allreduce(&n,&bivalue,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1180   }
1181 }
1182 
1183 /* ---------------------------------------------------------------------- */
1184 
compute_nreact()1185 void Stats::compute_nreact()
1186 {
1187   if (!collide) bivalue = 0;
1188   else {
1189     bigint n = collide->nreact_one;
1190     MPI_Allreduce(&n,&bivalue,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1191   }
1192 }
1193 
1194 /* ---------------------------------------------------------------------- */
1195 
compute_nsreact()1196 void Stats::compute_nsreact()
1197 {
1198   bigint n = surf->nreact_one;
1199   MPI_Allreduce(&n,&bivalue,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
1200 }
1201 
1202 /* ---------------------------------------------------------------------- */
1203 
compute_npave()1204 void Stats::compute_npave()
1205 {
1206   MPI_Allreduce(&update->nmove_running,&bivalue,1,MPI_SPARTA_BIGINT,
1207 		MPI_SUM,world);
1208   if (update->ntimestep == update->firststep) dvalue = 0.0;
1209   else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1210 }
1211 
1212 /* ---------------------------------------------------------------------- */
1213 
compute_ntouchave()1214 void Stats::compute_ntouchave()
1215 {
1216   MPI_Allreduce(&update->ntouch_running,&bivalue,1,MPI_SPARTA_BIGINT,
1217 		MPI_SUM,world);
1218   if (update->ntimestep == update->firststep) dvalue = 0.0;
1219   else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1220 }
1221 
1222 /* ---------------------------------------------------------------------- */
1223 
compute_ncommave()1224 void Stats::compute_ncommave()
1225 {
1226   MPI_Allreduce(&update->ncomm_running,&bivalue,1,MPI_SPARTA_BIGINT,
1227 		MPI_SUM,world);
1228   if (update->ntimestep == update->firststep) dvalue = 0.0;
1229   else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1230 }
1231 
1232 /* ---------------------------------------------------------------------- */
1233 
compute_nboundave()1234 void Stats::compute_nboundave()
1235 {
1236   MPI_Allreduce(&update->nboundary_running,&bivalue,1,MPI_SPARTA_BIGINT,
1237 		MPI_SUM,world);
1238   if (update->ntimestep == update->firststep) dvalue = 0.0;
1239   else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1240 }
1241 
1242 /* ---------------------------------------------------------------------- */
1243 
compute_nexitave()1244 void Stats::compute_nexitave()
1245 {
1246   MPI_Allreduce(&update->nexit_running,&bivalue,1,MPI_SPARTA_BIGINT,
1247 		MPI_SUM,world);
1248   if (update->ntimestep == update->firststep) dvalue = 0.0;
1249   else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1250 }
1251 
1252 /* ---------------------------------------------------------------------- */
1253 
compute_nscollave()1254 void Stats::compute_nscollave()
1255 {
1256   MPI_Allreduce(&update->nscollide_running,&bivalue,1,MPI_SPARTA_BIGINT,
1257 		MPI_SUM,world);
1258   if (update->ntimestep == update->firststep) dvalue = 0.0;
1259   else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1260 }
1261 
1262 /* ---------------------------------------------------------------------- */
1263 
compute_nscheckave()1264 void Stats::compute_nscheckave()
1265 {
1266   MPI_Allreduce(&update->nscheck_running,&bivalue,1,MPI_SPARTA_BIGINT,
1267 		MPI_SUM,world);
1268   if (update->ntimestep == update->firststep) dvalue = 0.0;
1269   else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1270 }
1271 
1272 /* ---------------------------------------------------------------------- */
1273 
compute_ncollave()1274 void Stats::compute_ncollave()
1275 {
1276   if (!collide) dvalue = 0.0;
1277   else {
1278     MPI_Allreduce(&collide->ncollide_running,&bivalue,1,MPI_SPARTA_BIGINT,
1279 		  MPI_SUM,world);
1280     if (update->ntimestep == update->firststep) dvalue = 0.0;
1281     else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1282   }
1283 }
1284 
1285 /* ---------------------------------------------------------------------- */
1286 
compute_nattemptave()1287 void Stats::compute_nattemptave()
1288 {
1289   if (!collide) dvalue = 0.0;
1290   else {
1291     MPI_Allreduce(&collide->nattempt_running,&bivalue,1,MPI_SPARTA_BIGINT,
1292 		  MPI_SUM,world);
1293     if (update->ntimestep == update->firststep) dvalue = 0.0;
1294     else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1295   }
1296 }
1297 
1298 /* ---------------------------------------------------------------------- */
1299 
compute_nreactave()1300 void Stats::compute_nreactave()
1301 {
1302   if (!collide) dvalue = 0.0;
1303   else {
1304     MPI_Allreduce(&collide->nreact_running,&bivalue,1,MPI_SPARTA_BIGINT,
1305 		  MPI_SUM,world);
1306     if (update->ntimestep == update->firststep) dvalue = 0.0;
1307     else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1308   }
1309 }
1310 
1311 /* ---------------------------------------------------------------------- */
1312 
compute_nsreactave()1313 void Stats::compute_nsreactave()
1314 {
1315   MPI_Allreduce(&surf->nreact_running,&bivalue,1,MPI_SPARTA_BIGINT,
1316                 MPI_SUM,world);
1317   if (update->ntimestep == update->firststep) dvalue = 0.0;
1318   else dvalue = 1.0*bivalue / (update->ntimestep - update->firststep);
1319 }
1320 
1321 /* ---------------------------------------------------------------------- */
1322 
compute_ngrid()1323 void Stats::compute_ngrid()
1324 {
1325   bivalue = grid->ncell;
1326 }
1327 
1328 /* ---------------------------------------------------------------------- */
1329 
compute_nsplit()1330 void Stats::compute_nsplit()
1331 {
1332   ivalue = grid->nsplit;
1333 }
1334 
1335 /* ---------------------------------------------------------------------- */
1336 
compute_maxlevel()1337 void Stats::compute_maxlevel()
1338 {
1339   ivalue = grid->maxlevel;
1340 }
1341 
1342 /* ---------------------------------------------------------------------- */
1343 
compute_vol()1344 void Stats::compute_vol()
1345 {
1346   if (domain->dimension == 3)
1347     dvalue = domain->xprd * domain->yprd * domain->zprd;
1348   else
1349     dvalue = domain->xprd * domain->yprd;
1350 }
1351 
1352 /* ---------------------------------------------------------------------- */
1353 
compute_lx()1354 void Stats::compute_lx()
1355 {
1356   dvalue = domain->xprd;
1357 }
1358 
1359 /* ---------------------------------------------------------------------- */
1360 
compute_ly()1361 void Stats::compute_ly()
1362 {
1363   dvalue = domain->yprd;
1364 }
1365 
1366 /* ---------------------------------------------------------------------- */
1367 
compute_lz()1368 void Stats::compute_lz()
1369 {
1370   dvalue = domain->zprd;
1371 }
1372 
1373 /* ---------------------------------------------------------------------- */
1374 
compute_xlo()1375 void Stats::compute_xlo()
1376 {
1377   dvalue = domain->boxlo[0];
1378 }
1379 
1380 /* ---------------------------------------------------------------------- */
1381 
compute_xhi()1382 void Stats::compute_xhi()
1383 {
1384   dvalue = domain->boxhi[0];
1385 }
1386 
1387 /* ---------------------------------------------------------------------- */
1388 
compute_ylo()1389 void Stats::compute_ylo()
1390 {
1391   dvalue = domain->boxlo[1];
1392 }
1393 
1394 /* ---------------------------------------------------------------------- */
1395 
compute_yhi()1396 void Stats::compute_yhi()
1397 {
1398   dvalue = domain->boxhi[1];
1399 }
1400 
1401 /* ---------------------------------------------------------------------- */
1402 
compute_zlo()1403 void Stats::compute_zlo()
1404 {
1405   dvalue = domain->boxlo[2];
1406 }
1407 
1408 /* ---------------------------------------------------------------------- */
1409 
compute_zhi()1410 void Stats::compute_zhi()
1411 {
1412   dvalue = domain->boxhi[2];
1413 }
1414