HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
vlasov.h File Reference

Vlasov Equation. More...

#include <stdbool.h>

Go to the source code of this file.

Data Structures

struct  Vlasov
 

Macros

#define _VLASOV_   "vlasov"
 
#define _MODEL_NDIMS_   2
 
#define _MODEL_NVARS_   1
 

Functions

int VlasovInitialize (void *, void *)
 
int VlasovCleanup (void *)
 

Detailed Description

Vlasov Equation.

Author
John Loffeld

1D-1V Vlasov Equation:

\begin{equation} \frac{\partial f}{\partial t} + v \frac{\partial f}{\partial x} + E \frac{\partial f}{\partial v} = 0, \end{equation}

where

  • \(f\) is the distribution function
  • \(x\) is the spatial coordinate
  • \(v\) is the velocity
  • \(E\) is the electric field

Reference:

  • Henon, "Vlasov equation?", Astronomy and Astrophysics, 114, 1982

Definition in file vlasov.h.


Data Structure Documentation

struct Vlasov

Definition at line 57 of file vlasov.h.

Data Fields
int self_consistent_electric_field

Use a self-consistent electric field? Requires FFTW library

int ndims_x

Number of spatial dimensions

int ndims_v

Number of velocity dimensions

int npts_local_x

Number of spatial grid points (local)

long npts_global_x

Number of spatial grid points (global)

int npts_local_x_wghosts

Number of spatial grid points with ghosts (local)

long npts_global_x_wghosts

Number of spatial grid points with ghosts (global)

double * e_field

electric field

double * potential

potential field

void * m

Pointer to MPI object of type MPIVariables

Macro Definition Documentation

#define _VLASOV_   "vlasov"

1D-1V Vlasov equation

Definition at line 37 of file vlasov.h.

#define _MODEL_NDIMS_   2

Number of spatial dimensions

Definition at line 43 of file vlasov.h.

#define _MODEL_NVARS_   1

Number of variables per grid point

Definition at line 45 of file vlasov.h.

Function Documentation

int VlasovInitialize ( void *  s,
void *  m 
)

Initialize the Vlasov physics module

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 }
MPI_Comm * comm
int ndims_x
Definition: vlasov.h:63
double(* ComputeCFL)(void *, void *, double, double)
Definition: hypar.h:269
double * e_field
Definition: vlasov.h:81
int * dim_global
Definition: hypar.h:33
double * u
Definition: hypar.h:116
double VlasovComputeCFL(void *s, void *m, double dt, double t)
#define _MODEL_NVARS_
Definition: euler1d.h:58
void * physics
Definition: hypar.h:266
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int VlasovUpwind(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: VlasovUpwind.c:17
int * dim_local
Definition: hypar.h:37
int(* PhysicsOutput)(void *, void *, double)
Definition: hypar.h:347
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
#define _MODEL_NDIMS_
Definition: euler1d.h:56
double * potential
Definition: vlasov.h:84
int ghosts
Definition: hypar.h:52
int VlasovPostStage(double *u, void *s, void *m, double waqt)
int VlasovEField(double *u, void *s, double t)
Definition: VlasovEField.c:258
MPI_Comm world
#define _MAX_STRING_SIZE_
Definition: basic.h:14
int(* PostStage)(double *, void *, void *, double)
Definition: hypar.h:336
long npts_global_x
Definition: vlasov.h:72
int npts_local_x_wghosts
Definition: vlasov.h:75
void * m
Definition: vlasov.h:87
int VlasovPreStep(double *u, void *s, void *m, double waqt)
Definition: vlasov.h:57
int npts_local_x
Definition: vlasov.h:69
int nvars
Definition: hypar.h:29
int VlasovAdvection(double *f, double *u, int dir, void *s, double t)
char SplitHyperbolicFlux[_MAX_STRING_SIZE_]
Definition: hypar.h:92
int ndims
Definition: hypar.h:26
int self_consistent_electric_field
Definition: vlasov.h:60
int(* PreStep)(double *, void *, void *, double)
Definition: hypar.h:339
Structure of MPI-related variables.
long npts_global_x_wghosts
Definition: vlasov.h:78
int ndims_v
Definition: vlasov.h:66
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
int MPIBroadcast_integer(int *, int, int, void *)
Definition: MPIBroadcast.c:23
int VlasovWriteEFieldAndPotential(void *, void *, double)
int VlasovCleanup ( void *  s)

Clean up the Vlasov physics module

Function to clean up all physics-related allocations for the Vlasov equations

Parameters
sSolver object of type HyPar

Definition at line 10 of file VlasovCleanup.c.

11 {
12  Vlasov *physics = (Vlasov*) s;
13 
14  free(physics->e_field);
15  free(physics->potential);
16 
17 #ifdef fftw
18  if(physics->self_consistent_electric_field) {
19  free(physics->sum_buffer);
20 
21  fftw_destroy_plan(physics->plan_forward_e);
22  fftw_destroy_plan(physics->plan_backward_e);
23 
24  fftw_free(physics->phys_buffer_e);
25  fftw_free(physics->fourier_buffer_e);
26 
27  fftw_destroy_plan(physics->plan_forward_phi);
28  fftw_destroy_plan(physics->plan_backward_phi);
29 
30  fftw_free(physics->phys_buffer_phi);
31  fftw_free(physics->fourier_buffer_phi);
32  }
33 #endif
34 
35  return(0);
36 }
double * e_field
Definition: vlasov.h:81
double * potential
Definition: vlasov.h:84
Definition: vlasov.h:57
int self_consistent_electric_field
Definition: vlasov.h:60