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