Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_MLNodeLaplacian.H
Go to the documentation of this file.
1#ifndef AMREX_ML_NODE_LAPLACIAN_H_
2#define AMREX_ML_NODE_LAPLACIAN_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_iMultiFab.H>
6#include <AMReX_MLNodeLinOp.H>
7#include <AMReX_MultiFab.H>
8#include <AMReX_Vector.H>
9
10// Feature macro: when defined, MLNodeLaplacian::updateVelocity and
11// MLNodeLaplacian::getFluxes correctly route through the anisotropic
12// (harmonic-average) mknewu kernel when setMapped(true) was used. Prior
13// to this, both functions only used sigma component 0 (sigma_x) for all
14// velocity/flux components, producing the wrong correction for
15// non-identity mesh mappings. EB is not covered by this feature (see
16// comments in AMReX_MLNodeLaplacian_misc.cpp).
17#define AMREX_MLNODELAP_HAS_MKNEWU_HA 1
18
19namespace amrex {
20
21// del dot (sigma grad phi) = rhs
22// where phi and rhs are nodal, and sigma is cell-centered.
23
26 : public MLNodeLinOp
27{
28public :
29
30 MLNodeLaplacian () = default;
31 MLNodeLaplacian (const Vector<Geometry>& a_geom,
32 const Vector<BoxArray>& a_grids,
33 const Vector<DistributionMapping>& a_dmap,
34 const LPInfo& a_info = LPInfo(),
35 const Vector<FabFactory<FArrayBox> const*>& a_factory = {},
36 Real a_const_sigma = Real(0.0));
37#ifdef AMREX_USE_EB
38 MLNodeLaplacian (const Vector<Geometry>& a_geom,
39 const Vector<BoxArray>& a_grids,
40 const Vector<DistributionMapping>& a_dmap,
41 const LPInfo& a_info,
42 const Vector<EBFArrayBoxFactory const*>& a_factory,
43 Real a_const_sigma = Real(0.0));
44#endif
45 ~MLNodeLaplacian () override = default;
46
51
52 void define (const Vector<Geometry>& a_geom,
53 const Vector<BoxArray>& a_grids,
54 const Vector<DistributionMapping>& a_dmap,
55 const LPInfo& a_info = LPInfo(),
56 const Vector<FabFactory<FArrayBox> const*>& a_factory = {},
57 Real a_const_sigma = Real(0.0));
58
59#ifdef AMREX_USE_EB
60 void define (const Vector<Geometry>& a_geom,
61 const Vector<BoxArray>& a_grids,
62 const Vector<DistributionMapping>& a_dmap,
63 const LPInfo& a_info,
64 const Vector<EBFArrayBoxFactory const*>& a_factory,
65 Real a_const_sigma = Real(0.0));
66#endif
67
68 std::string name () const override { return std::string("MLNodeLaplacian"); }
69
70 void setRZCorrection (bool rz) noexcept { m_is_rz = rz; }
71
72 void setNormalizationThreshold (Real t) noexcept { m_normalization_threshold = t; }
73
74 void setSigma (int amrlev, const MultiFab& a_sigma);
75
76 void compDivergence (const Vector<MultiFab*>& rhs, const Vector<MultiFab*>& vel);
77
78 void compRHS (const Vector<MultiFab*>& rhs, const Vector<MultiFab*>& vel,
79 const Vector<const MultiFab*>& rhnd,
80 const Vector<MultiFab*>& rhcc);
81
82 void updateVelocity (const Vector<MultiFab*>& vel, const Vector<MultiFab const*>& sol) const;
83
84 void compSyncResidualCoarse (MultiFab& sync_resid, const MultiFab& phi,
85 const MultiFab& vold, const MultiFab* rhcc,
86 const BoxArray& fine_grids, const IntVect& ref_ratio);
87
88 void compSyncResidualFine (MultiFab& sync_resid, const MultiFab& phi, const MultiFab& vold,
89 const MultiFab* rhcc);
90
91 void setGaussSeidel (bool flag) noexcept { m_use_gauss_seidel = flag; }
92 void setHarmonicAverage (bool flag) noexcept { m_use_harmonic_average = flag; }
93
94 void setMapped (bool flag) noexcept { m_use_mapped = flag; }
95
97 if (m_const_sigma == Real(0.0)) { m_coarsening_strategy = cs; }
98 }
99
104
105 void restriction (int amrlev, int cmglev, MultiFab& crse, MultiFab& fine) const final;
106 void interpolation (int amrlev, int fmglev, MultiFab& fine, const MultiFab& crse) const final;
107 void averageDownSolutionRHS (int camrlev, MultiFab& crse_sol, MultiFab& crse_rhs,
108 const MultiFab& fine_sol, const MultiFab& fine_rhs) final;
109
110 void reflux (int crse_amrlev,
111 MultiFab& res, const MultiFab& crse_sol, const MultiFab& crse_rhs,
112 MultiFab& fine_res, MultiFab& fine_sol, const MultiFab& fine_rhs) const final;
113
114 void prepareForSolve () final;
115 void Fapply (int amrlev, int mglev, MultiFab& out, const MultiFab& in) const final;
116 void Fsmooth (int amrlev, int mglev, MultiFab& sol, const MultiFab& rhs) const final;
117 void normalize (int amrlev, int mglev, MultiFab& mf) const final;
118
119 void fixUpResidualMask (int amrlev, iMultiFab& resmsk) final;
120
121 void getFluxes (const Vector<Array<MultiFab*,AMREX_SPACEDIM> >& /*a_flux*/,
122 const Vector<MultiFab*>& /*a_sol*/,
123 Location /*a_loc*/) const final {
124 amrex::Abort("MLNodeLaplacian::getFluxes: How did we get here?");
125 }
126 void getFluxes (const Vector<MultiFab*>& a_flux,
127 const Vector<MultiFab*>& a_sol) const final;
128 void unimposeNeumannBC (int amrlev, MultiFab& rhs) const final;
129 Vector<Real> getSolvabilityOffset (int amrlev, int mglev,
130 MultiFab const& rhs) const override;
131 void fixSolvabilityByOffset (int amrlev, int mglev, MultiFab& rhs,
132 Vector<Real> const& offset) const override;
133
134 void compGrad (int /*amrlev*/, const Array<MultiFab*,AMREX_SPACEDIM>& /*grad*/,
135 MultiFab& /*sol*/, Location /*loc*/) const final {
136 amrex::Abort("MLNodeLaplacian::compGrad: How did we get here?");
137 }
138 void compGrad (int amrlev, MultiFab& grad, MultiFab& sol) const;
139
140 void averageDownCoeffs ();
141 void averageDownCoeffsToCoarseAmrLevel (int flev);
142 void averageDownCoeffsSameAmrLevel (int amrlev);
143
144 void restrictInteriorNodes (int camrlev, MultiFab& crhs, MultiFab& frhs) const;
145
146 void FillBoundaryCoeff (MultiFab& sigma, const Geometry& geom);
147
148 void buildStencil ();
149
150#ifdef AMREX_USE_EB
151 void buildIntegral ();
152 void buildSurfaceIntegral ();
153
154 // Tells the solver that EB boundaries have inflow specified by "eb_vel"
155 void setEBInflowVelocity (int amrlev, const MultiFab& eb_vel);
156
157#endif
158
159#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
160 void fillIJMatrix (MFIter const& mfi,
162 Array4<int const> const& lid,
163 HypreNodeLap::Int* ncols,
164 HypreNodeLap::Int* cols,
165 Real* mat) const override;
166
167#ifdef AMREX_USE_GPU
168 void fillIJMatrix_gpu (MFIter const& mfi,
170 Array4<int const> const& lid,
171 HypreNodeLap::Int* ncols,
172 HypreNodeLap::Int* cols,
173 Real* mat) const;
174#endif
175
176 void fillIJMatrix_cpu (MFIter const& mfi,
178 Array4<int const> const& lid,
179 HypreNodeLap::Int* ncols,
180 HypreNodeLap::Int* cols,
181 Real* mat) const;
182
183 void fillRHS (MFIter const& mfi, Array4<int const> const& lid,
184 Real* rhs, Array4<Real const> const& bfab) const override;
185#endif
186
187protected:
188
189 void resizeMultiGrid (int new_size) final;
190
191private:
192
193 int m_is_rz = 0;
194
195 Real m_const_sigma = Real(0.0);
196 Vector<Vector<Array<std::unique_ptr<MultiFab>,AMREX_SPACEDIM> > > m_sigma;
198 Vector<std::unique_ptr<MultiFab> > m_nosigma_stencil;
199 Vector<Vector<Real> > m_s0_norm0;
200
201 Real m_normalization_threshold = Real(1.e-8);
202
203#ifdef AMREX_USE_EB
204 // they could be MultiCutFab
206 bool m_integral_built = false;
207
208 Vector<std::unique_ptr<MultiFab> > m_surface_integral;
209 bool m_surface_integral_built = false;
210 bool m_build_surface_integral = false;
211 Vector<std::unique_ptr<MultiFab> > m_eb_vel_dot_n;
212#endif
213
214 bool m_use_gauss_seidel = true;
215 bool m_use_harmonic_average = false;
216 bool m_use_mapped = false;
217
218 void checkPoint (std::string const& file_name) const final;
219};
220
221}
222
223#endif
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1139
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:564
Definition AMReX_FabFactory.H:50
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
HYPRE_Int Int
Definition AMReX_HypreNodeLap.H:59
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
Definition AMReX_MLNodeLaplacian.H:27
void getFluxes(const Vector< Array< MultiFab *, 3 > > &, const Vector< MultiFab * > &, Location) const final
Definition AMReX_MLNodeLaplacian.H:121
void restrictInteriorNodes(int camrlev, MultiFab &crhs, MultiFab &frhs) const
Definition AMReX_MLNodeLaplacian.cpp:751
void normalize(int amrlev, int mglev, MultiFab &mf) const final
Definition AMReX_MLNodeLaplacian.cpp:845
void FillBoundaryCoeff(MultiFab &sigma, const Geometry &geom)
Definition AMReX_MLNodeLaplacian.cpp:405
void define(const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const LPInfo &a_info=LPInfo(), const Vector< FabFactory< FArrayBox > const * > &a_factory={}, Real a_const_sigma=Real(0.0))
Definition AMReX_MLNodeLaplacian.cpp:39
void setHarmonicAverage(bool flag) noexcept
Definition AMReX_MLNodeLaplacian.H:92
void averageDownSolutionRHS(int camrlev, MultiFab &crse_sol, MultiFab &crse_rhs, const MultiFab &fine_sol, const MultiFab &fine_rhs) final
Definition AMReX_MLNodeLaplacian.cpp:736
std::string name() const override
Definition AMReX_MLNodeLaplacian.H:68
void prepareForSolve() final
Definition AMReX_MLNodeLaplacian.cpp:453
void fixUpResidualMask(int amrlev, iMultiFab &resmsk) final
Definition AMReX_MLNodeLaplacian.cpp:431
MLNodeLaplacian(const MLNodeLaplacian &)=delete
void compSyncResidualFine(MultiFab &sync_resid, const MultiFab &phi, const MultiFab &vold, const MultiFab *rhcc)
Definition AMReX_MLNodeLaplacian_sync.cpp:335
void Fapply(int amrlev, int mglev, MultiFab &out, const MultiFab &in) const final
Definition AMReX_MLNodeLaplacian_misc.cpp:199
void setSigma(int amrlev, const MultiFab &a_sigma)
Definition AMReX_MLNodeLaplacian.cpp:379
void setNormalizationThreshold(Real t) noexcept
Definition AMReX_MLNodeLaplacian.H:72
void buildStencil()
Definition AMReX_MLNodeLaplacian_sten.cpp:18
void compRHS(const Vector< MultiFab * > &rhs, const Vector< MultiFab * > &vel, const Vector< const MultiFab * > &rhnd, const Vector< MultiFab * > &rhcc)
Definition AMReX_MLNodeLaplacian_misc.cpp:1050
void setEBInflowVelocity(int amrlev, const MultiFab &eb_vel)
Definition AMReX_MLNodeLaplacian.cpp:1030
void buildIntegral()
Definition AMReX_MLNodeLaplacian_eb.cpp:15
void updateVelocity(const Vector< MultiFab * > &vel, const Vector< MultiFab const * > &sol) const
Definition AMReX_MLNodeLaplacian_misc.cpp:715
~MLNodeLaplacian() override=default
Vector< Real > getSolvabilityOffset(int amrlev, int mglev, MultiFab const &rhs) const override
Definition AMReX_MLNodeLaplacian.cpp:167
MLNodeLaplacian(MLNodeLaplacian &&)=delete
void setMapped(bool flag) noexcept
Definition AMReX_MLNodeLaplacian.H:94
void setGaussSeidel(bool flag) noexcept
Definition AMReX_MLNodeLaplacian.H:91
void interpolation(int amrlev, int fmglev, MultiFab &fine, const MultiFab &crse) const final
Definition AMReX_MLNodeLaplacian.cpp:586
void compGrad(int, const Array< MultiFab *, 3 > &, MultiFab &, Location) const final
Definition AMReX_MLNodeLaplacian.H:134
BottomSolver getDefaultBottomSolver() const final
Definition AMReX_MLNodeLaplacian.H:100
void restriction(int amrlev, int cmglev, MultiFab &crse, MultiFab &fine) const final
Definition AMReX_MLNodeLaplacian.cpp:472
void Fsmooth(int amrlev, int mglev, MultiFab &sol, const MultiFab &rhs) const final
Definition AMReX_MLNodeLaplacian_misc.cpp:347
MLNodeLaplacian & operator=(const MLNodeLaplacian &)=delete
void reflux(int crse_amrlev, MultiFab &res, const MultiFab &crse_sol, const MultiFab &crse_rhs, MultiFab &fine_res, MultiFab &fine_sol, const MultiFab &fine_rhs) const final
Definition AMReX_MLNodeLaplacian_sync.cpp:613
void averageDownCoeffs()
Definition AMReX_MLNodeLaplacian_misc.cpp:17
void buildSurfaceIntegral()
Definition AMReX_MLNodeLaplacian_eb.cpp:84
void setRZCorrection(bool rz) noexcept
Definition AMReX_MLNodeLaplacian.H:70
void compDivergence(const Vector< MultiFab * > &rhs, const Vector< MultiFab * > &vel)
Definition AMReX_MLNodeLaplacian_misc.cpp:1044
void fixSolvabilityByOffset(int amrlev, int mglev, MultiFab &rhs, Vector< Real > const &offset) const override
Definition AMReX_MLNodeLaplacian.cpp:282
void averageDownCoeffsToCoarseAmrLevel(int flev)
Definition AMReX_MLNodeLaplacian_misc.cpp:86
void unimposeNeumannBC(int amrlev, MultiFab &rhs) const final
Definition AMReX_MLNodeLaplacian.cpp:145
void resizeMultiGrid(int new_size) final
Definition AMReX_MLNodeLaplacian.cpp:121
void setCoarseningStrategy(CoarseningStrategy cs) noexcept
Definition AMReX_MLNodeLaplacian.H:96
void compSyncResidualCoarse(MultiFab &sync_resid, const MultiFab &phi, const MultiFab &vold, const MultiFab *rhcc, const BoxArray &fine_grids, const IntVect &ref_ratio)
Definition AMReX_MLNodeLaplacian_sync.cpp:16
void averageDownCoeffsSameAmrLevel(int amrlev)
Definition AMReX_MLNodeLaplacian_misc.cpp:102
Definition AMReX_MLNodeLinOp.H:17
CoarseningStrategy
Definition AMReX_MLNodeLinOp.H:20
CoarseningStrategy m_coarsening_strategy
Definition AMReX_MLNodeLinOp.H:149
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:50
BottomSolver
Definition AMReX_MLLinOp.H:32
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
A multidimensional array accessor.
Definition AMReX_Array4.H:283
Definition AMReX_MLLinOp.H:37
Location
Definition AMReX_MLLinOp.H:91