PLUTO
set_grid.c
Go to the documentation of this file.
1 /* ///////////////////////////////////////////////////////////////////// */
2 /*!
3  \file
4  \brief Set up the global and local grids in the three coordinate
5  directions.
6 
7  Collects functions for allocating memory and defining grid
8  coordinates.
9 
10  \author A. Mignone (mignone@ph.unito.it)
11  \date Aug 24, 2015
12 */
13 /* ///////////////////////////////////////////////////////////////////// */
14 #include "pluto.h"
15 
16 static void InitializeGrid (Runtime *, Grid *);
17 static void MakeGrid (int, Runtime *, double *, double *, double *);
18 static void stretch_fun (double, double *, double *, double *);
19 
20 /* ********************************************************************* */
21 void SetGrid (Runtime *rtime, Grid *GXYZ)
22 /*!
23  *
24  * \param [in] rtime pointer to a Runtime structure
25  * \param [out] GXYZ pointer to array of Grid structures
26  *
27  *********************************************************************** */
28 {
29  int i, idim;
30  int iL, iR, ngh;
31  char fname[512];
32  double *dx, *xlft, *xrgt;
33  Grid *G;
34  FILE *fg;
35  double xpatch_lft, xpatch_rgt;
36 
37  InitializeGrid(rtime, GXYZ);
38 
39  i = MAX(rtime->npoint[0], MAX(rtime->npoint[1], rtime->npoint[2]));
40  dx = ARRAY_1D(i + 2, double);
41  xlft = ARRAY_1D(i + 2, double);
42  xrgt = ARRAY_1D(i + 2, double);
43 
44  for (idim = 0; idim < 3; idim++) {
45 
46  G = GXYZ + idim;
47  ngh = G->nghost;
48 
49  iL = ngh - 1;
50  MakeGrid (idim, rtime, xlft, xrgt, dx);
51 /*MakeGrid (idim, rtime, G->xl_glob + iL, G->xr_glob + iL, G->dx_glob + iL); */
52 
53  /* ---- Assign values to grid structure members ---- */
54 
55  for (i = 1; i <= rtime->npoint[idim]; i++) {
56  G->dx_glob[i + ngh - 1] = dx[i];
57  G->xl_glob[i + ngh - 1] = xlft[i];
58  G->xr_glob[i + ngh - 1] = xrgt[i];
59  }
60 
61  iL = ngh;
62  iR = rtime->npoint[idim] + ngh - 1;
63  if (idim < DIMENSIONS){
64  G->xr_glob[iL - 1] = G->xl_glob[iL];
65  G->xl_glob[iR + 1] = G->xr_glob[iR];
66  }
67 /* ---- fill boundary values by copying adjacent cells ---- */
68 
69  for (i = 0; i < ngh; i++) {
70 
71  /* ---- left boundary ---- */
72 
73  G->dx_glob[i] = G->dx_glob[iL];
74  G->xl_glob[i] = G->xl_glob[iL] - (ngh - i)*G->dx_glob[iL];
75  G->xr_glob[i] = G->xl_glob[i] + G->dx_glob[iL];
76 
77  /* ---- right boundary ---- */
78 
79  G->dx_glob[iR + i + 1] = G->dx_glob[iR];
80  G->xl_glob[iR + i + 1] = G->xl_glob[iR] + (i + 1.0)*G->dx_glob[iR];
81  G->xr_glob[iR + i + 1] = G->xl_glob[iR] + (i + 2.0)*G->dx_glob[iR];
82  }
83 
84 /* ---- define geometrical cell center ---- */
85 
86  for (i = 0; i <= iR + ngh; i++) {
87  G->x_glob[i] = 0.5*(G->xl_glob[i] + G->xr_glob[i]);
88  }
89 
90 /* ---- define leftmost and rightmost domain extrema ---- */
91 
92  G->xi = G->xl_glob[iL];
93  G->xf = G->xr_glob[iR];
94 
95  g_domBeg[idim] = G->xl_glob[iL];
96  g_domEnd[idim] = G->xr_glob[iR];
97  }
98 
99 /* ---- free memory ---- */
100 
101  FreeArray1D(dx);
102  FreeArray1D(xlft);
103  FreeArray1D(xrgt);
104 
105 /* ---------------------------------------------------------------------
106  Write grid file
107  --------------------------------------------------------------------- */
108 
109  sprintf (fname,"%s/grid.out", rtime->output_dir);
110 #ifdef PLUTO3_Grid
111  fg = fopen(fname,"w");
112  for (idim = 0; idim < 3; idim++) {
113  G = GXYZ + idim;
114  ngh = G->nghost;
115  iL = G->gbeg;
116  iR = G->gend;
117 
118  fprintf (fg, "%d \n", iR - iL + 1);
119  for (i = iL; i <= iR; i++) {
120  fprintf (fg, " %d %12.6e %12.6e %12.6e %12.6e\n",
121  i-ngh+1, G->xl_glob[i], G->x_glob[i], G->xr_glob[i], G->dx_glob[i]);
122  }
123  }
124  fclose(fg);
125 #else
126  if (prank == 0){
127  time_t time_now;
128  time(&time_now);
129 
130  fg = fopen(fname,"w");
131  fprintf (fg, "# ******************************************************\n");
132  #ifdef CHOMBO
133  fprintf (fg, "# PLUTO-Chombo %s (Base) Grid File\n",PLUTO_VERSION);
134  #else
135  fprintf (fg, "# PLUTO %s Grid File\n",PLUTO_VERSION);
136  #endif
137 
138  fprintf (fg, "# Generated on %s",asctime(localtime(&time_now)));
139 
140 /* fprintf (fg, "# %s\n", getenv("PWD")); */
141 
142  fprintf (fg, "# \n");
143  fprintf (fg,"# DIMENSIONS: %d\n",DIMENSIONS);
144  #if GEOMETRY == CARTESIAN
145  fprintf (fg, "# GEOMETRY: CARTESIAN\n");
146  #elif GEOMETRY == CYLINDRICAL
147  fprintf (fg, "# GEOMETRY: CYLINDRICAL\n");
148  #elif GEOMETRY == POLAR
149  fprintf (fg, "# GEOMETRY: POLAR\n");
150  #elif GEOMETRY == SPHERICAL
151  fprintf (fg, "# GEOMETRY: SPHERICAL\n");
152  #endif
153  for (idim = 0; idim < DIMENSIONS; idim++){
154  fprintf (fg, "# X%d: [% f, % f], %d point(s), %d ghosts\n", idim+1,
155  g_domBeg[idim], g_domEnd[idim],
156  GXYZ[idim].np_int_glob, GXYZ[idim].nghost);
157  }
158  fprintf (fg, "# ******************************************************\n");
159  for (idim = 0; idim < 3; idim++) {
160  G = GXYZ + idim;
161  ngh = G->nghost;
162  iL = G->gbeg;
163  iR = G->gend;
164 
165  fprintf (fg, "%d \n", iR - iL + 1);
166  for (i = iL; i <= iR; i++) {
167  fprintf (fg, " %d %18.12e %18.12e\n",
168  i-ngh+1, G->xl_glob[i], G->xr_glob[i]);
169  }
170  }
171  fclose(fg);
172  }
173 #endif
174 
175 /* ---- define geometry factors, vol, area, etc... ---- */
176 
177  MakeGeometry(GXYZ);
178 
179 /* ----------------------------------------
180  print domain specifications
181  ---------------------------------------- */
182 
183  for (idim = 0; idim < DIMENSIONS; idim++){
184  print1 (" X%d: [% f, % f], %d point(s), %d ghosts\n", idim+1,
185  g_domBeg[idim], g_domEnd[idim],
186  GXYZ[idim].np_int_glob, GXYZ[idim].nghost);
187  }
188 }
189 
190 /* ********************************************************************* */
191 void FreeGrid (Grid *grid)
192 /*!
193  * Free array memory allocated previously.
194  *
195  *********************************************************************** */
196 {
197  int dir;
198 
199  for (dir = 0; dir < 3; dir++){
200  FreeArray1D(grid[dir].x);
201  FreeArray1D(grid[dir].xl);
202  FreeArray1D(grid[dir].xr);
203  FreeArray1D(grid[dir].dx);
204  FreeArray1D(grid[dir].xgc);
205  FreeArray1D(grid[dir].dV);
206  FreeArray1D(grid[dir].A-1);
207  FreeArray1D(grid[dir].r_1);
208  FreeArray1D(grid[dir].ct);
209  FreeArray1D(grid[dir].inv_dx);
210  FreeArray1D(grid[dir].inv_dxi);
211  }
212 }
213 
214 /* ********************************************************************* */
215 void InitializeGrid (Runtime *rtime, Grid *GXYZ)
216 /*!
217  * Allocate memory for grid arrays, set shortcut pointers
218  * for local grids.
219  *
220  *
221  *********************************************************************** */
222 {
223  int idim, ngh;
225  struct GRID *G;
226 
227  for (idim = 0; idim < 3; idim++) {
228 
229  G = GXYZ + idim;
230  if (GEOMETRY == CARTESIAN){
231  G->uniform = rtime->grid_is_uniform[idim];
232  }else{
233  G->uniform = 0;
234  }
235 
236 /* ----------------------------------------------------------------
237  Dimensions that are not used will have 1 grid point, no bound
238  ---------------------------------------------------------------- */
239 
240  if (idim >= DIMENSIONS) {
241  G->np_int = G->np_tot = G->np_int_glob = G->np_tot_glob = 1;
242  ngh = G->nghost = G->beg = G->end = G->gbeg = G->gend =
243  G->lbeg = G->lend = 0;
244  G->nproc = 1;
245  } else {
246  ngh = G->nghost;
247  }
248 
249  np_tot_glob = G->np_tot_glob;
250  np_int_glob = G->np_int_glob;
251  np_tot = G->np_tot;
252  np_int = G->np_int;
253 
254 /* -----------------------------------------------------------
255  Memory allocation.
256  ----------------------------------------------------------- */
257 
258  G->x_glob = ARRAY_1D(np_tot_glob, double);
259  G->xr_glob = ARRAY_1D(np_tot_glob, double);
260  G->xl_glob = ARRAY_1D(np_tot_glob, double);
261  G->dx_glob = ARRAY_1D(np_tot_glob, double);
262 
263 /* ---- define shortcuts for local grids ---- */
264 
265  G->x = G->x_glob + G->beg - ngh;
266  G->xr = G->xr_glob + G->beg - ngh;
267  G->xl = G->xl_glob + G->beg - ngh;
268  G->dx = G->dx_glob + G->beg - ngh;
269  }
270 }
271 
272 /* ********************************************************************* */
273 void MakeGrid (int idim, Runtime *rtime, double *xlft, double *xrgt, double *dx)
274 /*!
275  *
276  * Build grid nodes as defined by pluto.ini.
277  *
278  * Options are:
279  *
280  * 'u' = uniform grid, simply defined as
281  *
282  * dx = (xR - xL)/npoint, xleft(i) = xl + i*dx, xright(i) = xleft(i) + dx
283  *
284  * 's' = stretched grid; solve
285  *
286  * dx*( 1 + r + r^2 + r^3 + ... r^(N-1)) = xR - xL
287  *
288  * in the stretching ratio r, provided dx, N, xR and xL are known.
289  * dx is taken from the closest uniform grid.
290  *
291  * 'l+' = logarithmic grid, mesh size increases with x; it is defined as
292  *
293  *
294  * x + |xL| - xL
295  * y = Log[ ------------- ] , with uniform spacing y(i+1/2) - y(i-1/2) = dy
296  * |xL|
297  *
298  * dy = (yR - yL)/N and dx(i) becomes
299  * dx(i) = (x(i-) + fabs(xL) - xL)*(10^dy - 1.0);
300  *
301  * NOTE: xR must be positive and xL can take any value different from 0
302  *
303  * 'l-' = logarithmic grid, mesh size decreases with x; it is defined as
304  *
305  *
306  * xR + |xL| - x
307  * y = Log[ -------------- ] , with uniform spacing y(i+1/2) - y(i-1/2) = dy
308  * |xL|
309  *
310  * dy = -(yR - yL)/N
311  * dx(i) = (x(i-) - fabs(xL) - xR)*(10^dy - 1.0);
312  *
313  *********************************************************************** */
314 #define MAX_ITER 50
315 {
316  int i, iseg, n, nstart, nseg, npoint;
317  int log_inc, log_dec, next_seg_is;
318  int done_with_segment[16], all_segments_done;
319  int iL, iR;
320  int i_patch_lft[16], i_patch_rgt[16];
321  double x_patch_lft[16], x_patch_rgt[16];
322  double xR, xL;
323  double par[4];
324  double dalpha, alpha, f, df;
325  double dy;
326 
327  nseg = rtime->npatch[idim];
328 
329 /* ---------------------------------------------------
330  for each patch, find the leftmost and
331  rightmost indexes ilft and irgt
332  --------------------------------------------------- */
333 
334  i_patch_lft[1] = 1;
335  i_patch_rgt[1] = rtime->patch_npoint[idim][1];
336  x_patch_lft[1] = rtime->patch_left_node[idim][1];
337  x_patch_rgt[1] = rtime->patch_left_node[idim][2];
338 
339  for (iseg = 2; iseg <= nseg; iseg++) {
340  i_patch_lft[iseg] = i_patch_rgt[iseg - 1] + 1;
341  i_patch_rgt[iseg] = i_patch_lft[iseg] + rtime->patch_npoint[idim][iseg] - 1;
342  x_patch_lft[iseg] = rtime->patch_left_node[idim][iseg];
343  x_patch_rgt[iseg] = rtime->patch_left_node[idim][iseg + 1];
344  }
345 
346  done_with_segment[0] = 0;
347  done_with_segment[nseg + 1] = 0;
348 
349  for (iseg = 1; iseg <= nseg; iseg++) {
350 
351  done_with_segment[iseg] = 0;
352  xL = x_patch_lft[iseg];
353  xR = x_patch_rgt[iseg];
354  iL = i_patch_lft[iseg];
355  iR = i_patch_rgt[iseg];
356 
357  npoint = rtime->patch_npoint[idim][iseg];
358 
359 /* ---- first process only uniform or logarithmic grids ---- */
360 
361  if (rtime->patch_type[idim][iseg] == UNIFORM_GRID) {
362 
363  for (i = iL; i <= iR; i++) {
364  dx[i] = (xR - xL)/(double)(npoint);
365  xlft[i] = xL + (double)(i - iL)*dx[i];
366  xrgt[i] = xlft[i] + dx[i];
367  }
368  done_with_segment[iseg] = 1;
369 
370  } else if ( rtime->patch_type[idim][iseg] == LOGARITHMIC_INC_GRID ||
371  rtime->patch_type[idim][iseg] == LOGARITHMIC_DEC_GRID) {
372 
373  log_inc = rtime->patch_type[idim][iseg] == LOGARITHMIC_INC_GRID;
374  log_dec = rtime->patch_type[idim][iseg] == LOGARITHMIC_DEC_GRID;
375 
376  dy = log10( (xR + fabs(xL) - xL)/fabs(xL));
377  dy /= (double) (npoint);
378 
379  if (log_dec) dy *= -1.0;
380 
381  xlft[iL] = xL;
382 
383  for (i = iL; i <= iR; i++) {
384  if (log_inc) {
385  dx[i] = (xlft[i] + fabs(xL) - xL)*(pow(10.0,dy) - 1.0);
386  }else{
387  dx[i] = (xlft[i] - fabs(xL) - xR)*(pow(10.0,dy) - 1.0);
388  }
389  xrgt[i] = xlft[i] + dx[i];
390  xlft[i + 1] = xrgt[i];
391  }
392  done_with_segment[iseg] = 1;
393 
394  } else {
395  continue;
396  }
397  }
398 
399  all_segments_done = 1;
400  for (iseg = 1; iseg <= nseg; iseg++) {
401  all_segments_done = all_segments_done && done_with_segment[iseg];
402  }
403 
404 /* --------------------------------------------------------------
405  now do stretched grids ...
406  -------------------------------------------------------------- */
407 
408  while (!all_segments_done) { /* loop until all segments are processed... */
409 
410 /* ---------------------------------------------------
411  scan the entire grid and process only those
412  segments close to a segment that is already done
413  (i.e. done_with_segment[iseg] = 1)
414  --------------------------------------------------- */
415 
416  for (iseg = 1; iseg <= nseg; iseg++) {
417 
418  xL = x_patch_lft[iseg];
419  xR = x_patch_rgt[iseg];
420  iL = i_patch_lft[iseg];
421  iR = i_patch_rgt[iseg];
422  npoint = rtime->patch_npoint[idim][iseg];
423 
424 /* -----------------------------------------------------
425  now find whether the segment to right (iseg+1) or
426  to the left (iseg-1) is already done;
427  ----------------------------------------------------- */
428 
429  next_seg_is = 0;
430  if (done_with_segment[iseg + 1]) next_seg_is = iseg + 1;
431  if (done_with_segment[iseg - 1]) next_seg_is = iseg - 1;
432 
433 /* -----------------------------------------------------
434  if none of them is done, skip this segment and
435  move forward, otherwise process segment iseg,
436  otherwise process the grid
437  ----------------------------------------------------- */
438 
439  if (rtime->patch_type[idim][iseg] == STRETCHED_GRID &&
440  next_seg_is != 0 && !done_with_segment[iseg]) {
441 
442 /* ----------------------------------------------------------
443  nstart is:
444 
445  * the rightmost point if the next grid is done,
446  i.e., next_seg_is > iseg;
447  * the leftmost point if the previous grid is done,
448  i.e., next_seg_is < iseg;
449  ---------------------------------------------------------- */
450 
451  nstart = (next_seg_is > iseg ? iR:iL);
452 
453  alpha = 1.01; /* provide a first guess */
454 
455 /* --------------------------------------------------------------------
456  par[0] : contains the number of points
457  par[1] : Ratio L/dx : the first dx is the one that belongs to the
458  previous (iseg>1) or next (iseg = 1) uniform segment
459  --------------------------------------------------------------------- */
460 
461  par[0] = (double)npoint; /* Number of points */
462 
463  if (next_seg_is > iseg)
464  par[1] = (xR - xL)/dx[nstart + 1];
465  else
466  par[1] = (xR - xL)/dx[nstart - 1];
467 
468 /* ---- Use NEWTON algorithm to find the stretch factor ---- */
469 
470  for (n = 0; n <= MAX_ITER; n++) {
471  if (n == MAX_ITER) {
472  print1 ("Too many iterations during grid (%d) generation!\n",idim);
473  QUIT_PLUTO(1);
474  }
475 
476  stretch_fun (alpha, par, &f, &df);
477 
478  dalpha = f / df;
479  alpha -= dalpha;
480 
481  if (fabs (dalpha) < 1.e-14*alpha)
482  break;
483  }
484 
485  if (alpha > 1.2) {
486  print1 (" ! WARNING: alpha=%12.6e > 1.05 , dimension %d\n", alpha, idim);
487  print1 (" ! While stretching segment %d\n", iseg);
488  QUIT_PLUTO(1);
489  }
490 
491  print1 (" - Stretched grid on dim %d (ratio: %f)\n",idim, alpha);
492 
493  if (next_seg_is > iseg) {
494  xrgt[nstart] = xR;
495  for (i = nstart; i >= iL; i--) {
496  dx[i] = pow (alpha, nstart - i + 1)*dx[nstart + 1];
497  xlft[i] = xrgt[i] - dx[i];
498  xrgt[i - 1] = xlft[i];
499  }
500  } else {
501  xlft[nstart] = xL;
502  for (i = nstart; i <= iR; i++) {
503  dx[i] = pow (alpha, i - nstart + 1)*dx[nstart - 1];
504  xlft[i + 1] = xrgt[i] = xlft[i] + dx[i];
505  }
506  }
507 
508 /* ---- ok, this segment has been processed ---- */
509 
510  done_with_segment[iseg] = 1;
511 
512 /* ---- check whether there are other segments to process ---- */
513 
514  all_segments_done = 1;
515  for (i = 1; i <= nseg; i++)
516  all_segments_done = all_segments_done && done_with_segment[i];
517 
518 /* ---- exit from segment loop (iseg) and rescan the grid from the beginning ---- */
519 
520  break;
521  }
522  }
523  }
524 }
525 #undef MAX_ITER
526 
527 /* ############################################################## */
528 void stretch_fun (double x, double *par, double *f, double *dfdx)
529 /*
530  #
531  # S
532  #
533  ################################################################ */
534 {
535  double xN, L_dx, scrh;
536 
537  xN = par[0];
538  L_dx = par[1];
539 
540  scrh = x*(pow (x, xN) - 1.0)/(x - 1.0);
541  *f = log (scrh) - log (L_dx);
542  *dfdx = 1.0/x + xN*pow (x, xN - 1.0)/(pow (x, xN) - 1.0) -
543  1.0/(x - 1.0);
544 }
void FreeArray1D(void *v)
Definition: arrays.c:37
static double alpha
Definition: init.c:3
#define MAX(a, b)
Definition: macros.h:101
void MakeGeometry(Grid *)
Definition: set_geometry.c:16
int end
Global end index for the local array.
Definition: structs.h:116
static int n
Definition: analysis.c:3
int npoint[3]
Global number of zones in the interior domain.
Definition: structs.h:256
void print1(const char *fmt,...)
Definition: amrPluto.cpp:511
static void InitializeGrid(Runtime *, Grid *)
Definition: set_grid.c:215
char output_dir[256]
The name of the output directory (output_dir for static PLUTO, Output_dir for PLUTO-Chombo) ...
Definition: structs.h:269
double * xr
Definition: structs.h:81
int np_int_glob
Total number of points in the global domain (boundaries excluded).
Definition: structs.h:98
#define LOGARITHMIC_INC_GRID
Definition: pluto.h:40
tuple scrh
Definition: configure.py:200
#define UNIFORM_GRID
Definition: pluto.h:38
double * x_glob
Cell geometrical central points.
Definition: structs.h:80
double * dx_glob
Cell size.
Definition: structs.h:83
int lend
Local end index for the local array.
Definition: structs.h:118
int lbeg
Local start index for the local array.
Definition: structs.h:117
#define PLUTO_VERSION
Definition: pluto.h:16
double * dV
Cell volume.
Definition: structs.h:86
double * xr_glob
Cell right interface.
Definition: structs.h:81
double xf
Leftmost and rightmost point in the local domain.
Definition: structs.h:79
int gend
Global end index for the global array.
Definition: structs.h:114
int uniform
Definition: structs.h:119
double * dx
Definition: structs.h:83
int prank
Processor rank.
Definition: globals.h:33
int grid_is_uniform[3]
Definition: structs.h:259
int np_tot_glob
Total number of points in the global domain (boundaries included).
Definition: structs.h:96
#define MAX_ITER
#define LOGARITHMIC_DEC_GRID
Definition: pluto.h:41
double * inv_dx
Definition: structs.h:90
double patch_left_node[5][16]
Definition: structs.h:273
double * xl
Definition: structs.h:82
int beg
Global start index for the local array.
Definition: structs.h:115
static void MakeGrid(int, Runtime *, double *, double *, double *)
Definition: set_grid.c:273
Definition: structs.h:78
int gbeg
Global start index for the global array.
Definition: structs.h:113
void SetGrid(Runtime *rtime, Grid *GXYZ)
Definition: set_grid.c:21
int nghost
Number of ghost zones.
Definition: structs.h:104
double * xl_glob
Cell left interface.
Definition: structs.h:82
double g_domBeg[3]
Lower limits of the computational domain.
Definition: globals.h:125
double * x
Definition: structs.h:80
#define STRETCHED_GRID
Definition: pluto.h:39
double * xgc
Cell volumetric centroid (!= x when geometry != CARTESIAN).
Definition: structs.h:84
#define GEOMETRY
Definition: definitions_01.h:4
double * r_1
Geometrical factor 1/r.
Definition: structs.h:88
#define CARTESIAN
Definition: pluto.h:33
PLUTO main header file.
#define ARRAY_1D(nx, type)
Definition: prototypes.h:170
int i
Definition: analysis.c:2
int patch_npoint[5][16]
Definition: structs.h:261
double * ct
Geometrical factor cot(theta).
Definition: structs.h:89
int patch_type[5][16]
Definition: structs.h:262
int npatch[5]
The number of grid patches.
Definition: structs.h:260
int np_tot
Total number of points in the local domain (boundaries included).
Definition: structs.h:100
static void stretch_fun(double, double *, double *, double *)
Definition: set_grid.c:528
double g_domEnd[3]
Upper limits of the computational domain.
Definition: globals.h:126
void FreeGrid(Grid *grid)
Definition: set_grid.c:191
double * inv_dxi
inverse of the distance between the center of two cells, inv_dxi = .
Definition: structs.h:91
#define QUIT_PLUTO(e_code)
Definition: macros.h:125
int nproc
number of processors for this grid.
Definition: structs.h:120
double xi
Definition: structs.h:79
double * A
Right interface area, A[i] = .
Definition: structs.h:87
#define DIMENSIONS
Definition: definitions_01.h:2
int np_int
Total number of points in the local domain (boundaries excluded).
Definition: structs.h:102