Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_MLNodeLinOp.H
Go to the documentation of this file.
1#ifndef AMREX_ML_NODE_LINOP_H_H
2#define AMREX_ML_NODE_LINOP_H_H
3#include <AMReX_Config.H>
4
5#include <AMReX_MLLinOp.H>
6#include <AMReX_iMultiFab.H>
7
8#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
10#endif
11
12namespace amrex {
13
15 : public MLLinOp
16{
17public:
18
19 enum struct CoarseningStrategy : int { Sigma, RAP };
20
21 MLNodeLinOp ();
22 ~MLNodeLinOp () override = default;
23
24 MLNodeLinOp (const MLNodeLinOp&) = delete;
28
29 void define (const Vector<Geometry>& a_geom,
30 const Vector<BoxArray>& a_grids,
31 const Vector<DistributionMapping>& a_dmap,
32 const LPInfo& a_info = LPInfo(),
33 const Vector<FabFactory<FArrayBox> const*>& a_factory = {},
34 int a_eb_limit_coarsening = -1);
35
36 void setSmoothNumSweeps (int nsweeps) noexcept {
37 m_smooth_num_sweeps = nsweeps;
38 }
39
40 void setLevelBC (int /*amrlev*/, const MultiFab* /*levelbcdata*/,
41 const MultiFab* = nullptr, const MultiFab* = nullptr,
42 const MultiFab* = nullptr) final {}
43
44 void apply (int amrlev, int mglev, MultiFab& out, MultiFab& in, BCMode bc_mode,
45 StateMode s_mode, const MLMGBndry* bndry=nullptr) const final;
46
47 void smooth (int amrlev, int mglev, MultiFab& sol, const MultiFab& rhs,
48 bool skip_fillboundary=false) const override;
49
50 void solutionResidual (int amrlev, MultiFab& resid, MultiFab& x, const MultiFab& b,
51 const MultiFab* crse_bcdata=nullptr) override;
52 void correctionResidual (int amrlev, int mglev, MultiFab& resid, MultiFab& x, const MultiFab& b,
53 BCMode bc_mode, const MultiFab* crse_bcdata=nullptr) override;
54
55 Vector<Real> getSolvabilityOffset (int amrlev, int mglev,
56 MultiFab const& rhs) const override;
57 void fixSolvabilityByOffset (int amrlev, int mglev, MultiFab& rhs,
58 Vector<Real> const& offset) const override;
59
60 void prepareForSolve () override;
61
62 void preparePrecond () override;
63
64 void setDirichletNodesToZero (int amrlev, int mglev, MultiFab& mf) const override;
65
66 bool isSingular (int amrlev) const override
67 { return (amrlev == 0) ? m_is_bottom_singular : false; }
68 bool isBottomSingular () const override { return m_is_bottom_singular; }
69
70 Real xdoty (int amrlev, int mglev, const MultiFab& x, const MultiFab& y, bool local) const final;
71
73 Vector<MultiFab const*> const& y) const final;
74
75 Real norm2Precond (Vector<MultiFab const*> const& x) const final;
76
77 virtual void applyBC (int amrlev, int mglev, MultiFab& phi, BCMode bc_mode, StateMode state_mode,
78 bool skip_fillboundary=false) const;
79
80 virtual void Fapply (int amrlev, int mglev, MultiFab& out, const MultiFab& in) const = 0;
81 virtual void Fsmooth (int amrlev, int mglev, MultiFab& sol, const MultiFab& rhs) const = 0;
82
83 void nodalSync (int amrlev, int mglev, MultiFab& mf) const;
84
85 static std::unique_ptr<iMultiFab> makeOwnerMask (const BoxArray& ba,
86 const DistributionMapping& dm,
87 const Geometry& geom);
88
89 void buildMasks ();
90
91 // omask is either 0 or 1. 1 means the node is an unknown. 0 means it's known.
92 void setOversetMask (int amrlev, const iMultiFab& a_dmask);
93
94 virtual void fixUpResidualMask (int /*amrlev*/, iMultiFab& /*resmsk*/) { }
95
96 Real normInf (int amrlev, MultiFab const& mf, bool local) const override;
97
98 void avgDownResAmr (int, MultiFab&, MultiFab const&) const final { }
99
100 void interpolationAmr (int famrlev, MultiFab& fine, const MultiFab& crse,
101 IntVect const& nghost) const override;
102
103 void averageDownAndSync (Vector<MultiFab>& sol) const override;
104
105 void interpAssign (int amrlev, int fmglev, MultiFab& fine, MultiFab& crse) const override;
106
107#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
108 [[nodiscard]] std::unique_ptr<HypreNodeLap> makeHypreNodeLap(
109 int bottom_verbose,
110 const std::string& options_namespace) const override;
111
112 virtual void fillIJMatrix (MFIter const& /*mfi*/,
114 Array4<int const> const& /*lid*/,
115 HypreNodeLap::Int* /*ncols*/,
116 HypreNodeLap::Int* /*cols*/,
117 Real* /*mat*/) const
118 {
119 amrex::Abort("MLNodeLinOp::fillIJMatrix: how did we get here?");
120 }
121
122 virtual void fillRHS (MFIter const& /*mfi*/,
123 Array4<int const> const& /*lid*/,
124 Real* /*rhs*/,
125 Array4<Real const> const& /*bfab*/) const
126 {
127 amrex::Abort("MLNodeLinOp:fillRHS: how did we get here?");
128 }
129#endif
130
131protected:
132
133 void resizeMultiGrid (int new_size) override;
134
135 std::unique_ptr<iMultiFab> m_owner_mask_top; // ownership of nodes
136 std::unique_ptr<iMultiFab> m_owner_mask_bottom;
138 Vector<std::unique_ptr<iMultiFab> > m_cc_fine_mask; // cell-centered mask for cells covered by fine
139 Vector<std::unique_ptr<iMultiFab> > m_nd_fine_mask; // nodal mask: 0: this level node, 1: c/f boundary, 2: fine node
140 Vector<std::unique_ptr<LayoutData<int> > > m_has_fine_bndry; // does this fab contain c/f boundary?
144
146
147#ifdef AMREX_USE_EB
149#else
151#endif
152
153 bool m_masks_built = false;
155#ifdef AMREX_USE_GPU
157#else
158 int m_smooth_num_sweeps = 2;
159#endif
160 mutable bool m_in_solution_mode = true;
161
162private:
164};
165
166}
167
168#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:550
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:41
Definition AMReX_FabFactory.H:50
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:73
HYPRE_Int Int
Definition AMReX_HypreNodeLap.H:36
Definition AMReX_MFIter.H:57
Definition AMReX_MLLinOp.H:98
LinOpEnumType::StateMode StateMode
Definition AMReX_MLLinOp.H:113
LinOpEnumType::BCMode BCMode
Definition AMReX_MLLinOp.H:112
Definition AMReX_MLMGBndry.H:12
Definition AMReX_MLNodeLinOp.H:16
MLNodeLinOp(MLNodeLinOp &&)=delete
Vector< std::unique_ptr< iMultiFab > > m_nd_fine_mask
Definition AMReX_MLNodeLinOp.H:139
void nodalSync(int amrlev, int mglev, MultiFab &mf) const
Definition AMReX_MLNodeLinOp.cpp:117
MultiFab m_coarse_dot_mask
Definition AMReX_MLNodeLinOp.H:142
Real dotProductPrecond(Vector< MultiFab const * > const &x, Vector< MultiFab const * > const &y) const final
Definition AMReX_MLNodeLinOp.cpp:195
void correctionResidual(int amrlev, int mglev, MultiFab &resid, MultiFab &x, const MultiFab &b, BCMode bc_mode, const MultiFab *crse_bcdata=nullptr) override
Definition AMReX_MLNodeLinOp.cpp:152
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={}, int a_eb_limit_coarsening=-1)
Definition AMReX_MLNodeLinOp.cpp:23
bool isBottomSingular() const override
Is the bottom of MG singular?
Definition AMReX_MLNodeLinOp.H:68
bool m_in_solution_mode
Definition AMReX_MLNodeLinOp.H:160
virtual void fixUpResidualMask(int, iMultiFab &)
Definition AMReX_MLNodeLinOp.H:94
bool m_masks_built
Definition AMReX_MLNodeLinOp.H:153
CoarseningStrategy
Definition AMReX_MLNodeLinOp.H:19
std::unique_ptr< iMultiFab > m_owner_mask_bottom
Definition AMReX_MLNodeLinOp.H:136
virtual void Fsmooth(int amrlev, int mglev, MultiFab &sol, const MultiFab &rhs) const =0
MultiFab m_bottom_dot_mask
Definition AMReX_MLNodeLinOp.H:141
void interpAssign(int amrlev, int fmglev, MultiFab &fine, MultiFab &crse) const override
Definition AMReX_MLNodeLinOp.cpp:626
bool m_is_bottom_singular
Definition AMReX_MLNodeLinOp.H:163
MLNodeLinOp & operator=(const MLNodeLinOp &)=delete
void interpolationAmr(int famrlev, MultiFab &fine, const MultiFab &crse, IntVect const &nghost) const override
Definition AMReX_MLNodeLinOp.cpp:573
Real normInf(int amrlev, MultiFab const &mf, bool local) const override
Definition AMReX_MLNodeLinOp.cpp:561
MLNodeLinOp(const MLNodeLinOp &)=delete
void avgDownResAmr(int, MultiFab &, MultiFab const &) const final
Definition AMReX_MLNodeLinOp.H:98
Vector< std::unique_ptr< LayoutData< int > > > m_has_fine_bndry
Definition AMReX_MLNodeLinOp.H:140
bool m_overset_dirichlet_mask
Definition AMReX_MLNodeLinOp.H:154
Vector< Real > getSolvabilityOffset(int amrlev, int mglev, MultiFab const &rhs) const override
Definition AMReX_MLNodeLinOp.cpp:222
void setLevelBC(int, const MultiFab *, const MultiFab *=nullptr, const MultiFab *=nullptr, const MultiFab *=nullptr) final
Definition AMReX_MLNodeLinOp.H:40
void preparePrecond() override
Definition AMReX_MLNodeLinOp.cpp:409
Real norm2Precond(Vector< MultiFab const * > const &x) const final
Definition AMReX_MLNodeLinOp.cpp:209
void fixSolvabilityByOffset(int amrlev, int mglev, MultiFab &rhs, Vector< Real > const &offset) const override
Definition AMReX_MLNodeLinOp.cpp:247
Vector< std::unique_ptr< iMultiFab > > m_cc_fine_mask
Definition AMReX_MLNodeLinOp.H:138
Real xdoty(int amrlev, int mglev, const MultiFab &x, const MultiFab &y, bool local) const final
Definition AMReX_MLNodeLinOp.cpp:179
void solutionResidual(int amrlev, MultiFab &resid, MultiFab &x, const MultiFab &b, const MultiFab *crse_bcdata=nullptr) override
Definition AMReX_MLNodeLinOp.cpp:123
void averageDownAndSync(Vector< MultiFab > &sol) const override
Definition AMReX_MLNodeLinOp.cpp:605
void prepareForSolve() override
Definition AMReX_MLNodeLinOp.cpp:100
void buildMasks()
Definition AMReX_MLNodeLinOp.cpp:284
std::unique_ptr< iMultiFab > m_owner_mask_top
Definition AMReX_MLNodeLinOp.H:135
void smooth(int amrlev, int mglev, MultiFab &sol, const MultiFab &rhs, bool skip_fillboundary=false) const override
Definition AMReX_MLNodeLinOp.cpp:169
Vector< Vector< std::unique_ptr< iMultiFab > > > m_dirichlet_mask
Definition AMReX_MLNodeLinOp.H:137
Vector< std::unique_ptr< iMultiFab > > m_norm_fine_mask
Definition AMReX_MLNodeLinOp.H:145
MLNodeLinOp()
Definition AMReX_MLNodeLinOp.cpp:17
int m_smooth_num_sweeps
Definition AMReX_MLNodeLinOp.H:156
void apply(int amrlev, int mglev, MultiFab &out, MultiFab &in, BCMode bc_mode, StateMode s_mode, const MLMGBndry *bndry=nullptr) const final
Definition AMReX_MLNodeLinOp.cpp:161
void resizeMultiGrid(int new_size) override
Definition AMReX_MLNodeLinOp.cpp:528
void setSmoothNumSweeps(int nsweeps) noexcept
Definition AMReX_MLNodeLinOp.H:36
void setDirichletNodesToZero(int amrlev, int mglev, MultiFab &mf) const override
Definition AMReX_MLNodeLinOp.cpp:463
static std::unique_ptr< iMultiFab > makeOwnerMask(const BoxArray &ba, const DistributionMapping &dm, const Geometry &geom)
Definition AMReX_MLNodeLinOp.cpp:108
Vector< MultiFab > m_precond_weight_mask
Definition AMReX_MLNodeLinOp.H:143
bool isSingular(int amrlev) const override
Is it singular on given AMR level?
Definition AMReX_MLNodeLinOp.H:66
~MLNodeLinOp() override=default
CoarseningStrategy m_coarsening_strategy
Definition AMReX_MLNodeLinOp.H:150
virtual void applyBC(int amrlev, int mglev, MultiFab &phi, BCMode bc_mode, StateMode state_mode, bool skip_fillboundary=false) const
Definition AMReX_MLNodeLinOp.cpp:498
void setOversetMask(int amrlev, const iMultiFab &a_dmask)
Definition AMReX_MLNodeLinOp.cpp:480
virtual void Fapply(int amrlev, int mglev, MultiFab &out, const MultiFab &in) const =0
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
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:225
Definition AMReX_Array4.H:61
Definition AMReX_MLLinOp.H:35