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
10namespace amrex {
11
12// del dot (sigma grad phi) = rhs
13// where phi and rhs are nodal, and sigma is cell-centered.
14
17 : public MLNodeLinOp
18{
19public :
20
21 MLNodeLaplacian () = default;
22 MLNodeLaplacian (const Vector<Geometry>& a_geom,
23 const Vector<BoxArray>& a_grids,
24 const Vector<DistributionMapping>& a_dmap,
25 const LPInfo& a_info = LPInfo(),
26 const Vector<FabFactory<FArrayBox> const*>& a_factory = {},
27 Real a_const_sigma = Real(0.0));
28#ifdef AMREX_USE_EB
29 MLNodeLaplacian (const Vector<Geometry>& a_geom,
30 const Vector<BoxArray>& a_grids,
31 const Vector<DistributionMapping>& a_dmap,
32 const LPInfo& a_info,
33 const Vector<EBFArrayBoxFactory const*>& a_factory,
34 Real a_const_sigma = Real(0.0));
35#endif
36 ~MLNodeLaplacian () override = default;
37
42
43 void define (const Vector<Geometry>& a_geom,
44 const Vector<BoxArray>& a_grids,
45 const Vector<DistributionMapping>& a_dmap,
46 const LPInfo& a_info = LPInfo(),
47 const Vector<FabFactory<FArrayBox> const*>& a_factory = {},
48 Real a_const_sigma = Real(0.0));
49
50#ifdef AMREX_USE_EB
51 void define (const Vector<Geometry>& a_geom,
52 const Vector<BoxArray>& a_grids,
53 const Vector<DistributionMapping>& a_dmap,
54 const LPInfo& a_info,
55 const Vector<EBFArrayBoxFactory const*>& a_factory,
56 Real a_const_sigma = Real(0.0));
57#endif
58
59 std::string name () const override { return std::string("MLNodeLaplacian"); }
60
61 void setRZCorrection (bool rz) noexcept { m_is_rz = rz; }
62
64
65 void setSigma (int amrlev, const MultiFab& a_sigma);
66
67 void compDivergence (const Vector<MultiFab*>& rhs, const Vector<MultiFab*>& vel);
68
69 void compRHS (const Vector<MultiFab*>& rhs, const Vector<MultiFab*>& vel,
70 const Vector<const MultiFab*>& rhnd,
71 const Vector<MultiFab*>& rhcc);
72
73 void updateVelocity (const Vector<MultiFab*>& vel, const Vector<MultiFab const*>& sol) const;
74
75 void compSyncResidualCoarse (MultiFab& sync_resid, const MultiFab& phi,
76 const MultiFab& vold, const MultiFab* rhcc,
77 const BoxArray& fine_grids, const IntVect& ref_ratio);
78
79 void compSyncResidualFine (MultiFab& sync_resid, const MultiFab& phi, const MultiFab& vold,
80 const MultiFab* rhcc);
81
82 void setGaussSeidel (bool flag) noexcept { m_use_gauss_seidel = flag; }
83 void setHarmonicAverage (bool flag) noexcept { m_use_harmonic_average = flag; }
84
85 void setMapped (bool flag) noexcept { m_use_mapped = flag; }
86
88 if (m_const_sigma == Real(0.0)) { m_coarsening_strategy = cs; }
89 }
90
95
96 void restriction (int amrlev, int cmglev, MultiFab& crse, MultiFab& fine) const final;
97 void interpolation (int amrlev, int fmglev, MultiFab& fine, const MultiFab& crse) const final;
98 void averageDownSolutionRHS (int camrlev, MultiFab& crse_sol, MultiFab& crse_rhs,
99 const MultiFab& fine_sol, const MultiFab& fine_rhs) final;
100
101 void reflux (int crse_amrlev,
102 MultiFab& res, const MultiFab& crse_sol, const MultiFab& crse_rhs,
103 MultiFab& fine_res, MultiFab& fine_sol, const MultiFab& fine_rhs) const final;
104
105 void prepareForSolve () final;
106 void Fapply (int amrlev, int mglev, MultiFab& out, const MultiFab& in) const final;
107 void Fsmooth (int amrlev, int mglev, MultiFab& sol, const MultiFab& rhs) const final;
108 void normalize (int amrlev, int mglev, MultiFab& mf) const final;
109
110 void fixUpResidualMask (int amrlev, iMultiFab& resmsk) final;
111
112 void getFluxes (const Vector<Array<MultiFab*,AMREX_SPACEDIM> >& /*a_flux*/,
113 const Vector<MultiFab*>& /*a_sol*/,
114 Location /*a_loc*/) const final {
115 amrex::Abort("MLNodeLaplacian::getFluxes: How did we get here?");
116 }
117 void getFluxes (const Vector<MultiFab*>& a_flux,
118 const Vector<MultiFab*>& a_sol) const final;
119 void unimposeNeumannBC (int amrlev, MultiFab& rhs) const final;
120 Vector<Real> getSolvabilityOffset (int amrlev, int mglev,
121 MultiFab const& rhs) const override;
122 void fixSolvabilityByOffset (int amrlev, int mglev, MultiFab& rhs,
123 Vector<Real> const& offset) const override;
124
125 void compGrad (int /*amrlev*/, const Array<MultiFab*,AMREX_SPACEDIM>& /*grad*/,
126 MultiFab& /*sol*/, Location /*loc*/) const final {
127 amrex::Abort("MLNodeLaplacian::compGrad: How did we get here?");
128 }
129 void compGrad (int amrlev, MultiFab& grad, MultiFab& sol) const;
130
131 void averageDownCoeffs ();
132 void averageDownCoeffsToCoarseAmrLevel (int flev);
133 void averageDownCoeffsSameAmrLevel (int amrlev);
134
135 void restrictInteriorNodes (int camrlev, MultiFab& crhs, MultiFab& frhs) const;
136
137 void FillBoundaryCoeff (MultiFab& sigma, const Geometry& geom);
138
139 void buildStencil ();
140
141#ifdef AMREX_USE_EB
142 void buildIntegral ();
143 void buildSurfaceIntegral ();
144
145 // Tells the solver that EB boundaries have inflow specified by "eb_vel"
146 void setEBInflowVelocity (int amrlev, const MultiFab& eb_vel);
147
148#endif
149
150#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
151 void fillIJMatrix (MFIter const& mfi,
153 Array4<int const> const& lid,
154 HypreNodeLap::Int* ncols,
155 HypreNodeLap::Int* cols,
156 Real* mat) const override;
157
158#ifdef AMREX_USE_GPU
159 void fillIJMatrix_gpu (MFIter const& mfi,
161 Array4<int const> const& lid,
162 HypreNodeLap::Int* ncols,
163 HypreNodeLap::Int* cols,
164 Real* mat) const;
165#endif
166
167 void fillIJMatrix_cpu (MFIter const& mfi,
169 Array4<int const> const& lid,
170 HypreNodeLap::Int* ncols,
171 HypreNodeLap::Int* cols,
172 Real* mat) const;
173
174 void fillRHS (MFIter const& mfi, Array4<int const> const& lid,
175 Real* rhs, Array4<Real const> const& bfab) const override;
176#endif
177
178protected:
179
180 void resizeMultiGrid (int new_size) final;
181
182private:
183
184 int m_is_rz = 0;
185
191
193
194#ifdef AMREX_USE_EB
195 // they could be MultiCutFab
197 bool m_integral_built = false;
198
203#endif
204
207 bool m_use_mapped = false;
208
209 void checkPoint (std::string const& file_name) const final;
210};
211
212}
213
214#endif
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1089
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:567
Definition AMReX_FabFactory.H:50
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
HYPRE_Int Int
Definition AMReX_HypreNodeLap.H:36
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:63
Definition AMReX_MLNodeLaplacian.H:18
void getFluxes(const Vector< Array< MultiFab *, 3 > > &, const Vector< MultiFab * > &, Location) const final
Definition AMReX_MLNodeLaplacian.H:112
Real m_normalization_threshold
Definition AMReX_MLNodeLaplacian.H:192
void restrictInteriorNodes(int camrlev, MultiFab &crhs, MultiFab &frhs) const
Definition AMReX_MLNodeLaplacian.cpp:751
bool m_use_gauss_seidel
Definition AMReX_MLNodeLaplacian.H:205
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
bool m_use_mapped
Definition AMReX_MLNodeLaplacian.H:207
void setHarmonicAverage(bool flag) noexcept
Definition AMReX_MLNodeLaplacian.H:83
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:59
void prepareForSolve() final
Definition AMReX_MLNodeLaplacian.cpp:453
void fixUpResidualMask(int amrlev, iMultiFab &resmsk) final
Definition AMReX_MLNodeLaplacian.cpp:431
bool m_surface_integral_built
Definition AMReX_MLNodeLaplacian.H:200
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
Real m_const_sigma
Definition AMReX_MLNodeLaplacian.H:186
void setSigma(int amrlev, const MultiFab &a_sigma)
Definition AMReX_MLNodeLaplacian.cpp:379
void setNormalizationThreshold(Real t) noexcept
Definition AMReX_MLNodeLaplacian.H:63
void buildStencil()
Definition AMReX_MLNodeLaplacian_sten.cpp:17
Vector< Vector< Real > > m_s0_norm0
Definition AMReX_MLNodeLaplacian.H:190
void compRHS(const Vector< MultiFab * > &rhs, const Vector< MultiFab * > &vel, const Vector< const MultiFab * > &rhnd, const Vector< MultiFab * > &rhcc)
Definition AMReX_MLNodeLaplacian_misc.cpp:974
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
int m_is_rz
Definition AMReX_MLNodeLaplacian.H:184
Vector< Vector< Array< std::unique_ptr< MultiFab >, 3 > > > m_sigma
Definition AMReX_MLNodeLaplacian.H:187
~MLNodeLaplacian() override=default
Vector< Real > getSolvabilityOffset(int amrlev, int mglev, MultiFab const &rhs) const override
Definition AMReX_MLNodeLaplacian.cpp:167
bool m_build_surface_integral
Definition AMReX_MLNodeLaplacian.H:201
MLNodeLaplacian(MLNodeLaplacian &&)=delete
bool m_integral_built
Definition AMReX_MLNodeLaplacian.H:197
void setMapped(bool flag) noexcept
Definition AMReX_MLNodeLaplacian.H:85
void setGaussSeidel(bool flag) noexcept
Definition AMReX_MLNodeLaplacian.H:82
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:125
BottomSolver getDefaultBottomSolver() const final
Definition AMReX_MLNodeLaplacian.H:91
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
bool m_use_harmonic_average
Definition AMReX_MLNodeLaplacian.H:206
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
Vector< std::unique_ptr< MultiFab > > m_nosigma_stencil
Definition AMReX_MLNodeLaplacian.H:189
void setRZCorrection(bool rz) noexcept
Definition AMReX_MLNodeLaplacian.H:61
void checkPoint(std::string const &file_name) const final
Definition AMReX_MLNodeLaplacian.cpp:941
void compDivergence(const Vector< MultiFab * > &rhs, const Vector< MultiFab * > &vel)
Definition AMReX_MLNodeLaplacian_misc.cpp:968
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
Vector< std::unique_ptr< MultiFab > > m_integral
Definition AMReX_MLNodeLaplacian.H:196
Vector< std::unique_ptr< MultiFab > > m_eb_vel_dot_n
Definition AMReX_MLNodeLaplacian.H:202
Vector< Vector< std::unique_ptr< MultiFab > > > m_stencil
Definition AMReX_MLNodeLaplacian.H:188
void resizeMultiGrid(int new_size) final
Definition AMReX_MLNodeLaplacian.cpp:121
void setCoarseningStrategy(CoarseningStrategy cs) noexcept
Definition AMReX_MLNodeLaplacian.H:87
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
Vector< std::unique_ptr< MultiFab > > m_surface_integral
Definition AMReX_MLNodeLaplacian.H:199
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:25
Definition AMReX_Amr.cpp:49
BottomSolver
Definition AMReX_MLLinOp.H:31
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
Definition AMReX_Array4.H:61
Definition AMReX_MLLinOp.H:36
Location
Definition AMReX_MLLinOp.H:90