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

Contains the function to compute the electric field. More...

#include <stdlib.h>
#include <stdbool.h>
#include <math.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <physicalmodels/vlasov.h>
#include <mpivars.h>
#include <hypar.h>

Go to the source code of this file.

Functions

static int FFTFreqNum (int bin, int N)
 
static int SetEFieldPrescribed (double *u, void *s, double t)
 
static int SetEFieldSelfConsistent (double *u, void *s, double t)
 
int VlasovEField (double *u, void *s, double t)
 

Detailed Description

Contains the function to compute the electric field.

Author
John Loffeld, Ping-Hsuan Tsai

Definition in file VlasovEField.c.

Function Documentation

◆ FFTFreqNum()

static int FFTFreqNum ( int  bin,
int  N 
)
static

Compute the Fourier transform sample frequency for bin. Note that bin is 0-based. Calling this routine one bin at a time is not efficient, but we do not care about speed here.

Parameters
binThe bin number for the Fourier frequency to be returned
NThe total number of Fourier frequencies

Definition at line 20 of file VlasovEField.c.

23 {
24  double remainder;
25  int last_positive = 0;
26 
27  remainder = N % 2;
28  // Note that last_positive is an integer
29  if (remainder) {
30  // N is odd
31  last_positive = (N-1) / 2;
32  } else {
33  // N is even
34  last_positive = N/2 - 1;
35  }
36 
37  if (bin <= last_positive) {
38  return bin;
39  } else {
40  return -(N - bin);
41  }
42 }

◆ SetEFieldPrescribed()

static int SetEFieldPrescribed ( double *  u,
void *  s,
double  t 
)
static

Sets the prescribed electric field.

Parameters
uConserved solution
sSolver object of type HyPar
tCurrent time

Definition at line 45 of file VlasovEField.c.

49 {
50 
51  HyPar *solver = (HyPar*) s;
52  Vlasov *param = (Vlasov*) solver->physics;
53  MPIVariables *mpi = (MPIVariables *) param->m;
54 
55  int* dim_local = solver->dim_local;
56  int ghosts = solver->ghosts;
57 
58  int ndims_x = param->ndims_x;
59 
60  int dim_x[ndims_x];
61  _ArrayCopy1D_(dim_local,dim_x,ndims_x);
62 
63  int bounds[ndims_x];
64  _ArrayCopy1D_(dim_x,bounds,ndims_x);
65  for (int d=0; d<ndims_x; d++) bounds[d] += 2*ghosts;
66 
67  int offset[ndims_x];
68  _ArraySetValue_(offset,ndims_x,-ghosts);
69 
70  int done = 0;
71  int index[ndims_x]; _ArraySetValue_(index, ndims_x, 0);
72 
73  while (!done) {
74 
75  double xvec[ndims_x];
76  for (int d=0; d<ndims_x; d++) {
77  _GetCoordinate_(d,index[d]-ghosts,dim_local,ghosts,solver->x,xvec[d]);
78  }
79 
80  int p;
81  _ArrayIndex1DWO_(ndims_x,dim_x,index,offset,ghosts,p);
82  param->e_field[p*ndims_x+0] = 0.1 * cos(xvec[0]);
83 
84  _ArrayIncrementIndex_(ndims_x,bounds,index,done);
85  }
86 
87  return 0;
88 }
Definition: vlasov.h:57
double * x
Definition: hypar.h:107
double * e_field
Definition: vlasov.h:81
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
void * m
Definition: vlasov.h:87
#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
Structure of MPI-related variables.
#define _ArrayCopy1D_(x, y, size)
int ndims_x
Definition: vlasov.h:63

◆ SetEFieldSelfConsistent()

static int SetEFieldSelfConsistent ( double *  u,
void *  s,
double  t 
)
static

Compute the self-consistent electric field over the local domain: The field is solved from the solution values using a Poisson solve in Fourier space.

Parameters
uConserved solution
sSolver object of type HyPar
tCurrent time

Definition at line 92 of file VlasovEField.c.

96 {
97  HyPar *solver = (HyPar*) s;
98  Vlasov *param = (Vlasov*) solver->physics;
99  MPIVariables *mpi = (MPIVariables *) param->m;
100 
101  if (param->ndims_x > 1) {
102  fprintf(stderr,"Error in SetEFieldSelfConsistent():\n");
103  fprintf(stderr," Implemented for 1 spatial dimension only.\n");
104  return 1;
105  }
106 
107 #ifndef fftw
108 
109  fprintf(stderr,"Error in SetEFieldSelfConsistent():\n");
110  fprintf(stderr," Using a self-consistent electric field requires FFTW.\n");
111  exit(1);
112 
113 #else
114 
115  int *dim = solver->dim_local;
116  int N = solver->dim_global[0];
117  int ghosts = solver->ghosts;
118  int ndims = solver->ndims;
119 
120  double *sum_buffer = param->sum_buffer;
121  double *field = param->e_field;
122  fftw_complex *phys_buffer_e = param->phys_buffer_e;
123  fftw_complex *fourier_buffer_e = param->fourier_buffer_e;
124  fftw_plan plan_forward_e = param->plan_forward_e;
125  fftw_plan plan_backward_e = param->plan_backward_e;
126  ptrdiff_t local_ni = param->local_ni;
127  ptrdiff_t local_i_start = param->local_i_start;
128  ptrdiff_t local_no = param->local_no;
129  ptrdiff_t local_o_start = param->local_o_start;
130 
131  fftw_complex *phys_buffer_phi = param->phys_buffer_phi;
132  fftw_complex *fourier_buffer_phi = param->fourier_buffer_phi;
133  fftw_plan plan_forward_phi = param->plan_forward_phi;
134  fftw_plan plan_backward_phi = param->plan_backward_phi;
135 
136  int index[ndims], bounds[ndims], bounds_noghost[ndims], offset[ndims];
137 
138  // set bounds for array index to include ghost points
139  _ArrayCopy1D_(dim,bounds,ndims);
140  for (int i = 0; i < ndims; i++) bounds[i] += 2*ghosts;
141 
142  // set bounds for array index to NOT include ghost points
143  _ArrayCopy1D_(dim,bounds_noghost,ndims);
144 
145  // set offset such that index is compatible with ghost point arrangement
146  _ArraySetValue_(offset,ndims,-ghosts);
147 
148  // First, integrate the particle distribution over velocity.
149  // Since the array dimension we want is not unit stride,
150  // first manually add up the local part of the array.
151  int done = 0; _ArraySetValue_(index,ndims,0);
152  _ArraySetValue_(sum_buffer,dim[0],0);
153  while (!done) {
154  //int p; _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,p);
155  int p; _ArrayIndex1D_(ndims,dim,index,ghosts,p);
156 
157  // accumulate f at this spatial location
158  double dvinv; _GetCoordinate_(1,index[1],dim,ghosts,solver->dxinv,dvinv);
159  double x; _GetCoordinate_(0,index[0],dim,ghosts,solver->x,x);
160  double v; _GetCoordinate_(1,index[1],dim,ghosts,solver->x,v);
161 
162  sum_buffer[index[0]] += u[p] / dvinv;
163 
164  _ArrayIncrementIndex_(ndims,bounds_noghost,index,done);
165  }
166 
167  // Now we can add up globally using MPI reduction
168  for (int i = 0; i < dim[0]; i++) {
169  MPISum_double(&sum_buffer[i], &sum_buffer[i], 1, &mpi->comm[1]);
170  }
171 
172  // Find the average density over all x
173  double average_velocity = 0.0;
174  for (int i = 0; i < dim[0]; i++) {
175  average_velocity += sum_buffer[i];
176  }
177  MPISum_double(&average_velocity, &average_velocity, 1, &mpi->comm[0]);
178  average_velocity /= (double) N;
179 
180  // Copy velocity-integrated values into complex-valued FFTW buffer
181  for (int i = 0; i < dim[0]; i++) {
182  phys_buffer_e[i][0] = sum_buffer[i] - average_velocity;
183  phys_buffer_e[i][1] = 0.0;
184  }
185 
186  // Execute the FFT
187  MPI_Barrier(mpi->comm[0]);
188  fftw_execute(plan_forward_e);
189  MPI_Barrier(mpi->comm[0]);
190 
191  // Simultaneously do a Poisson solve and take derivative in frequency space
192  int freq_start = 0;
193  if (local_o_start == 0) {
194  freq_start = 1;
195  fourier_buffer_e[0][0] = 0.0;
196  fourier_buffer_e[0][1] = 0.0;
197  }
198 
199  // Simultaneously do a Poisson solve in frequency space
200  int freq_start_phi = 0;
201  if (local_o_start == 0) {
202  freq_start_phi = 1;
203  fourier_buffer_phi[0][0] = 0.0;
204  fourier_buffer_phi[0][1] = 0.0;
205  }
206 
207  for (int i = freq_start_phi; i < local_no; i++) {
208  double freq_num = FFTFreqNum(i + local_o_start, N);
209  double thek = freq_num;
210  // Swapping values is due to multiplication by i
211  fourier_buffer_phi[i][0] = fourier_buffer_e[i][0] / (thek*thek);
212  fourier_buffer_phi[i][1] = fourier_buffer_e[i][1]/ (thek*thek);
213  }
214 
215 
216  for (int i = freq_start; i < local_no; i++) {
217  double freq_num = FFTFreqNum(i + local_o_start, N);
218  double thek = freq_num;
219  double temp = fourier_buffer_e[i][0];
220  // Swapping values is due to multiplication by i
221  fourier_buffer_e[i][0] = - fourier_buffer_e[i][1] / thek;
222  fourier_buffer_e[i][1] = temp / thek;
223  }
224 
225  // Do an inverse Fourier transform to get back physical solved values
226  MPI_Barrier(mpi->comm[0]);
227  fftw_execute(plan_backward_e);
228  MPI_Barrier(mpi->comm[0]);
229 
230  // Do an inverse Fourier transform to get back physical solved values
231  MPI_Barrier(mpi->comm[0]);
232  fftw_execute(plan_backward_phi);
233  MPI_Barrier(mpi->comm[0]);
234 
235  // copy the solved electric field into the e buffer
236  for (int i = 0; i < dim[0]; i++) {
237  param->e_field[i + ghosts] = - phys_buffer_e[i][0] / (double) N;
238  }
239 
240  // copy the solved potential field into the potential buffer
241  for (int i = 0; i < dim[0]; i++) {
242  param->potential[i + ghosts] = phys_buffer_phi[i][0] / (double) N;
243  }
244 
245  // Do halo exchange on the e
246  MPIExchangeBoundaries1D(mpi, param->e_field, dim[0], ghosts, 0, ndims);
247 
248  // Do halo exchange on the potential
249  MPIExchangeBoundaries1D(mpi, param->potential, dim[0], ghosts, 0, ndims);
250 
251 #endif
252 
253  return 0;
254 }
Definition: vlasov.h:57
double * potential
Definition: vlasov.h:84
static int FFTFreqNum(int bin, int N)
Definition: VlasovEField.c:20
int MPISum_double(double *, double *, int, void *)
Definition: MPISum.c:39
double * x
Definition: hypar.h:107
double * e_field
Definition: vlasov.h:81
int MPIExchangeBoundaries1D(void *, double *, int, int, int, int)
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
void * m
Definition: vlasov.h:87
#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
Structure of MPI-related variables.
#define _ArrayCopy1D_(x, y, size)
int * dim_global
Definition: hypar.h:33
int ndims_x
Definition: vlasov.h:63
double * dxinv
Definition: hypar.h:110

◆ 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