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
|
/* Copyright 2020 Remi Lehe
*
* This file is part of WarpX.
*
* License: BSD-3-Clause-LBNL
*/
#include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H"
#include "BoundaryConditions/PMLComponent.H"
#ifndef WARPX_DIM_RZ
# include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianYeeAlgorithm.H"
# include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianCKCAlgorithm.H"
# include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianNodalAlgorithm.H"
#else
# include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CylindricalYeeAlgorithm.H"
#endif
#include "Utils/WarpXAlgorithmSelection.H"
#include <AMReX.H>
#include <AMReX_Array4.H>
#include <AMReX_Config.H>
#include <AMReX_Extension.H>
#include <AMReX_GpuContainers.H>
#include <AMReX_GpuControl.H>
#include <AMReX_GpuLaunch.H>
#include <AMReX_GpuQualifiers.H>
#include <AMReX_IndexType.H>
#include <AMReX_MFIter.H>
#include <AMReX_MultiFab.H>
#include <AMReX_REAL.H>
#include <AMReX_BaseFwd.H>
#include <array>
using namespace amrex;
/**
* \brief Update the E field, over one timestep
*/
void FiniteDifferenceSolver::EvolveFPML (
amrex::MultiFab* Ffield,
std::array< amrex::MultiFab*, 3 > const Efield,
amrex::Real const dt ) {
// Select algorithm (The choice of algorithm is a runtime option,
// but we compile code for each algorithm, using templates)
#ifdef WARPX_DIM_RZ
amrex::ignore_unused(Ffield, Efield, dt);
amrex::Abort("PML are not implemented in cylindrical geometry.");
#else
if (m_do_nodal) {
EvolveFPMLCartesian <CartesianNodalAlgorithm> ( Ffield, Efield, dt );
} else if (m_fdtd_algo == MaxwellSolverAlgo::Yee) {
EvolveFPMLCartesian <CartesianYeeAlgorithm> ( Ffield, Efield, dt );
} else if (m_fdtd_algo == MaxwellSolverAlgo::CKC) {
EvolveFPMLCartesian <CartesianCKCAlgorithm> ( Ffield, Efield, dt );
} else {
amrex::Abort("EvolveFPML: Unknown algorithm");
}
#endif
}
#ifndef WARPX_DIM_RZ
template<typename T_Algo>
void FiniteDifferenceSolver::EvolveFPMLCartesian (
amrex::MultiFab* Ffield,
std::array< amrex::MultiFab*, 3 > const Efield,
amrex::Real const dt ) {
// Loop through the grids, and over the tiles within each grid
#ifdef AMREX_USE_OMP
#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
#endif
for ( MFIter mfi(*Ffield, TilingIfNotGPU()); mfi.isValid(); ++mfi ) {
// Extract field data for this grid/tile
Array4<Real> const& F = Ffield->array(mfi);
Array4<Real> const& Ex = Efield[0]->array(mfi);
Array4<Real> const& Ey = Efield[1]->array(mfi);
Array4<Real> const& Ez = Efield[2]->array(mfi);
// Extract stencil coefficients
Real const * const AMREX_RESTRICT coefs_x = m_stencil_coefs_x.dataPtr();
int const n_coefs_x = m_stencil_coefs_x.size();
Real const * const AMREX_RESTRICT coefs_y = m_stencil_coefs_y.dataPtr();
int const n_coefs_y = m_stencil_coefs_y.size();
Real const * const AMREX_RESTRICT coefs_z = m_stencil_coefs_z.dataPtr();
int const n_coefs_z = m_stencil_coefs_z.size();
// Extract tileboxes for which to loop
Box const& tf = mfi.tilebox(Ffield->ixType().ixType());
// Loop over the cells and update the fields
amrex::ParallelFor(tf,
[=] AMREX_GPU_DEVICE (int i, int j, int k){
F(i, j, k, PMLComp::x) += dt * (
T_Algo::DownwardDx(Ex, coefs_x, n_coefs_x, i, j, k, PMLComp::xx)
+ T_Algo::DownwardDx(Ex, coefs_x, n_coefs_x, i, j, k, PMLComp::xy)
+ T_Algo::DownwardDx(Ex, coefs_x, n_coefs_x, i, j, k, PMLComp::xz) );
F(i, j, k, PMLComp::y) += dt * (
T_Algo::DownwardDy(Ey, coefs_y, n_coefs_y, i, j, k, PMLComp::yx)
+ T_Algo::DownwardDy(Ey, coefs_y, n_coefs_y, i, j, k, PMLComp::yy)
+ T_Algo::DownwardDy(Ey, coefs_y, n_coefs_y, i, j, k, PMLComp::yz) );
F(i, j, k, PMLComp::z) += dt * (
T_Algo::DownwardDz(Ez, coefs_z, n_coefs_z, i, j, k, PMLComp::zx)
+ T_Algo::DownwardDz(Ez, coefs_z, n_coefs_z, i, j, k, PMLComp::zy)
+ T_Algo::DownwardDz(Ez, coefs_z, n_coefs_z, i, j, k, PMLComp::zz) );
}
);
}
}
#endif // corresponds to ifndef WARPX_DIM_RZ
|