Line data Source code
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)
9 : #include <AMReX_HypreNodeLap.H>
10 : #endif
11 :
12 : namespace amrex {
13 :
14 : class MLNodeLinOp
15 : : public MLLinOp
16 : {
17 : public:
18 :
19 : enum struct CoarseningStrategy : int { Sigma, RAP };
20 :
21 : MLNodeLinOp ();
22 126 : ~MLNodeLinOp () override = default;
23 :
24 : MLNodeLinOp (const MLNodeLinOp&) = delete;
25 : MLNodeLinOp (MLNodeLinOp&&) = delete;
26 : MLNodeLinOp& operator= (const MLNodeLinOp&) = delete;
27 : MLNodeLinOp& operator= (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 0 : void setLevelBC (int /*amrlev*/, const MultiFab* /*levelbcdata*/,
41 : const MultiFab* = nullptr, const MultiFab* = nullptr,
42 0 : 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 prepareForGMRES () 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 :
72 : virtual void applyBC (int amrlev, int mglev, MultiFab& phi, BCMode bc_mode, StateMode state_mode,
73 : bool skip_fillboundary=false) const;
74 :
75 : virtual void Fapply (int amrlev, int mglev, MultiFab& out, const MultiFab& in) const = 0;
76 : virtual void Fsmooth (int amrlev, int mglev, MultiFab& sol, const MultiFab& rhs) const = 0;
77 :
78 : void nodalSync (int amrlev, int mglev, MultiFab& mf) const;
79 :
80 : static std::unique_ptr<iMultiFab> makeOwnerMask (const BoxArray& ba,
81 : const DistributionMapping& dm,
82 : const Geometry& geom);
83 :
84 : void buildMasks ();
85 :
86 : // omask is either 0 or 1. 1 means the node is an unknown. 0 means it's known.
87 : void setOversetMask (int amrlev, const iMultiFab& a_dmask);
88 :
89 : virtual void fixUpResidualMask (int /*amrlev*/, iMultiFab& /*resmsk*/) { }
90 :
91 : Real normInf (int amrlev, MultiFab const& mf, bool local) const override;
92 :
93 260 : void avgDownResAmr (int, MultiFab&, MultiFab const&) const final { }
94 :
95 : void interpolationAmr (int famrlev, MultiFab& fine, const MultiFab& crse,
96 : IntVect const& nghost) const override;
97 :
98 : void averageDownAndSync (Vector<MultiFab>& sol) const override;
99 :
100 : void interpAssign (int amrlev, int fmglev, MultiFab& fine, MultiFab& crse) const override;
101 :
102 : #if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
103 : [[nodiscard]] std::unique_ptr<HypreNodeLap> makeHypreNodeLap(
104 : int bottom_verbose,
105 : const std::string& options_namespace) const override;
106 :
107 : virtual void fillIJMatrix (MFIter const& /*mfi*/,
108 : Array4<HypreNodeLap::AtomicInt const> const& /*gid*/,
109 : Array4<int const> const& /*lid*/,
110 : HypreNodeLap::Int* /*ncols*/,
111 : HypreNodeLap::Int* /*cols*/,
112 : Real* /*mat*/) const
113 : {
114 : amrex::Abort("MLNodeLinOp::fillIJMatrix: how did we get here?");
115 : }
116 :
117 : virtual void fillRHS (MFIter const& /*mfi*/,
118 : Array4<int const> const& /*lid*/,
119 : Real* /*rhs*/,
120 : Array4<Real const> const& /*bfab*/) const
121 : {
122 : amrex::Abort("MLNodeLinOp:fillRHS: how did we get here?");
123 : }
124 : #endif
125 :
126 : protected:
127 :
128 : void resizeMultiGrid (int new_size) override;
129 :
130 : std::unique_ptr<iMultiFab> m_owner_mask_top; // ownership of nodes
131 : std::unique_ptr<iMultiFab> m_owner_mask_bottom;
132 : Vector<Vector<std::unique_ptr<iMultiFab> > > m_dirichlet_mask; // dirichlet?
133 : Vector<std::unique_ptr<iMultiFab> > m_cc_fine_mask; // cell-centered mask for cells covered by fine
134 : Vector<std::unique_ptr<iMultiFab> > m_nd_fine_mask; // nodal mask: 0: this level node, 1: c/f boundary, 2: fine node
135 : Vector<std::unique_ptr<LayoutData<int> > > m_has_fine_bndry; // does this fab contain c/f boundary?
136 : MultiFab m_bottom_dot_mask;
137 : MultiFab m_coarse_dot_mask;
138 :
139 : Vector<std::unique_ptr<iMultiFab> > m_norm_fine_mask;
140 :
141 : #ifdef AMREX_USE_EB
142 : CoarseningStrategy m_coarsening_strategy = CoarseningStrategy::RAP;
143 : #else
144 : CoarseningStrategy m_coarsening_strategy = CoarseningStrategy::Sigma;
145 : #endif
146 :
147 : bool m_masks_built = false;
148 : bool m_overset_dirichlet_mask = false;
149 : #ifdef AMREX_USE_GPU
150 : int m_smooth_num_sweeps = 4;
151 : #else
152 : int m_smooth_num_sweeps = 2;
153 : #endif
154 : mutable bool m_in_solution_mode = true;
155 :
156 : private:
157 : bool m_is_bottom_singular = false;
158 : };
159 :
160 : }
161 :
162 : #endif
|