Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_MLEBNodeFDLaplacian.H
Go to the documentation of this file.
1#ifndef AMREX_MLEBNODEFDLAPLACIAN_H_
2#define AMREX_MLEBNODEFDLAPLACIAN_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_Array.H>
6#ifdef AMREX_USE_EB
8#endif
9#include <AMReX_MLNodeLinOp.H>
10
11#include <limits>
12
13namespace amrex {
14
15// Although the class has EB in the name, it works for non-EB build too.
16//
17// del dot (sigma grad phi) = rhs, for non-RZ
18// where phi and rhs are nodal multifab, and sigma is a tensor constant
19// with only diagonal components. The EB is assumed to be Dirichlet.
20//
21// del dot (sigma grad phi) - alpha/r^2 phi = rhs, for RZ where alpha is a
22// scalar constant that is zero by default.
23//
24// New feature: for non-RZ, sigma can also be a single-component
25// cell-centered multifab.
26
29 : public MLNodeLinOp
30{
31public:
32
33 MLEBNodeFDLaplacian () = default;
34
35#ifdef AMREX_USE_EB
37 const Vector<BoxArray>& a_grids,
38 const Vector<DistributionMapping>& a_dmap,
39 const LPInfo& a_info,
40 const Vector<EBFArrayBoxFactory const*>& a_factory);
41#endif
42
44 const Vector<BoxArray>& a_grids,
45 const Vector<DistributionMapping>& a_dmap,
46 const LPInfo& a_info);
47
48 ~MLEBNodeFDLaplacian () override = default;
49
54
55 void setSigma (Array<Real,AMREX_SPACEDIM> const& a_sigma) noexcept;
56
57 void setSigma (int amrlev, MultiFab const& a_sigma);
58
59 void setRZ (bool flag);
60
61 void setAlpha (Real a_alpha);
62
63#ifdef AMREX_USE_EB
64
65 // Phi on EB
66 void setEBDirichlet (Real a_phi_eb);
67 //
68 template <typename F>
69 std::enable_if_t<IsCallableR<Real,F,AMREX_D_DECL(Real,Real,Real)>::value>
70 setEBDirichlet (F const& f);
71
72 void define (const Vector<Geometry>& a_geom,
73 const Vector<BoxArray>& a_grids,
74 const Vector<DistributionMapping>& a_dmap,
75 const LPInfo& a_info,
76 const Vector<EBFArrayBoxFactory const*>& a_factory);
77
78 [[nodiscard]] std::unique_ptr<FabFactory<FArrayBox> > makeFactory (int amrlev, int mglev) const final;
79
80 [[nodiscard]] bool scaleRHS (int amrlev, MultiFab* rhs) const final;
81
82#endif
83
84 void define (const Vector<Geometry>& a_geom,
85 const Vector<BoxArray>& a_grids,
86 const Vector<DistributionMapping>& a_dmap,
87 const LPInfo& a_info);
88
89 [[nodiscard]] std::string name () const override { return std::string("MLEBNodeFDLaplacian"); }
90
91 void restriction (int amrlev, int cmglev, MultiFab& crse, MultiFab& fine) const final;
92 void interpolation (int amrlev, int fmglev, MultiFab& fine, const MultiFab& crse) const final;
93
94 [[nodiscard]] bool needsUpdate () const override {
95 return (m_needs_update || MLNodeLinOp::needsUpdate());
96 }
97 void update () override;
98
99 void prepareForSolve () final;
100 void Fapply (int amrlev, int mglev, MultiFab& out, const MultiFab& in) const final;
101 void Fsmooth (int amrlev, int mglev, MultiFab& sol, const MultiFab& rhs) const final;
102 void normalize (int amrlev, int mglev, MultiFab& mf) const final;
103
104 void fixUpResidualMask (int amrlev, iMultiFab& resmsk) final;
105
106 [[nodiscard]] bool isSingular (int) const final { return false; }
107 [[nodiscard]] bool isBottomSingular () const final { return false; }
108
109 void compGrad (int amrlev, const Array<MultiFab*,AMREX_SPACEDIM>& grad,
110 MultiFab& sol, Location /*loc*/) const override;
111
112#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
113 void fillIJMatrix (MFIter const& mfi,
115 Array4<int const> const& lid,
116 HypreNodeLap::Int* ncols,
117 HypreNodeLap::Int* cols,
118 Real* mat) const override;
119
120 void fillRHS (MFIter const& mfi,
121 Array4<int const> const& lid,
122 Real* rhs,
123 Array4<Real const> const& bfab) const override;
124#endif
125
126 void postSolve (Vector<MultiFab*> const& sol) const override;
127
128private:
129 GpuArray<Real,AMREX_SPACEDIM> m_sigma{{AMREX_D_DECL(1_rt,1_rt,1_rt)}};
130 Vector<Vector<std::unique_ptr<MultiFab>>> m_sigma_mf;
131 bool m_has_sigma_mf = false;
132 bool m_needs_update = true;
133 Real m_s_phi_eb = std::numeric_limits<Real>::lowest();
134 Vector<MultiFab> m_phi_eb;
135 int m_rz = false;
136 Real m_rz_alpha = 0._rt;
137
138 void update_sigma ();
139};
140
141#ifdef AMREX_USE_EB
142
143template <typename F>
144std::enable_if_t<IsCallableR<Real,F,AMREX_D_DECL(Real,Real,Real)>::value>
146{
147 m_phi_eb.resize(m_num_amr_levels);
148 for (int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev) {
149 auto const* factory = dynamic_cast<EBFArrayBoxFactory const*>(m_factory[amrlev][0].get());
150 if (factory) {
151 Geometry const& geom = m_geom[amrlev][0];
152 auto const problo = geom.ProbLoArray();
153 auto const cellsize = geom.CellSizeArray();
154 if (m_phi_eb[amrlev].empty()) {
155 m_phi_eb[amrlev].define(amrex::convert(m_grids[amrlev][0],IntVect(1)),
156 m_dmap[amrlev][0], 1, 1);
157 m_phi_eb[amrlev].setVal(0.0);
158 }
159 auto const& flags = factory->getMultiEBCellFlagFab();
160 auto const& levset = factory->getLevelSet();
161#ifdef AMREX_USE_OMP
162#pragma omp parallel if (Gpu::notInLaunchRegion())
163#endif
164 for (MFIter mfi(m_phi_eb[amrlev],TilingIfNotGPU()); mfi.isValid(); ++mfi)
165 {
166 const Box& ndbx = mfi.growntilebox();
167 const auto& flag = flags[mfi];
168 if (flag.getType() != FabType::regular) {
169 Array4<Real const> const lstarr = levset.const_array(mfi);
170 Array4<Real> const& phi = m_phi_eb[amrlev].array(mfi);
171 AMREX_HOST_DEVICE_FOR_3D(ndbx, i, j, k,
172 {
173 if (lstarr(i,j,k) >= Real(0.0)) {
174 phi(i,j,k) = f(AMREX_D_DECL(problo[0]+i*cellsize[0],
175 problo[1]+j*cellsize[1],
176 problo[2]+k*cellsize[2]));
177 }
178 });
179 }
180 }
181 }
182 }
183}
184
185#endif
186
187}
188
189#endif
#define AMREX_HOST_DEVICE_FOR_3D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:106
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
GpuArray< Real, 3 > CellSizeArray() const noexcept
Definition AMReX_CoordSys.H:76
Definition AMReX_EBFabFactory.H:25
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
GpuArray< Real, 3 > ProbLoArray() const noexcept
Definition AMReX_Geometry.H:187
HYPRE_Int Int
Definition AMReX_HypreNodeLap.H:36
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:85
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
Definition AMReX_MLEBNodeFDLaplacian.H:30
void restriction(int amrlev, int cmglev, MultiFab &crse, MultiFab &fine) const final
Definition AMReX_MLEBNodeFDLaplacian.cpp:168
std::string name() const override
Definition AMReX_MLEBNodeFDLaplacian.H:89
bool isBottomSingular() const final
Is the bottom of MG singular?
Definition AMReX_MLEBNodeFDLaplacian.H:107
MLEBNodeFDLaplacian & operator=(const MLEBNodeFDLaplacian &)=delete
bool needsUpdate() const override
Does it need update if it's reused?
Definition AMReX_MLEBNodeFDLaplacian.H:94
void compGrad(int amrlev, const Array< MultiFab *, 3 > &grad, MultiFab &sol, Location) const override
Definition AMReX_MLEBNodeFDLaplacian.cpp:604
void Fsmooth(int amrlev, int mglev, MultiFab &sol, const MultiFab &rhs) const final
Definition AMReX_MLEBNodeFDLaplacian.cpp:487
void interpolation(int amrlev, int fmglev, MultiFab &fine, const MultiFab &crse) const final
Definition AMReX_MLEBNodeFDLaplacian.cpp:215
void setAlpha(Real a_alpha)
Definition AMReX_MLEBNodeFDLaplacian.cpp:68
void setEBDirichlet(Real a_phi_eb)
Definition AMReX_MLEBNodeFDLaplacian.cpp:80
MLEBNodeFDLaplacian(const MLEBNodeFDLaplacian &)=delete
bool scaleRHS(int amrlev, MultiFab *rhs) const final
Definition AMReX_MLEBNodeFDLaplacian.cpp:328
void setRZ(bool flag)
Definition AMReX_MLEBNodeFDLaplacian.cpp:58
void define(const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const LPInfo &a_info, const Vector< EBFArrayBoxFactory const * > &a_factory)
Definition AMReX_MLEBNodeFDLaplacian.cpp:86
void Fapply(int amrlev, int mglev, MultiFab &out, const MultiFab &in) const final
Definition AMReX_MLEBNodeFDLaplacian.cpp:364
void normalize(int amrlev, int mglev, MultiFab &mf) const final
Definition AMReX_MLEBNodeFDLaplacian.cpp:592
void setSigma(Array< Real, 3 > const &a_sigma) noexcept
Definition AMReX_MLEBNodeFDLaplacian.cpp:36
void update() override
Update for reuse.
Definition AMReX_MLEBNodeFDLaplacian.cpp:757
MLEBNodeFDLaplacian(MLEBNodeFDLaplacian &&)=delete
void prepareForSolve() final
Definition AMReX_MLEBNodeFDLaplacian.cpp:259
std::unique_ptr< FabFactory< FArrayBox > > makeFactory(int amrlev, int mglev) const final
Definition AMReX_MLEBNodeFDLaplacian.cpp:154
~MLEBNodeFDLaplacian() override=default
void fixUpResidualMask(int amrlev, iMultiFab &resmsk) final
Definition AMReX_MLEBNodeFDLaplacian.cpp:598
void postSolve(Vector< MultiFab * > const &sol) const override
Definition AMReX_MLEBNodeFDLaplacian.cpp:720
bool isSingular(int) const final
Is it singular on given AMR level?
Definition AMReX_MLEBNodeFDLaplacian.H:106
Vector< Vector< std::unique_ptr< FabFactory< FAB > > > > m_factory
Definition AMReX_MLLinOp.H:618
virtual bool needsUpdate() const
Does it need update if it's reused?
Definition AMReX_MLLinOp.H:271
Vector< Vector< BoxArray > > m_grids
Definition AMReX_MLLinOp.H:616
Vector< Vector< DistributionMapping > > m_dmap
Definition AMReX_MLLinOp.H:617
Vector< Vector< Geometry > > m_geom
first Vector is for amr level and second is mg level
Definition AMReX_MLLinOp.H:615
LinOpEnumType::Location Location
Definition AMReX_MLLinOp.H:118
int m_num_amr_levels
Definition AMReX_MLLinOp.H:600
Definition AMReX_MLNodeLinOp.H:17
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
A Collection of IArrayBoxes.
Definition AMReX_iMultiFab.H:34
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
std::array< T, N > Array
Definition AMReX_Array.H:26
Definition AMReX_Amr.cpp:49
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
bool TilingIfNotGPU() noexcept
Definition AMReX_MFIter.H:12
A multidimensional array accessor.
Definition AMReX_Array4.H:224
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:41
Test if a given type T is callable with arguments of type Args...
Definition AMReX_TypeTraits.H:219
Definition AMReX_MLLinOp.H:36