WarpX
Loading...
Searching...
No Matches
FiniteDifferenceSolver.H
Go to the documentation of this file.
1/* Copyright 2020-2024 The WarpX Community
2 *
3 * This file is part of WarpX.
4 *
5 * Authors: Remi Lehe (LBNL)
6 * S. Eric Clark (Helion Energy)
7 *
8 * License: BSD-3-Clause-LBNL
9 */
10
11#ifndef WARPX_FINITE_DIFFERENCE_SOLVER_H_
12#define WARPX_FINITE_DIFFERENCE_SOLVER_H_
13
17
21
22#include <ablastr/utils/Enums.H>
24
25#include <AMReX_GpuContainers.H>
26#include <AMReX_REAL.H>
27
28#include <AMReX_BaseFwd.H>
29
30#include <array>
31#include <memory>
32
40{
41 public:
42
43 // Constructor
54 std::array<amrex::Real,3> cell_size,
56
58 int lev,
59 PatchType patch_type,
60 std::array< std::unique_ptr<amrex::iMultiFab>, 3 >& flag_info_cell,
61 std::array< std::unique_ptr<amrex::LayoutData<FaceInfoBox> >, 3 >& borrowing,
62 amrex::Real dt );
63
65 int lev,
66 PatchType patch_type,
67 ablastr::fields::VectorField const& Efield,
68 std::array< std::unique_ptr<amrex::iMultiFab>,3 > const& eb_update_E,
69 amrex::Real dt );
70
71 void EvolveF ( amrex::MultiFab* Ffield,
72 ablastr::fields::VectorField const& Efield,
73 amrex::MultiFab* rhofield,
74 int rho_comp,
75 amrex::Real dt );
76
77 void EvolveG (amrex::MultiFab* Gfield,
78 ablastr::fields::VectorField const& Bfield,
79 amrex::Real dt);
80
81 void EvolveECTRho ( ablastr::fields::VectorField const& Efield,
82 ablastr::fields::VectorField const& edge_lengths,
83 ablastr::fields::VectorField const& face_areas,
84 ablastr::fields::VectorField const& ECTRhofield,
85 int lev );
86
90 amrex::Box domain_box,
91 amrex::Real dt,
94
95 void ComputeDivE (
96 ablastr::fields::VectorField const & Efield,
97 amrex::MultiFab& divE
98 );
99
113 void MacroscopicEvolveE (
114 MacroscopicSolverAlgo macroscopic_solver_algo,
115 ablastr::fields::VectorField const& Efield,
116 ablastr::fields::VectorField const& Bfield,
117 ablastr::fields::VectorField const& Jfield,
118 std::array< std::unique_ptr<amrex::iMultiFab>,3 > const& eb_update_E,
119 amrex::Real dt,
120 std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
121
122 void EvolveBPML (
124 PatchType patch_type,
125 int level,
126 amrex::Real dt,
127 bool dive_cleaning
128 );
129
130 void EvolveEPML (
132 PatchType patch_type,
133 int level,
134 MultiSigmaBox const& sigba,
135 amrex::Real dt,
136 bool pml_has_particles
137 );
138
139 void EvolveFPML ( amrex::MultiFab* Ffield,
141 amrex::Real dt );
142
161 ablastr::fields::VectorField const& Jifield,
162 ablastr::fields::VectorField const& Bfield,
163 amrex::MultiFab const& rhofield,
164 amrex::MultiFab const& Pefield,
165 std::array< std::unique_ptr<amrex::iMultiFab>,3> const& eb_update_E,
166 int lev, HybridPICModel const* hybrid_model,
167 bool solve_for_Faraday );
168
180 ablastr::fields::VectorField const& Bfield,
181 std::array< std::unique_ptr<amrex::iMultiFab>,3> const& eb_update_E,
182 int lev );
183
193 void ComputeCurlA (
195 ablastr::fields::VectorField const& Afield,
196 std::array< std::unique_ptr<amrex::iMultiFab>,3> const& eb_update_B,
197 int lev );
198
199 private:
200
203
204#if defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER)
205 amrex::Real m_dr, m_rmin;
207 // host-only
209 // device copy after init
212#elif defined(WARPX_DIM_RSPHERE)
213 amrex::Real m_dr, m_rmin;
214 // host-only
216 // device copy after init
218#else
219 // host-only
220 amrex::Vector<amrex::Real> m_h_stencil_coefs_x, m_h_stencil_coefs_y, m_h_stencil_coefs_z;
221 // device copy after init
225#endif
226
227 public:
228 // The member functions below contain extended __device__ lambda.
229 // In order to compile with nvcc, they need to be public.
230
231#if defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER)
232 template< typename T_Algo >
233 void EvolveBCylindrical (
234 ablastr::fields::VectorField const& Bfield,
235 ablastr::fields::VectorField const& Efield,
236 int lev,
237 amrex::Real dt );
238
239 template< typename T_Algo >
240 void EvolveECylindrical (
241 ablastr::fields::VectorField const& Efield,
242 ablastr::fields::VectorField const& Bfield,
243 ablastr::fields::VectorField const& Jfield,
244 std::array< std::unique_ptr<amrex::iMultiFab>,3 > const& eb_update_E,
245 amrex::MultiFab const* Ffield,
246 int lev,
247 amrex::Real dt );
248
249 template< typename T_Algo >
250 void EvolveFCylindrical (
251 amrex::MultiFab* Ffield,
252 ablastr::fields::VectorField const & Efield,
253 amrex::MultiFab* rhofield,
254 int rho_comp,
255 amrex::Real dt );
256
257 template< typename T_Algo >
259 ablastr::fields::VectorField const & Efield,
260 amrex::MultiFab& divE
261 );
262
263 template<typename T_Algo>
265 ablastr::fields::VectorField const& Efield,
266 ablastr::fields::VectorField const& Jfield,
267 ablastr::fields::VectorField const& Jifield,
268 ablastr::fields::VectorField const& Bfield,
269 amrex::MultiFab const& rhofield,
270 amrex::MultiFab const& Pefield,
271 std::array< std::unique_ptr<amrex::iMultiFab>,3> const& eb_update_E,
272 int lev, HybridPICModel const* hybrid_model,
273 bool solve_for_Faraday );
274
275 template<typename T_Algo>
278 ablastr::fields::VectorField const& Bfield,
279 std::array< std::unique_ptr<amrex::iMultiFab>,3> const& eb_update_E,
280 int lev
281 );
282
283 template<typename T_Algo>
286 ablastr::fields::VectorField const& Afield,
287 std::array< std::unique_ptr<amrex::iMultiFab>,3> const& eb_update_B,
288 int lev
289 );
290
291#elif defined(WARPX_DIM_RSPHERE)
292 template< typename T_Algo >
293 void EvolveBSpherical (
294 ablastr::fields::VectorField const& Bfield,
295 ablastr::fields::VectorField const& Efield,
296 int lev,
297 amrex::Real dt );
298
299 template< typename T_Algo >
300 void EvolveESpherical (
301 ablastr::fields::VectorField const& Efield,
302 ablastr::fields::VectorField const& Bfield,
303 ablastr::fields::VectorField const& Jfield,
304 amrex::MultiFab const* Ffield,
305 int lev,
306 amrex::Real dt );
307
308 template< typename T_Algo >
309 void EvolveFSpherical (
310 amrex::MultiFab* Ffield,
311 ablastr::fields::VectorField const & Efield,
312 amrex::MultiFab* rhofield,
313 int rho_comp,
314 amrex::Real dt );
315
316 template< typename T_Algo >
317 void ComputeDivESpherical (
318 ablastr::fields::VectorField const & Efield,
319 amrex::MultiFab& divE
320 );
321
322 template<typename T_Algo>
323 void HybridPICSolveESpherical (
324 ablastr::fields::VectorField const& Efield,
325 ablastr::fields::VectorField const& Jfield,
326 ablastr::fields::VectorField const& Jifield,
327 ablastr::fields::VectorField const& Bfield,
328 amrex::MultiFab const& rhofield,
329 amrex::MultiFab const& Pefield,
330 int lev, HybridPICModel const* hybrid_model,
331 bool solve_for_Faraday );
332
333 template<typename T_Algo>
334 void CalculateCurrentAmpereSpherical (
336 ablastr::fields::VectorField const& Bfield,
337 int lev
338 );
339
340 template<typename T_Algo>
341 void ComputeCurlASpherical (
343 ablastr::fields::VectorField const& Afield,
344 std::array< std::unique_ptr<amrex::iMultiFab>,3> const& eb_update_B,
345 int lev
346 );
347
348#else
349 template< typename T_Algo >
350 void EvolveBCartesian (
351 ablastr::fields::VectorField const& Bfield,
352 ablastr::fields::VectorField const& Efield,
353 amrex::MultiFab const * Gfield,
354 int lev, amrex::Real dt );
355
356 template< typename T_Algo >
357 void EvolveECartesian (
358 ablastr::fields::VectorField const& Efield,
359 ablastr::fields::VectorField const& Bfield,
360 ablastr::fields::VectorField const& Jfield,
361 std::array< std::unique_ptr<amrex::iMultiFab>,3 > const& eb_update_E,
362 amrex::MultiFab const* Ffield,
363 int lev, amrex::Real dt );
364
365 template< typename T_Algo >
366 void EvolveFCartesian (
367 amrex::MultiFab* Ffield,
369 amrex::MultiFab* rhofield,
370 int rho_comp,
371 amrex::Real dt );
372
373 template< typename T_Algo >
374 void EvolveGCartesian (
375 amrex::MultiFab* Gfield,
376 ablastr::fields::VectorField const& Bfield,
377 amrex::Real dt);
378
379 void EvolveRhoCartesianECT (
380 ablastr::fields::VectorField const& Efield,
381 ablastr::fields::VectorField const& edge_lengths,
382 ablastr::fields::VectorField const& face_areas,
383 ablastr::fields::VectorField const& ECTRhofield, int lev);
384
385 void EvolveBCartesianECT (
386 ablastr::fields::VectorField const& Bfield,
387 ablastr::fields::VectorField const& face_areas,
388 ablastr::fields::VectorField const& area_mod,
389 ablastr::fields::VectorField const& ECTRhofield,
391 std::array< std::unique_ptr<amrex::iMultiFab>, 3 >& flag_info_cell,
392 std::array< std::unique_ptr<amrex::LayoutData<FaceInfoBox> >, 3 >& borrowing,
393 int lev, amrex::Real dt
394 );
395
396 template< typename T_Algo >
397 void ComputeDivECartesian (
398 ablastr::fields::VectorField const & Efield,
399 amrex::MultiFab& divE );
400
401 template< typename T_Algo, typename T_MacroAlgo >
402 void MacroscopicEvolveECartesian (
403 ablastr::fields::VectorField const& Efield,
404 ablastr::fields::VectorField const& Bfield,
405 ablastr::fields::VectorField const& Jfield,
406 std::array< std::unique_ptr<amrex::iMultiFab>,3 > const& eb_update_E,
407 amrex::Real dt,
408 std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
409
410 template< typename T_Algo >
411 void EvolveBPMLCartesian (
412 std::array< amrex::MultiFab*, 3 > Bfield,
414 amrex::Real dt,
415 bool dive_cleaning);
416
417 template< typename T_Algo >
418 void EvolveEPMLCartesian (
420 std::array< amrex::MultiFab*, 3 > Bfield,
421 std::array< amrex::MultiFab*, 3 > Jfield,
422 std::array< amrex::MultiFab*, 3 > edge_lengths,
423 amrex::MultiFab* Ffield,
424 MultiSigmaBox const& sigba,
425 amrex::Real dt, bool pml_has_particles );
426
427 template< typename T_Algo >
428 void EvolveFPMLCartesian ( amrex::MultiFab* Ffield,
430 amrex::Real dt );
431
432 template<typename T_Algo>
433 void HybridPICSolveECartesian (
434 ablastr::fields::VectorField const& Efield,
435 ablastr::fields::VectorField const& Jfield,
436 ablastr::fields::VectorField const& Jifield,
437 ablastr::fields::VectorField const& Bfield,
438 amrex::MultiFab const& rhofield,
439 amrex::MultiFab const& Pefield,
440 std::array< std::unique_ptr<amrex::iMultiFab>,3> const& eb_update_E,
441 int lev, HybridPICModel const* hybrid_model,
442 bool solve_for_Faraday );
443
444 template<typename T_Algo>
445 void CalculateCurrentAmpereCartesian (
447 ablastr::fields::VectorField const& Bfield,
448 std::array< std::unique_ptr<amrex::iMultiFab>,3> const& eb_update_E,
449 int lev
450 );
451
452 template<typename T_Algo>
453 void ComputeCurlACartesian (
455 ablastr::fields::VectorField const& Afield,
456 std::array< std::unique_ptr<amrex::iMultiFab>,3> const& eb_update_B,
457 int lev
458 );
459#endif
460
461};
462
463#endif // WARPX_FINITE_DIFFERENCE_SOLVER_H_
ElectromagneticSolverAlgo
Definition WarpXAlgorithmSelection.H:58
MacroscopicSolverAlgo
struct to select algorithm for macroscopic Maxwell solver LaxWendroff (semi-implicit) represents sigm...
Definition WarpXAlgorithmSelection.H:48
void ComputeDivE(ablastr::fields::VectorField const &Efield, amrex::MultiFab &divE)
Update the F field, over one timestep.
Definition ComputeDivE.cpp:46
int m_nmodes
Definition FiniteDifferenceSolver.H:206
void ComputeCurlA(ablastr::fields::VectorField &Bfield, ablastr::fields::VectorField const &Afield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_B, int lev)
Calculation of B field from the vector potential A B = (curl x A) / mu0.
Definition ComputeCurlA.cpp:27
void EvolveB(ablastr::fields::MultiFabRegister &fields, int lev, PatchType patch_type, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > &flag_info_cell, std::array< std::unique_ptr< amrex::LayoutData< FaceInfoBox > >, 3 > &borrowing, amrex::Real dt)
Update the B field, over one timestep.
Definition EvolveB.cpp:53
ablastr::utils::enums::GridType m_grid_type
Definition FiniteDifferenceSolver.H:202
FiniteDifferenceSolver(ElectromagneticSolverAlgo fdtd_algo, std::array< amrex::Real, 3 > cell_size, ablastr::utils::enums::GridType grid_type)
Initialize the finite-difference Maxwell solver (for a given refinement level)
Definition FiniteDifferenceSolver.cpp:32
void EvolveECTRho(ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField const &edge_lengths, ablastr::fields::VectorField const &face_areas, ablastr::fields::VectorField const &ECTRhofield, int lev)
Update the B field, over one timestep.
Definition EvolveECTRho.cpp:48
void EvolveE(ablastr::fields::MultiFabRegister &fields, int lev, PatchType patch_type, ablastr::fields::VectorField const &Efield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, amrex::Real dt)
Update the E field, over one timestep.
Definition EvolveE.cpp:55
void EvolveBCylindrical(ablastr::fields::VectorField const &Bfield, ablastr::fields::VectorField const &Efield, int lev, amrex::Real dt)
Definition EvolveB.cpp:394
void EvolveF(amrex::MultiFab *Ffield, ablastr::fields::VectorField const &Efield, amrex::MultiFab *rhofield, int rho_comp, amrex::Real dt)
Update the F field, over one timestep.
Definition EvolveF.cpp:48
void HybridPICSolveE(ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField &Jfield, ablastr::fields::VectorField const &Jifield, ablastr::fields::VectorField const &Bfield, amrex::MultiFab const &rhofield, amrex::MultiFab const &Pefield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, int lev, HybridPICModel const *hybrid_model, bool solve_for_Faraday)
E-update in the hybrid PIC algorithm as described in Winske et al. (2003) Eq. 10. https://link....
Definition HybridPICSolveE.cpp:484
void EvolveEPML(ablastr::fields::MultiFabRegister &fields, PatchType patch_type, int level, MultiSigmaBox const &sigba, amrex::Real dt, bool pml_has_particles)
Update the E field, over one timestep.
Definition EvolveEPML.cpp:46
void ComputeCurlACylindrical(ablastr::fields::VectorField &Bfield, ablastr::fields::VectorField const &Afield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_B, int lev)
Calculate B from the curl of A i.e. B = curl(A) output field on B field mesh staggering.
Definition ComputeCurlA.cpp:77
void EvolveFPML(amrex::MultiFab *Ffield, ablastr::fields::VectorField Efield, amrex::Real dt)
Update the E field, over one timestep.
Definition EvolveFPML.cpp:40
void EvolveFCylindrical(amrex::MultiFab *Ffield, ablastr::fields::VectorField const &Efield, amrex::MultiFab *rhofield, int rho_comp, amrex::Real dt)
Definition EvolveF.cpp:144
void MacroscopicEvolveE(MacroscopicSolverAlgo macroscopic_solver_algo, ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField const &Bfield, ablastr::fields::VectorField const &Jfield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, amrex::Real dt, std::unique_ptr< MacroscopicProperties > const &macroscopic_properties)
Macroscopic E-update for non-vacuum medium using the user-selected finite-difference algorithm and ma...
Definition MacroscopicEvolveE.cpp:38
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_z
Definition FiniteDifferenceSolver.H:211
void HybridPICSolveECylindrical(ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField const &Jfield, ablastr::fields::VectorField const &Jifield, ablastr::fields::VectorField const &Bfield, amrex::MultiFab const &rhofield, amrex::MultiFab const &Pefield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, int lev, HybridPICModel const *hybrid_model, bool solve_for_Faraday)
Definition HybridPICSolveE.cpp:534
void EvolveBPML(ablastr::fields::MultiFabRegister &fields, PatchType patch_type, int level, amrex::Real dt, bool dive_cleaning)
Update the B field, over one timestep.
Definition EvolveBPML.cpp:42
amrex::Vector< amrex::Real > m_h_stencil_coefs_z
Definition FiniteDifferenceSolver.H:208
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_r
Definition FiniteDifferenceSolver.H:210
void EvolveG(amrex::MultiFab *Gfield, ablastr::fields::VectorField const &Bfield, amrex::Real dt)
Definition EvolveG.cpp:38
amrex::Real m_dr
Definition FiniteDifferenceSolver.H:205
void ComputeDivECylindrical(ablastr::fields::VectorField const &Efield, amrex::MultiFab &divE)
Definition ComputeDivE.cpp:136
void ApplySilverMuellerBoundary(ablastr::fields::VectorField &Efield, ablastr::fields::VectorField &Bfield, amrex::Box domain_box, amrex::Real dt, amrex::Array< FieldBoundaryType, 3 > field_boundary_lo, amrex::Array< FieldBoundaryType, 3 > field_boundary_hi)
Update the B field at the boundary, using the Silver-Mueller condition.
Definition ApplySilverMuellerBoundary.cpp:37
amrex::Vector< amrex::Real > m_h_stencil_coefs_r
Definition FiniteDifferenceSolver.H:208
amrex::Real m_rmin
Definition FiniteDifferenceSolver.H:205
ElectromagneticSolverAlgo m_fdtd_algo
Definition FiniteDifferenceSolver.H:201
void CalculateCurrentAmpereCylindrical(ablastr::fields::VectorField &Jfield, ablastr::fields::VectorField const &Bfield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, int lev)
Calculate total current from Ampere's law without displacement current i.e. J = 1/mu_0 curl x B.
Definition HybridPICSolveE.cpp:80
void EvolveECylindrical(ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField const &Bfield, ablastr::fields::VectorField const &Jfield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, amrex::MultiFab const *Ffield, int lev, amrex::Real dt)
Definition EvolveE.cpp:240
void CalculateCurrentAmpere(ablastr::fields::VectorField &Jfield, ablastr::fields::VectorField const &Bfield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, int lev)
Calculation of total current using Ampere's law (without displacement current): J = (curl x B) / mu0.
Definition HybridPICSolveE.cpp:31
This class contains the parameters needed to evaluate hybrid field solutions (kinetic ions with fluid...
Definition HybridPICModel.H:41
Definition PML.H:124
Definition EffectivePotentialPoissonSolver.H:63
std::array< amrex::MultiFab *, 3 > VectorField
Definition MultiFabRegister.H:191
GridType
Definition Enums.H:23
PatchType
Definition Enums.H:30
PODVector< T, ArenaAllocator< T > > DeviceVector
BoxND< 3 > Box
std::array< T, N > Array
Definition MultiFabRegister.H:262