HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
VlasovInitialize.c File Reference

Initialize the Vlasov module. More...

#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <bandedmatrix.h>
#include <physicalmodels/vlasov.h>
#include <mpivars.h>
#include <hypar.h>

Go to the source code of this file.

Functions

double VlasovComputeCFL (void *, void *, double, double)
 
int VlasovAdvection (double *, double *, int, void *, double)
 
int VlasovUpwind (double *, double *, double *, double *, double *, double *, int, void *, double)
 
int VlasovWriteEFieldAndPotential (void *, void *, double)
 
int VlasovPreStep (double *, void *, void *, double)
 
int VlasovPostStage (double *, void *, void *, double)
 
int VlasovEField (double *, void *, double)
 
int VlasovInitialize (void *s, void *m)
 

Detailed Description

Initialize the Vlasov module.

Author
John Loffeld

Definition in file VlasovInitialize.c.

Function Documentation

◆ VlasovComputeCFL()

double VlasovComputeCFL ( void *  s,
void *  m,
double  dt,
double  t 
)

Compute the maximum CFL number

Computes the maximum CFL number over the domain. Note that the CFL is computed over the local domain on this processor only.

Parameters
sSolver object of type HyPar
mMPI object of type MPIVariables
dtTime step size for which to compute the CFL
tTime

Definition at line 19 of file VlasovComputeCFL.c.

24 {
25  HyPar *solver = (HyPar*) s;
26  Vlasov *params = (Vlasov*) solver->physics;
27 
28  int ndims = solver->ndims;
29  int ghosts = solver->ghosts;
30  int *dim = solver->dim_local;
31  double *u = solver->u;
32 
33  double max_cfl = 0;
34  int done = 0;
35  int index[ndims];
36  _ArraySetValue_(index,ndims,0);
37 
38  while (!done) {
39 
40  for (int dir=0; dir<ndims; dir++) {
41  double dxinv;
42  _GetCoordinate_(dir,index[dir],dim,ghosts,solver->dxinv,dxinv);
43  double eig = VlasovAdvectionCoeff(index, dir, solver);
44  double local_cfl = eig*dt*dxinv;
45  if (local_cfl > max_cfl) max_cfl = local_cfl;
46  }
47 
48  _ArrayIncrementIndex_(ndims,dim,index,done);
49  }
50 
51  return(max_cfl);
52 }
Definition: vlasov.h:57
double * u
Definition: hypar.h:116
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
#define _ArrayIncrementIndex_(N, imax, i, done)
#define _GetCoordinate_(dir, i, dim, ghosts, x, coord)
Definition: basic.h:31
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
double VlasovAdvectionCoeff(int *, int, void *)
double * dxinv
Definition: hypar.h:110

◆ VlasovAdvection()

int VlasovAdvection ( double *  f,
double *  u,
int  dir,
void *  s,
double  t 
)

Compute the advection term

Compute the hyperbolic flux over the local domain.

\begin{equation} {\bf F}\left({\bf u}(x,v)\right) = c {\bf u}(x,v), \end{equation}

where the advection coefficient \(c\) is computed by VlasovAdvectionCoeff().

Parameters
fArray to hold the computed flux (same size and layout as u)
uArray containing the conserved solution
dirSpatial dimension
sSolver object of type HyPar
tCurrent time

Definition at line 22 of file VlasovAdvection.c.

28 {
29 
30  HyPar *solver = (HyPar*) s;
31  Vlasov *param = (Vlasov*) solver->physics;
32 
33  int* dim = solver->dim_local;
34  int ghosts = solver->ghosts;
35  int ndims = solver->ndims;
36 
37  // set bounds for array index to include ghost points
38  int bounds[ndims];
39  _ArrayCopy1D_(dim,bounds,ndims);
40  for (int i=0; i<ndims; i++) bounds[i] += 2*ghosts;
41 
42  // set offset such that index is compatible with ghost point arrangement
43  int offset[ndims];
44  _ArraySetValue_(offset,ndims,-ghosts);
45 
46  int done = 0;
47  int index_wog[ndims];
48  int index[ndims]; _ArraySetValue_(index,ndims,0);
49  while (!done) {
50 
51  _ArrayCopy1D_(index, index_wog, ndims);
52  for (int i = 0; i < ndims; i++) index_wog[i] -= ghosts;
53  double adv_coeff = VlasovAdvectionCoeff(index_wog, dir, solver);
54 
55  int p; _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,p);
56  f[p] = adv_coeff * u[p];
57  _ArrayIncrementIndex_(ndims,bounds,index,done);
58 
59  }
60 
61  return 0;
62 }
Definition: vlasov.h:57
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
#define _ArrayIncrementIndex_(N, imax, i, done)
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
#define _ArrayCopy1D_(x, y, size)
double VlasovAdvectionCoeff(int *, int, void *)

◆ VlasovUpwind()

int VlasovUpwind ( double *  fI,
double *  fL,
double *  fR,
double *  uL,
double *  uR,
double *  u,
int  dir,
void *  s,
double  t 
)

Compute the upwind flux at interfaces

Upwinding scheme for the Vlasov equations

Parameters
fIComputed upwind interface flux
fLLeft-biased reconstructed interface flux
fRRight-biased reconstructed interface flux
uLLeft-biased reconstructed interface solution
uRRight-biased reconstructed interface solution
uCell-centered solution
dirSpatial dimension
sSolver object of type HyPar
tCurrent solution time

Definition at line 17 of file VlasovUpwind.c.

27 {
28  HyPar *solver = (HyPar*) s;
29  Vlasov *param = (Vlasov*) solver->physics;
30  int done;
31 
32  int ndims = solver->ndims;
33  int *dim = solver->dim_local;
34  int ghosts = solver->ghosts;
35 
36  int index_outer[ndims],
37  index_inter[ndims],
38  bounds_outer[ndims],
39  bounds_inter[ndims];
40  _ArrayCopy1D_(dim,bounds_outer,ndims); bounds_outer[dir] = 1;
41  _ArrayCopy1D_(dim,bounds_inter,ndims); bounds_inter[dir] += 1;
42 
43  done = 0; _ArraySetValue_(index_outer,ndims,0);
44  while (!done) {
45  _ArrayCopy1D_(index_outer,index_inter,ndims);
46  for (index_inter[dir] = 0; index_inter[dir] < bounds_inter[dir]; index_inter[dir]++) {
47 
48  int p; _ArrayIndex1D_(ndims,bounds_inter,index_inter,0,p);
49  int indexL[ndims]; _ArrayCopy1D_(index_inter,indexL,ndims); indexL[dir]--;
50  int indexR[ndims]; _ArrayCopy1D_(index_inter,indexR,ndims);
51  int idxL[ndims], idxR[ndims];
52 
53  int neig = 2+4*(ndims-1);
54  double eig[neig];
55  int count = 0;
56  eig[count] = VlasovAdvectionCoeff(indexL, dir, solver); count++;
57  eig[count] = VlasovAdvectionCoeff(indexR, dir, solver); count++;
58  for (int tdir = 0; tdir < ndims; tdir++ ) {
59  if (tdir != dir) {
60  _ArrayCopy1D_(indexL, idxL, ndims); idxL[tdir]--;
61  _ArrayCopy1D_(indexR, idxR, ndims); idxR[tdir]--;
62  eig[count] = VlasovAdvectionCoeff(idxL, dir, solver); count++;
63  eig[count] = VlasovAdvectionCoeff(idxR, dir, solver); count++;
64  _ArrayCopy1D_(indexL, idxL, ndims); idxL[tdir]++;
65  _ArrayCopy1D_(indexR, idxR, ndims); idxR[tdir]++;
66  eig[count] = VlasovAdvectionCoeff(idxL, dir, solver); count++;
67  eig[count] = VlasovAdvectionCoeff(idxR, dir, solver); count++;
68  }
69  }
70  if (count != neig) {
71  fprintf(stderr, "Error in VlasovUpwind(): count != neig for dir=%d\n",dir);
72  return 1;
73  }
74 
75  int all_positive = 1, all_negative = 1;
76  double maxabs_eig = 0.0;
77  for (int n = 0; n<neig; n++) {
78  if (eig[n] <= 0) all_positive = 0;
79  if (eig[n] >= 0) all_negative = 0;
80  if (absolute(eig[n]) > maxabs_eig) {
81  maxabs_eig = absolute(eig[n]);
82  }
83  }
84 
85  if (all_positive) {
86  fI[p] = fL[p];
87  } else if (all_negative) {
88  fI[p] = fR[p];
89  } else {
90  fI[p] = 0.5 * (fL[p] + fR[p] - maxabs_eig * (uR[p] - uL[p]));
91  }
92 
93  }
94  _ArrayIncrementIndex_(ndims,bounds_outer,index_outer,done);
95 
96  }
97 
98  return 0;
99 }
Definition: vlasov.h:57
#define absolute(a)
Definition: math_ops.h:32
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _ArrayIndex1D_(N, imax, i, ghost, index)
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
#define _ArrayIncrementIndex_(N, imax, i, done)
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
double VlasovAdvectionCoeff(int *, int, void *)
#define _ArrayCopy1D_(x, y, size)

◆ VlasovWriteEFieldAndPotential()

int VlasovWriteEFieldAndPotential ( void *  s,
void *  m,
double  a_t 
)

Write E-field and potential to file

Write out the electric field and potential to file

Parameters
sSolver object of type HyPar
mMPI object of type MPIVariables
a_tCurrent simulation time

Definition at line 23 of file VlasovWriteEFieldAndPotential.c.

26 {
27  HyPar *solver = (HyPar*) s;
28  MPIVariables *mpi = (MPIVariables*) m;
29  Vlasov *param = (Vlasov*) solver->physics;
30 
31  {
32  char fname_root[_MAX_STRING_SIZE_] = "efield";
33  VlasovWriteSpatialField( solver, mpi, param->e_field, fname_root );
34  if (!strcmp(solver->plot_solution,"yes")) {
35  VlasovPlotSpatialField( solver, mpi, param->e_field, a_t, fname_root );
36  }
37  }
38 
39  if (param->self_consistent_electric_field) {
40  char fname_root[_MAX_STRING_SIZE_] = "potential";
41  VlasovWriteSpatialField( solver, mpi, param->potential, fname_root );
42  if (!strcmp(solver->plot_solution,"yes")) {
43  VlasovPlotSpatialField( solver, mpi, param->potential, a_t, fname_root );
44  }
45  }
46 
47  return 0;
48 }
Definition: vlasov.h:57
char plot_solution[_MAX_STRING_SIZE_]
Definition: hypar.h:194
double * potential
Definition: vlasov.h:84
int self_consistent_electric_field
Definition: vlasov.h:60
double * e_field
Definition: vlasov.h:81
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _MAX_STRING_SIZE_
Definition: basic.h:14
int VlasovPlotSpatialField(void *, void *, double *, double, char *)
int VlasovWriteSpatialField(void *, void *, double *, char *)
void * physics
Definition: hypar.h:266
Structure of MPI-related variables.

◆ VlasovPreStep()

int VlasovPreStep ( double *  u,
void *  s,
void *  m,
double  waqt 
)

Vlasov-specific function called at the beginning of each time-step: Calls the function to set the electric field

Parameters
uSolution (conserved variables)
sSolver object of type HyPar
mMPI object of type MPIVariables
waqtCurrent solution time

Definition at line 14 of file VlasovFunctions.c.

19 {
20  HyPar *solver = (HyPar*) s;
21  Vlasov *param = (Vlasov*) solver->physics;
22 
23  int ierr = VlasovEField(u, solver, waqt);
24  if (ierr) return ierr;
25 
26  return 0;
27 }
Definition: vlasov.h:57
int VlasovEField(double *, void *, double)
Definition: VlasovEField.c:258
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
void * physics
Definition: hypar.h:266

◆ VlasovPostStage()

int VlasovPostStage ( double *  u,
void *  s,
void *  m,
double  waqt 
)

Vlasov-specific function called at the end of each stage in a multistage time integrator: Calls the function to set the electric field

Parameters
uSolution (conserved variables)
sSolver object of type HyPar
mMPI object of type MPIVariables
waqtCurrent solution time

Definition at line 33 of file VlasovFunctions.c.

38 {
39  HyPar *solver = (HyPar*) s;
40  Vlasov *param = (Vlasov*) solver->physics;
41 
42  int ierr = VlasovEField(u, solver, waqt);
43  if (ierr) return ierr;
44 
45  return 0;
46 }
Definition: vlasov.h:57
int VlasovEField(double *, void *, double)
Definition: VlasovEField.c:258
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
void * physics
Definition: hypar.h:266

◆ VlasovEField()

int VlasovEField ( double *  u,
void *  s,
double  t 
)

Compute electric field.

Parameters
uConserved solution
sSolver object of type HyPar
tCurrent time

Definition at line 258 of file VlasovEField.c.

262 {
263  HyPar *solver = (HyPar*) s;
264  Vlasov *param = (Vlasov*) solver->physics;
265 
266  int ierr;
267 
268  if (param->self_consistent_electric_field) {
269  ierr = SetEFieldSelfConsistent(u, s, t);
270  } else {
271  ierr = SetEFieldPrescribed(u, s, t);
272  }
273 
274  return ierr;
275 }
Definition: vlasov.h:57
int self_consistent_electric_field
Definition: vlasov.h:60
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
static int SetEFieldSelfConsistent(double *u, void *s, double t)
Definition: VlasovEField.c:92
void * physics
Definition: hypar.h:266
static int SetEFieldPrescribed(double *u, void *s, double t)
Definition: VlasovEField.c:45

◆ VlasovInitialize()

int VlasovInitialize ( void *  s,
void *  m 
)

Initialize the Vlasov physics module - allocate and set physics-related parameters, read physics-related inputs from file, and set the physics-related function pointers in HyPar

Parameters
sSolver object of type HyPar
mObject of type MPIVariables containing MPI-related info

Definition at line 43 of file VlasovInitialize.c.

46 {
47  HyPar *solver = (HyPar*) s;
48  MPIVariables *mpi = (MPIVariables*) m;
49  Vlasov *physics = (Vlasov*) solver->physics;
50 
51  int *dim_global = solver->dim_global;
52  int *dim_local = solver->dim_local;
53  int ghosts = solver->ghosts;
54 
55  if (solver->nvars != _MODEL_NVARS_) {
56  if (!mpi->rank) {
57  fprintf(stderr,"Error in VlasovInitialize(): nvars has to be %d.\n",_MODEL_NVARS_);
58  }
59  return(1);
60  }
61  if (solver->ndims != _MODEL_NDIMS_) {
62  if (!mpi->rank) {
63  fprintf(stderr,"Error in VlasovInitialize(): ndims has to be %d.\n",_MODEL_NDIMS_);
64  }
65  return(1);
66  }
67 
68  /* default is prescribed electric field */
69  physics->self_consistent_electric_field = 0;
70  physics->ndims_x = 1;
71  physics->ndims_v = 1;
72 
73  /* reading physical model specific inputs - all processes */
74  if (!mpi->rank) {
75  FILE *in;
76  in = fopen("physics.inp","r");
77  if (in) {
78  printf("Reading physical model inputs from file \"physics.inp\".\n");
79  char word[_MAX_STRING_SIZE_];
80  int ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
81  if (!strcmp(word, "begin")){
82  while (strcmp(word, "end")){
83  int ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
84  if (!strcmp(word, "self_consistent_electric_field")) {
85  /* read whether electric field is self-consistent or prescribed */
86  int ferr = fscanf(in,"%d", &physics->self_consistent_electric_field);
87  if (ferr != 1) return(1);
88  } else if (!strcmp(word, "x_ndims")) {
89  /* read number of spatial dimensions */
90  int ferr = fscanf(in,"%d", &physics->ndims_x);
91  if (ferr != 1) return(1);
92  } else if (!strcmp(word, "v_ndims")) {
93  /* read number of velocity dimensions */
94  int ferr = fscanf(in,"%d", &physics->ndims_v);
95  if (ferr != 1) return(1);
96  } else if (strcmp(word,"end")) {
97  char useless[_MAX_STRING_SIZE_];
98  int ferr = fscanf(in,"%s",useless); if (ferr != 1) return(ferr);
99  printf("Warning: keyword %s in file \"physics.inp\" with value %s not ",
100  word, useless);
101  printf("recognized or extraneous. Ignoring.\n");
102  }
103  }
104  } else {
105  fprintf(stderr,"Error: Illegal format in file \"physics.inp\".\n");
106  return(1);
107  }
108  }
109  fclose(in);
110  }
111 
112  if ((physics->ndims_x+physics->ndims_v) != solver->ndims) {
113  if (!mpi->rank) {
114  fprintf(stderr,"Error in VlasovInitialize:\n");
115  fprintf(stderr, " space + vel dims not equal to ndims!\n");
116  }
117  return(1);
118  }
119 
120  if (!strcmp(solver->SplitHyperbolicFlux,"yes")) {
121  if (!mpi->rank) {
122  fprintf(stderr,"Error in VlasovInitialize: This physical model does not have a splitting ");
123  fprintf(stderr,"of the hyperbolic term defined.\n");
124  }
125  return(1);
126  }
127 
128 #ifndef serial
129  /* Broadcast parsed problem data */
130  MPIBroadcast_integer(&physics->ndims_x,1,0,&mpi->world);
131  MPIBroadcast_integer(&physics->ndims_v,1,0,&mpi->world);
133  1,0,&mpi->world);
134 #endif
135 
136  /* compute local number of x-space points with ghosts */
137  physics->npts_local_x_wghosts = 1;
138  physics->npts_local_x = 1;
139  physics->npts_global_x_wghosts = 1;
140  physics->npts_global_x = 1;
141  for (int d=0; d<physics->ndims_x; d++) {
142  physics->npts_local_x_wghosts *= (dim_local[d]+2*ghosts);
143  physics->npts_local_x *= dim_local[d];
144  physics->npts_global_x_wghosts *= (dim_global[d]+2*ghosts);
145  physics->npts_global_x *= dim_global[d];
146  }
147 
148  /* allocate array to hold the electric field (needs to have ghosts);
149  note that number of electric field components is the number of
150  spatial dimensions */
151  physics->e_field = (double*) calloc( physics->npts_local_x_wghosts
152  * physics->ndims_x,
153  sizeof(double) );
154  physics->potential = (double*) calloc( physics->npts_local_x_wghosts
155  * physics->ndims_x,
156  sizeof(double) );
157 
158  /* Put the mpi object in the params for access in other functions */
159  physics->m = m;
160 
161  if (physics->self_consistent_electric_field) {
162 #ifdef fftw
163  /* If using FFTW, make sure MPI is enabled
164  Currently we only support the distrubuted memory version
165  */
166 #ifdef serial
167  fprintf(stderr,"Error in VlasovInitialize(): Using FFTW requires MPI to be enabled.\n");
168  return(1);
169 #else
170 
171  if (physics->ndims_x > 1) {
172  fprintf(stderr,"Error in VlasovInitialize():\n");
173  fprintf(stderr," Self-consistent electric field is implemented for only 1 space dimension.\n");
174  return 1;
175  }
176 
177  /* Create a scratch buffer for moving between real and complex values */
178  physics->sum_buffer = (double*) calloc(dim_local[0], sizeof(double));
179 
180  /* Initialize FFTW and set up data buffers used for the transforms */
181  fftw_mpi_init();
182  physics->alloc_local = fftw_mpi_local_size_1d(dim_global[0], mpi->comm[0],
183  FFTW_FORWARD, 0,
184  &physics->local_ni,
185  &physics->local_i_start,
186  &physics->local_no,
187  &physics->local_o_start);
188  if (dim_local[0] != physics->local_ni) {
189  fprintf(stderr,"Error in VlasovInitialize(): The FFTW data distribution is incompatible with the HyPar one.\n");
190  fprintf(stderr,"Decompose the spatial dimension so that the degrees of freedom are evenly divided.\n");
191  return(1);
192  }
193 
194  physics->phys_buffer_e = fftw_alloc_complex(physics->alloc_local);
195  physics->fourier_buffer_e = fftw_alloc_complex(physics->alloc_local);
196 
197  physics->plan_forward_e = fftw_mpi_plan_dft_1d(dim_global[0],
198  physics->phys_buffer_e,
199  physics->fourier_buffer_e,
200  mpi->comm[0],
201  FFTW_FORWARD,
202  FFTW_ESTIMATE);
203 
204  physics->plan_backward_e = fftw_mpi_plan_dft_1d(dim_global[0],
205  physics->fourier_buffer_e,
206  physics->phys_buffer_e,
207  mpi->comm[0],
208  FFTW_BACKWARD,
209  FFTW_ESTIMATE);
210 
211  physics->phys_buffer_phi = fftw_alloc_complex(physics->alloc_local);
212  physics->fourier_buffer_phi = fftw_alloc_complex(physics->alloc_local);
213 
214  physics->plan_forward_phi = fftw_mpi_plan_dft_1d(dim_global[0],
215  physics->phys_buffer_phi,
216  physics->fourier_buffer_phi,
217  mpi->comm[0],
218  FFTW_FORWARD,
219  FFTW_ESTIMATE);
220 
221  physics->plan_backward_phi = fftw_mpi_plan_dft_1d(dim_global[0],
222  physics->fourier_buffer_phi,
223  physics->phys_buffer_phi,
224  mpi->comm[0],
225  FFTW_BACKWARD,
226  FFTW_ESTIMATE);
227 #endif
228 #else
229  fprintf(stderr,"Error in VlasovInitialize():\n");
230  fprintf(stderr," Self-consistent electric field requires FFTW library.\n");
231  return(1);
232 #endif
233  }
234 
235  /* initializing physical model-specific functions */
236  solver->PreStep = VlasovPreStep;
237  solver->ComputeCFL = VlasovComputeCFL;
238  solver->FFunction = VlasovAdvection;
239  solver->Upwind = VlasovUpwind;
241  solver->PostStage = VlasovPostStage;
242 
243  int ierr = VlasovEField(solver->u, solver, 0.0);
244  if (ierr) return ierr;
245 
246  return 0;
247 }
Definition: vlasov.h:57
int(* PhysicsOutput)(void *, void *, double)
Definition: hypar.h:347
int nvars
Definition: hypar.h:29
double * potential
Definition: vlasov.h:84
char SplitHyperbolicFlux[_MAX_STRING_SIZE_]
Definition: hypar.h:92
int self_consistent_electric_field
Definition: vlasov.h:60
int MPIBroadcast_integer(int *, int, int, void *)
Definition: MPIBroadcast.c:23
double * u
Definition: hypar.h:116
int VlasovUpwind(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: VlasovUpwind.c:17
double * e_field
Definition: vlasov.h:81
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int ndims
Definition: hypar.h:26
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
int npts_local_x
Definition: vlasov.h:69
double(* ComputeCFL)(void *, void *, double, double)
Definition: hypar.h:269
long npts_global_x
Definition: vlasov.h:72
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _MAX_STRING_SIZE_
Definition: basic.h:14
int VlasovWriteEFieldAndPotential(void *, void *, double)
int VlasovPreStep(double *, void *, void *, double)
long npts_global_x_wghosts
Definition: vlasov.h:78
#define _MODEL_NDIMS_
Definition: euler1d.h:56
MPI_Comm world
int VlasovEField(double *, void *, double)
Definition: VlasovEField.c:258
MPI_Comm * comm
int(* PreStep)(double *, void *, void *, double)
Definition: hypar.h:339
int * dim_local
Definition: hypar.h:37
void * m
Definition: vlasov.h:87
int ndims_v
Definition: vlasov.h:66
int npts_local_x_wghosts
Definition: vlasov.h:75
void * physics
Definition: hypar.h:266
int(* PostStage)(double *, void *, void *, double)
Definition: hypar.h:336
int VlasovPostStage(double *, void *, void *, double)
int ghosts
Definition: hypar.h:52
int VlasovAdvection(double *, double *, int, void *, double)
Structure of MPI-related variables.
double VlasovComputeCFL(void *, void *, double, double)
#define _MODEL_NVARS_
Definition: euler1d.h:58
int * dim_global
Definition: hypar.h:33
int ndims_x
Definition: vlasov.h:63