1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
|
/* Copyright 2020 Remi Lehe
*
* This file is part of WarpX.
*
* License: BSD-3-Clause-LBNL
*/
#ifndef WARPX_FINITE_DIFFERENCE_SOLVER_H_
#define WARPX_FINITE_DIFFERENCE_SOLVER_H_
#include <AMReX_MultiFab.H>
/**
* \brief Top-level class for the electromagnetic finite-difference solver
*
* Stores the coefficients of the finite-difference stencils,
* and has member functions to update fields over one time step.
*/
class FiniteDifferenceSolver
{
public:
// Constructor
/** \brief Initialize the finite-difference Maxwell solver (for a given refinement level)
*
* This function initializes the stencil coefficients for the chosen finite-difference algorithm
*
* \param fdtd_algo Identifies the chosen algorithm, as defined in WarpXAlgorithmSelection.H
* \param cell_size Cell size along each dimension, for the chosen refinement level
* \param do_nodal Whether the solver is applied to a nodal or staggered grid
*/
FiniteDifferenceSolver (
int const fdtd_algo,
std::array<amrex::Real,3> cell_size,
bool const do_nodal );
void EvolveB ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
amrex::Real const dt );
void EvolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
std::unique_ptr<amrex::MultiFab> const& Ffield,
amrex::Real const dt );
void EvolveF ( std::unique_ptr<amrex::MultiFab>& Ffield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
std::unique_ptr<amrex::MultiFab> const& rhofield,
int const rhocomp,
amrex::Real const dt );
void ComputeDivE ( const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
amrex::MultiFab& divE );
private:
int m_fdtd_algo;
bool m_do_nodal;
#ifdef WARPX_DIM_RZ
amrex::Real m_dr, m_rmin;
amrex::Real m_nmodes;
amrex::Gpu::ManagedVector<amrex::Real> m_stencil_coefs_r;
amrex::Gpu::ManagedVector<amrex::Real> m_stencil_coefs_z;
#else
amrex::Gpu::ManagedVector<amrex::Real> m_stencil_coefs_x;
amrex::Gpu::ManagedVector<amrex::Real> m_stencil_coefs_y;
amrex::Gpu::ManagedVector<amrex::Real> m_stencil_coefs_z;
#endif
public:
// The member functions below contain extended __device__ lambda.
// In order to compile with nvcc, they need to be public.
#ifdef WARPX_DIM_RZ
template< typename T_Algo >
void EvolveBCylindrical (
std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
amrex::Real const dt );
template< typename T_Algo >
void EvolveECylindrical (
std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
std::unique_ptr<amrex::MultiFab> const& Ffield,
amrex::Real const dt );
template< typename T_Algo >
void EvolveFCylindrical (
std::unique_ptr<amrex::MultiFab>& Ffield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
std::unique_ptr<amrex::MultiFab> const& rhofield,
int const rhocomp,
amrex::Real const dt );
template< typename T_Algo >
void ComputeDivECylindrical (
const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
amrex::MultiFab& divE );
#else
template< typename T_Algo >
void EvolveBCartesian (
std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
amrex::Real const dt );
template< typename T_Algo >
void EvolveECartesian (
std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
std::unique_ptr<amrex::MultiFab> const& Ffield,
amrex::Real const dt );
template< typename T_Algo >
void EvolveFCartesian (
std::unique_ptr<amrex::MultiFab>& Ffield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
std::unique_ptr<amrex::MultiFab> const& rhofield,
int const rhocomp,
amrex::Real const dt );
template< typename T_Algo >
void ComputeDivECartesian (
const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
amrex::MultiFab& divE );
#endif
};
#endif // WARPX_FINITE_DIFFERENCE_SOLVER_H_
|