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
|
#include "DivEFunctor.H"
#include "Utils/CoarsenIO.H"
#include "Utils/TextMsg.H"
#ifdef WARPX_DIM_RZ
# include "Utils/WarpXAlgorithmSelection.H"
#endif
#include "WarpX.H"
#include <AMReX.H>
#include <AMReX_BoxArray.H>
#include <AMReX_IntVect.H>
#include <AMReX_MultiFab.H>
DivEFunctor::DivEFunctor(const std::array<const amrex::MultiFab* const, 3> arr_mf_src, const int lev,
const amrex::IntVect crse_ratio,
bool convertRZmodes2cartesian, const int ncomp)
: ComputeDiagFunctor(ncomp, crse_ratio), m_arr_mf_src(arr_mf_src), m_lev(lev),
m_convertRZmodes2cartesian(convertRZmodes2cartesian)
{
amrex::ignore_unused(m_arr_mf_src);
}
void
DivEFunctor::operator()(amrex::MultiFab& mf_dst, const int dcomp, const int /*i_buffer*/) const
{
auto& warpx = WarpX::GetInstance();
// Guard cell is set to 1 for generality. However, for a cell-centered
// output Multifab, mf_dst, the guard-cell data is not needed especially considering
// the operations performend in the CoarsenAndInterpolate function.
constexpr int ng = 1;
// For staggered and nodal calculations, divE is computed on the nodes.
// The temporary divE MultiFab is generated to comply with the location of divE.
amrex::IntVect cell_type = amrex::IntVect::TheNodeVector();
#ifdef WARPX_DIM_RZ
// For RZ spectral, all quantities are cell centered.
if (WarpX::maxwell_solver_id == MaxwellSolverAlgo::PSATD)
cell_type = amrex::IntVect::TheCellVector();
#endif
const amrex::BoxArray& ba = amrex::convert(warpx.boxArray(m_lev), cell_type);
amrex::MultiFab divE(ba, warpx.DistributionMap(m_lev), warpx.ncomps, ng );
warpx.ComputeDivE(divE, m_lev);
#ifdef WARPX_DIM_RZ
if (m_convertRZmodes2cartesian) {
// In cylindrical geometry, sum real part of all modes of divE in
// temporary multifab mf_dst_stag, and cell-center it to mf_dst.
WARPX_ALWAYS_ASSERT_WITH_MESSAGE(
nComp()==1,
"The RZ averaging over modes must write into 1 single component");
amrex::MultiFab mf_dst_stag(divE.boxArray(), warpx.DistributionMap(m_lev), 1, divE.nGrowVect());
// Mode 0
amrex::MultiFab::Copy(mf_dst_stag, divE, 0, 0, 1, divE.nGrowVect());
for (int ic=1 ; ic < divE.nComp() ; ic += 2) {
// Real part of all modes > 0
amrex::MultiFab::Add(mf_dst_stag, divE, ic, 0, 1, divE.nGrowVect());
}
CoarsenIO::Coarsen( mf_dst, mf_dst_stag, dcomp, 0, nComp(), 0, m_crse_ratio);
} else {
CoarsenIO::Coarsen( mf_dst, divE, dcomp, 0, nComp(), 0, m_crse_ratio);
}
#else
// In cartesian geometry, coarsen and interpolate from simulation MultiFab, divE,
// to output diagnostic MultiFab, mf_dst.
CoarsenIO::Coarsen( mf_dst, divE, dcomp, 0, nComp(), 0, m_crse_ratio);
amrex::ignore_unused(m_convertRZmodes2cartesian);
#endif
}
|