Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_MLEBABecLap.H
Go to the documentation of this file.
1#ifndef AMREX_MLEBABECLAP_H_
2#define AMREX_MLEBABECLAP_H_
3#include <AMReX_Config.H>
4
7#include <AMReX_Array.H>
8#include <limits>
9
10namespace amrex {
11
12// (alpha * a - beta * (del dot b grad)) phi
13
16 : public MLCellABecLap
17{
18public:
19
20 MLEBABecLap () = default;
21 MLEBABecLap (const Vector<Geometry>& a_geom,
22 const Vector<BoxArray>& a_grids,
23 const Vector<DistributionMapping>& a_dmap,
24 const LPInfo& a_info,
25 const Vector<EBFArrayBoxFactory const*>& a_factory,
26 int a_ncomp = 1);
27
28 ~MLEBABecLap () override;
29
30 MLEBABecLap (const MLEBABecLap&) = delete;
34
35 void define (const Vector<Geometry>& a_geom,
36 const Vector<BoxArray>& a_grids,
37 const Vector<DistributionMapping>& a_dmap,
38 const LPInfo& a_info,
39 const Vector<EBFArrayBoxFactory const*>& a_factory,
40 int a_ncomp = 1);
41
42 void setPhiOnCentroid ();
43
44 void setScalars (Real a, Real b);
45 void setACoeffs (int amrlev, const MultiFab& alpha);
46 void setACoeffs (int amrlev, Real alpha);
47
48 void setBCoeffs (int amrlev, const Array<MultiFab const*,AMREX_SPACEDIM>& beta,
49 Location a_beta_loc);
51 {setBCoeffs (amrlev, beta, Location::FaceCenter);}
52
53 void setBCoeffs (int amrlev, Real beta);
54 void setBCoeffs (int amrlev, Vector<Real> const& beta);
55
56 // Tells the solver that EB boundaries have Dirichlet bc's specified by "phi"
57 void setEBDirichlet (int amrlev, const MultiFab& phi, const MultiFab& beta);
58 void setEBDirichlet (int amrlev, const MultiFab& phi, Real beta);
59 void setEBDirichlet (int amrlev, const MultiFab& phi, Vector<Real> const& beta);
60
61 // Tells the solver that EB boundaries have homogeneous Dirichlet bc's
62 void setEBHomogDirichlet (int amrlev, const MultiFab& beta);
63 void setEBHomogDirichlet (int amrlev, Real beta);
64 void setEBHomogDirichlet (int amrlev, Vector<Real> const& beta);
65
66 int getNComp () const override { return m_ncomp; }
67
68 bool needsUpdate () const override {
70 }
71 void update () override;
72
73 std::unique_ptr<FabFactory<FArrayBox> > makeFactory (int amrlev, int mglev) const final;
74
75 bool isCrossStencil () const override { return false; }
76
77 void applyBC (int amrlev, int mglev, MultiFab& in, BCMode bc_mode, StateMode s_mode,
78 const MLMGBndry* bndry=nullptr, bool skip_fillboundary=false) const final;
79 void apply (int amrlev, int mglev, MultiFab& out, MultiFab& in, BCMode bc_mode,
80 StateMode s_mode, const MLMGBndry* bndry=nullptr) const override;
81 void compGrad (int amrlev, const Array<MultiFab*,AMREX_SPACEDIM>& grad,
82 MultiFab& sol, Location loc) const final;
83
84 void prepareForSolve () override;
85 bool isSingular (int amrlev) const override { return m_is_singular[amrlev]; }
86 bool isBottomSingular () const override { return m_is_singular[0]; }
87 void Fapply (int amrlev, int mglev, MultiFab& out, const MultiFab& in) const final;
88 void Fsmooth (int amrlev, int mglev, MultiFab& sol, const MultiFab& rhs, int redblack) const final;
89 void FFlux (int amrlev, const MFIter& mfi,
91 const FArrayBox& sol, Location loc,
92 int face_only=0) const final;
93
94 void normalize (int amrlev, int mglev, MultiFab& mf) const final;
95
96 Real getAScalar () const final { return m_a_scalar; }
97 Real getBScalar () const final { return m_b_scalar; }
98 MultiFab const* getACoeffs (int amrlev, int mglev) const final
99 { return &(m_a_coeffs[amrlev][mglev]); }
100 Array<MultiFab const*,AMREX_SPACEDIM> getBCoeffs (int amrlev, int mglev) const final
101 { return amrex::GetArrOfConstPtrs(m_b_coeffs[amrlev][mglev]); }
102
103 std::unique_ptr<MLLinOp> makeNLinOp (int /*grid_size*/) const final {
104 amrex::Abort("MLABecLaplacian::makeNLinOp: Not implemented");
105 return std::unique_ptr<MLLinOp>{};
106 }
107
108 void restriction (int, int, MultiFab& crse, MultiFab& fine) const final;
109
110 void interpolation (int amrlev, int fmglev, MultiFab& fine, const MultiFab& crse) const final;
111
112 void averageDownSolutionRHS (int camrlev, MultiFab& crse_sol, MultiFab& crse_rhs,
113 const MultiFab& fine_sol, const MultiFab& fine_rhs) final;
114
115 void getEBFluxes (const Vector<MultiFab*>& a_flux,
116 const Vector<MultiFab*>& a_sol) const override;
117
119
120#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
121 [[nodiscard]] std::unique_ptr<Hypre> makeHypre (Hypre::Interface hypre_interface) const override;
122#endif
123
124#if defined(AMREX_USE_PETSC) && (AMREX_SPACEDIM > 1)
125 [[nodiscard]] std::unique_ptr<PETScABecLap> makePETSc () const override;
126#endif
127
128 Real m_a_scalar = std::numeric_limits<Real>::quiet_NaN();
129 Real m_b_scalar = std::numeric_limits<Real>::quiet_NaN();
132
133 bool m_scalars_set = false;
134 bool m_acoef_set = false;
135
136protected:
137
138 int m_ncomp = 1;
139
140 bool m_needs_update = true;
141
142 Location m_beta_loc; // Location of coefficients: face centers or face centroids
143 Location m_phi_loc; // Location of solution variable: cell centers or cell centroids
144
146
149
151
152 mutable int m_is_eb_inhomog;
153
154 //
155 // functions
156 //
157 bool isEBDirichlet () const noexcept { return m_eb_phi[0] != nullptr; }
158
161 const Vector<MultiFab*>& b_eb);
162 void averageDownCoeffs ();
163 void averageDownCoeffsToCoarseAmrLevel (int flev);
164
165 [[nodiscard]] bool supportRobinBC () const noexcept override { return true; }
166};
167
169
170// Tag for applyBC EB version
171template <typename T>
172struct MLMGABCEBTag {
177 T bcloc;
178 amrex::Box bx;
179 amrex::BoundCond bctype;
180 int blen;
181 int comp;
182 int dir;
183 int is_eb;
184 int side;
185
187 amrex::Box const& box() const noexcept { return bx; }
188};
189
191
192}
193
194#endif
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< int const > mask
Definition AMReX_InterpFaceRegister.cpp:93
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
Maintain an identifier for boundary condition types.
Definition AMReX_BoundCond.H:20
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
Interface
Definition AMReX_Hypre.H:21
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:63
Definition AMReX_MLCellABecLap.H:14
bool needsUpdate() const override
Does it need update if it's reused?
Definition AMReX_MLCellABecLap.H:47
typename MLLinOpT< MF >::Location Location
Definition AMReX_MLCellABecLap.H:20
typename MLLinOpT< MF >::BCMode BCMode
Definition AMReX_MLCellLinOp.H:28
typename MLLinOpT< MF >::StateMode StateMode
Definition AMReX_MLCellLinOp.H:29
Definition AMReX_MLEBABecLap.H:17
~MLEBABecLap() override
bool isCrossStencil() const override
Definition AMReX_MLEBABecLap.H:75
MLEBABecLap(MLEBABecLap &&)=delete
void prepareForSolve() override
Definition AMReX_MLEBABecLap.cpp:686
int m_is_eb_inhomog
Definition AMReX_MLEBABecLap.H:152
MLEBABecLap(const MLEBABecLap &)=delete
bool m_acoef_set
Definition AMReX_MLEBABecLap.H:134
Vector< Vector< std::unique_ptr< MultiFab > > > m_eb_b_coeffs
Definition AMReX_MLEBABecLap.H:148
void setPhiOnCentroid()
Definition AMReX_MLEBABecLap.cpp:103
Real getBScalar() const final
Definition AMReX_MLEBABecLap.H:97
Real m_b_scalar
Definition AMReX_MLEBABecLap.H:129
std::unique_ptr< MLLinOp > makeNLinOp(int) const final
Definition AMReX_MLEBABecLap.H:103
Vector< std::unique_ptr< MultiFab > > m_eb_phi
Definition AMReX_MLEBABecLap.H:147
bool supportRobinBC() const noexcept override
Definition AMReX_MLEBABecLap.H:165
void Fapply(int amrlev, int mglev, MultiFab &out, const MultiFab &in) const final
Definition AMReX_MLEBABecLap_F.cpp:12
void setBCoeffs(int amrlev, const Array< MultiFab const *, 3 > &beta, Location a_beta_loc)
Definition AMReX_MLEBABecLap.cpp:141
void restriction(int, int, MultiFab &crse, MultiFab &fine) const final
Definition AMReX_MLEBABecLap.cpp:950
void FFlux(int amrlev, const MFIter &mfi, const Array< FArrayBox *, 3 > &flux, const FArrayBox &sol, Location loc, int face_only=0) const final
Definition AMReX_MLEBABecLap_F.cpp:304
void applyRobinBCTermsCoeffs()
Definition AMReX_MLEBABecLap.cpp:1406
void setACoeffs(int amrlev, const MultiFab &alpha)
Definition AMReX_MLEBABecLap.cpp:125
MultiFab const * getACoeffs(int amrlev, int mglev) const final
Definition AMReX_MLEBABecLap.H:98
void normalize(int amrlev, int mglev, MultiFab &mf) const final
Definition AMReX_MLEBABecLap.cpp:857
Vector< Vector< MultiFab > > m_a_coeffs
Definition AMReX_MLEBABecLap.H:130
void averageDownSolutionRHS(int camrlev, MultiFab &crse_sol, MultiFab &crse_rhs, const MultiFab &fine_sol, const MultiFab &fine_rhs) final
Definition AMReX_MLEBABecLap.cpp:1005
void interpolation(int amrlev, int fmglev, MultiFab &fine, const MultiFab &crse) const final
Definition AMReX_MLEBABecLap.cpp:958
void applyBC(int amrlev, int mglev, MultiFab &in, BCMode bc_mode, StateMode s_mode, const MLMGBndry *bndry=nullptr, bool skip_fillboundary=false) const final
Definition AMReX_MLEBABecLap.cpp:1015
void averageDownCoeffs()
Definition AMReX_MLEBABecLap.cpp:607
bool isEBDirichlet() const noexcept
Definition AMReX_MLEBABecLap.H:157
void Fsmooth(int amrlev, int mglev, MultiFab &sol, const MultiFab &rhs, int redblack) const final
Definition AMReX_MLEBABecLap_F.cpp:157
void getEBFluxes(const Vector< MultiFab * > &a_flux, const Vector< MultiFab * > &a_sol) const override
Definition AMReX_MLEBABecLap.cpp:1333
MLEBABecLap & operator=(const MLEBABecLap &)=delete
Vector< int > m_is_singular
Definition AMReX_MLEBABecLap.H:150
std::unique_ptr< FabFactory< FArrayBox > > makeFactory(int amrlev, int mglev) const final
Definition AMReX_MLEBABecLap.cpp:34
void setEBDirichlet(int amrlev, const MultiFab &phi, const MultiFab &beta)
Definition AMReX_MLEBABecLap.cpp:190
Vector< Vector< iMultiFab > > m_cc_mask
Definition AMReX_MLEBABecLap.H:145
int m_ncomp
Definition AMReX_MLEBABecLap.H:138
void apply(int amrlev, int mglev, MultiFab &out, MultiFab &in, BCMode bc_mode, StateMode s_mode, const MLMGBndry *bndry=nullptr) const override
Definition AMReX_MLEBABecLap.cpp:1288
void averageDownCoeffsToCoarseAmrLevel(int flev)
Definition AMReX_MLEBABecLap.cpp:663
bool needsUpdate() const override
Does it need update if it's reused?
Definition AMReX_MLEBABecLap.H:68
void setEBHomogDirichlet(int amrlev, const MultiFab &beta)
Definition AMReX_MLEBABecLap.cpp:400
void compGrad(int amrlev, const Array< MultiFab *, 3 > &grad, MultiFab &sol, Location loc) const final
Definition AMReX_MLEBABecLap.cpp:731
bool m_needs_update
Definition AMReX_MLEBABecLap.H:140
Real getAScalar() const final
Definition AMReX_MLEBABecLap.H:96
MLEBABecLap()=default
Real m_a_scalar
Definition AMReX_MLEBABecLap.H:128
int getNComp() const override
Return number of components.
Definition AMReX_MLEBABecLap.H:66
bool isSingular(int amrlev) const override
Is it singular on given AMR level?
Definition AMReX_MLEBABecLap.H:85
void setScalars(Real a, Real b)
Definition AMReX_MLEBABecLap.cpp:109
void update() override
Update for reuse.
Definition AMReX_MLEBABecLap.cpp:1297
void averageDownCoeffsSameAmrLevel(int amrlev, Vector< MultiFab > &a, Vector< Array< MultiFab, 3 > > &b, const Vector< MultiFab * > &b_eb)
Definition AMReX_MLEBABecLap.cpp:632
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, int a_ncomp=1)
Definition AMReX_MLEBABecLap.cpp:44
bool isBottomSingular() const override
Is the bottom of MG singular?
Definition AMReX_MLEBABecLap.H:86
Location m_phi_loc
Definition AMReX_MLEBABecLap.H:143
void setBCoeffs(int amrlev, const Array< MultiFab const *, 3 > &beta)
Definition AMReX_MLEBABecLap.H:50
Array< MultiFab const *, 3 > getBCoeffs(int amrlev, int mglev) const final
Definition AMReX_MLEBABecLap.H:100
bool m_scalars_set
Definition AMReX_MLEBABecLap.H:133
Vector< Vector< Array< MultiFab, 3 > > > m_b_coeffs
Definition AMReX_MLEBABecLap.H:131
Location m_beta_loc
Definition AMReX_MLEBABecLap.H:142
Definition AMReX_MLMGBndry.H:12
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
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
std::array< T, N > Array
Definition AMReX_Array.H:25
Definition AMReX_Amr.cpp:49
std::array< T const *, 3 > GetArrOfConstPtrs(const std::array< T, 3 > &a) noexcept
Definition AMReX_Array.H:1016
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
std::unique_ptr< Hypre > makeHypre(const BoxArray &grids, const DistributionMapping &dmap, const Geometry &geom, MPI_Comm comm_, Hypre::Interface interface, const iMultiFab *overset_mask)
Definition AMReX_Hypre.cpp:12
Definition AMReX_Array4.H:61
Definition AMReX_MLLinOp.H:36