Block-Structured AMR Software Framework
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Loading...
Searching...
No Matches
AMReX_MLNodeABecLaplacian.H
Go to the documentation of this file.
1#ifndef AMREX_MLNODEABECLAPLACIAN_H_
2#define AMREX_MLNODEABECLAPLACIAN_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_MLNodeLinOp.H>
6
7namespace amrex {
8
9// (alpha * a - beta * (del dot b grad)) phi = rhs
10// a, phi and rhs are nodal. b is cell-centered.
11
13 : public MLNodeLinOp
14{
15public:
16
17 MLNodeABecLaplacian () = default;
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 ~MLNodeABecLaplacian () override = default;
24
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
36 std::string name () const override { return std::string("MLNodeABecLaplacian"); }
37
38 void setScalars (Real a, Real b) {
39 m_a_scalar = a;
40 m_b_scalar = b;
41 }
42
43 void setACoeffs (int amrlev, Real a_acoef);
44 void setACoeffs (int amrlev, const MultiFab& a_acoef);
45
46 void setBCoeffs (int amrlev, Real a_bcoef);
47 void setBCoeffs (int amrlev, const MultiFab& a_bcoef);
48
49 void Fapply (int amrlev, int mglev, MultiFab& out, const MultiFab& in) const final;
50 void Fsmooth (int amrlev, int mglev, MultiFab& sol, const MultiFab& rhs) const final;
51
52 void fixUpResidualMask (int amrlev, iMultiFab& resmsk) final;
53
54 bool isSingular (int /*amrlev*/) const final { return false; }
55 bool isBottomSingular () const final { return false; }
56
57 void restriction (int amrlev, int cmglev, MultiFab& crse, MultiFab& fine) const final;
58 void interpolation (int amrlev, int fmglev, MultiFab& fine, const MultiFab& crse) const final;
59 void averageDownSolutionRHS (int camrlev, MultiFab& crse_sol, MultiFab& crse_rhs,
60 const MultiFab& fine_sol, const MultiFab& fine_rhs) final;
61
62 void reflux (int crse_amrlev,
63 MultiFab& res, const MultiFab& crse_sol, const MultiFab& crse_rhs,
64 MultiFab& fine_res, MultiFab& fine_sol, const MultiFab& fine_rhs) const final;
65
66 void prepareForSolve () final;
67
68 [[nodiscard]] bool needsUpdate () const final { return m_needs_update; }
69
70 void update () final;
71
72 void averageDownCoeffs ();
74 void averageDownCoeffsSameAmrLevel (int amrlev);
75
76private:
77
78 bool m_needs_update = true;
79
80 Real m_a_scalar = std::numeric_limits<Real>::quiet_NaN();
81 Real m_b_scalar = std::numeric_limits<Real>::quiet_NaN();
84};
85
86}
87
88#endif
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
Definition AMReX_FabFactory.H:50
Definition AMReX_MLNodeABecLaplacian.H:14
void setACoeffs(int amrlev, Real a_acoef)
Definition AMReX_MLNodeABecLaplacian.cpp:55
void averageDownCoeffsSameAmrLevel(int amrlev)
Definition AMReX_MLNodeABecLaplacian.cpp:336
void Fapply(int amrlev, int mglev, MultiFab &out, const MultiFab &in) const final
Definition AMReX_MLNodeABecLaplacian.cpp:85
bool m_needs_update
Definition AMReX_MLNodeABecLaplacian.H:78
Vector< Vector< MultiFab > > m_a_coeffs
Definition AMReX_MLNodeABecLaplacian.H:82
void fixUpResidualMask(int amrlev, iMultiFab &resmsk) final
Definition AMReX_MLNodeABecLaplacian.cpp:271
std::string name() const override
Definition AMReX_MLNodeABecLaplacian.H:36
Vector< Vector< MultiFab > > m_b_coeffs
Definition AMReX_MLNodeABecLaplacian.H:83
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={})
Definition AMReX_MLNodeABecLaplacian.cpp:17
bool isSingular(int) const final
Is it singular on given AMR level?
Definition AMReX_MLNodeABecLaplacian.H:54
~MLNodeABecLaplacian() override=default
void setBCoeffs(int amrlev, Real a_bcoef)
Definition AMReX_MLNodeABecLaplacian.cpp:70
void update() final
Update for reuse.
Definition AMReX_MLNodeABecLaplacian.cpp:263
void Fsmooth(int amrlev, int mglev, MultiFab &sol, const MultiFab &rhs) const final
Definition AMReX_MLNodeABecLaplacian.cpp:116
bool isBottomSingular() const final
Is the bottom of MG singular?
Definition AMReX_MLNodeABecLaplacian.H:55
MLNodeABecLaplacian & operator=(const MLNodeABecLaplacian &)=delete
void averageDownCoeffsToCoarseAmrLevel(int flev)
Definition AMReX_MLNodeABecLaplacian.cpp:321
MLNodeABecLaplacian(MLNodeABecLaplacian &&)=delete
void setScalars(Real a, Real b)
Definition AMReX_MLNodeABecLaplacian.H:38
void prepareForSolve() final
Definition AMReX_MLNodeABecLaplacian.cpp:249
Real m_b_scalar
Definition AMReX_MLNodeABecLaplacian.H:81
MLNodeABecLaplacian(const MLNodeABecLaplacian &)=delete
void averageDownSolutionRHS(int camrlev, MultiFab &crse_sol, MultiFab &crse_rhs, const MultiFab &fine_sol, const MultiFab &fine_rhs) final
Definition AMReX_MLNodeABecLaplacian.cpp:232
Real m_a_scalar
Definition AMReX_MLNodeABecLaplacian.H:80
bool needsUpdate() const final
Does it need update if it's reused?
Definition AMReX_MLNodeABecLaplacian.H:68
void averageDownCoeffs()
Definition AMReX_MLNodeABecLaplacian.cpp:287
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_MLNodeABecLaplacian.cpp:240
void interpolation(int amrlev, int fmglev, MultiFab &fine, const MultiFab &crse) const final
Definition AMReX_MLNodeABecLaplacian.cpp:204
void restriction(int amrlev, int cmglev, MultiFab &crse, MultiFab &fine) const final
Definition AMReX_MLNodeABecLaplacian.cpp:173
Definition AMReX_MLNodeLinOp.H:16
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
Definition AMReX_MLLinOp.H:35