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