Block-Structured AMR Software Framework
AMReX_EBFluxRegister.H
Go to the documentation of this file.
1 #ifndef AMREX_EBFLUXREGISTER_H_
2 #define AMREX_EBFLUXREGISTER_H_
3 #include <AMReX_Config.H>
4 
5 #include <AMReX_YAFluxRegister.H>
6 #include <AMReX_EBCellFlag.H>
7 
8 extern "C" {
11 }
12 
13 namespace amrex {
14 
57  : public YAFluxRegister
58 {
59 public:
60 
61  EBFluxRegister () = default;
62 
63  EBFluxRegister (const BoxArray& fba, const BoxArray& cba,
64  const DistributionMapping& fdm, const DistributionMapping& cdm,
65  const Geometry& fgeom, const Geometry& cgeom,
66  const IntVect& ref_ratio, int fine_lev, int nvar);
67 
68  void define (const BoxArray& fba, const BoxArray& cba,
69  const DistributionMapping& fdm, const DistributionMapping& cdm,
70  const Geometry& fgeom, const Geometry& cgeom,
71  const IntVect& ref_ratio, int fine_lev, int nvar);
72 
74  void CrseAdd (const MFIter& mfi,
75  const std::array<FArrayBox const*, AMREX_SPACEDIM>& flux,
76  const Real* dx, Real dt,
77  const FArrayBox& volfrac,
78  const std::array<FArrayBox const*, AMREX_SPACEDIM>& areafrac,
79  RunOn runon);
80  void CrseAdd (const MFIter& mfi,
81  const std::array<FArrayBox const*, AMREX_SPACEDIM>& flux,
82  const Real* dx, Real dt,
83  const FArrayBox& volfrac,
84  const std::array<FArrayBox const*, AMREX_SPACEDIM>& areafrac,
85  int srccomp, int destcomp, int numcomp, RunOn runon);
86 
88  void FineAdd (const MFIter& mfi,
89  const std::array<FArrayBox const*, AMREX_SPACEDIM>& flux,
90  const Real* dx, Real dt,
91  const FArrayBox& volfrac,
92  const std::array<FArrayBox const*, AMREX_SPACEDIM>& areafrac,
93  const FArrayBox& dm,
94  RunOn runon);
95  void FineAdd (const MFIter& mfi,
96  const std::array<FArrayBox const*, AMREX_SPACEDIM>& flux,
97  const Real* dx, Real dt,
98  const FArrayBox& volfrac,
99  const std::array<FArrayBox const*, AMREX_SPACEDIM>& areafrac,
100  const FArrayBox& dm,
101  int srccomp, int destcomp, int numcomp, RunOn runon);
103  void FineAdd (const MFIter& mfi,
104  const std::array<FArrayBox const*, AMREX_SPACEDIM>& flux,
105  const Real* dx, Real dt,
106  const FArrayBox& volfrac,
107  const std::array<FArrayBox const*, AMREX_SPACEDIM>& areafrac,
108  int srccomp, int destcomp, int numcomp, RunOn runon);
109 
110  void Reflux (MultiFab& crse_state, const amrex::MultiFab& crse_vfrac,
111  MultiFab& fine_state, const amrex::MultiFab& fine_vfrac);
112  void Reflux (MultiFab& crse_state, const amrex::MultiFab& crse_vfrac,
113  MultiFab& fine_state, const amrex::MultiFab& fine_vfrac,
114  int srccomp, int destcomp, int numcomp);
116  void Reflux (MultiFab& crse_state, const amrex::MultiFab& crse_vfrac,
117  int srccomp, int destcomp, int numcomp);
118 
119  FArrayBox* getCrseData (const MFIter& mfi) {
120  return &(m_crse_data[mfi]);
121  }
122 
123  const IArrayBox* getCrseFlag (const MFIter& mfi) const {
124  return &(m_crse_flag[mfi]);
125  }
126 
127 private:
128 
130 
131 public: // for cuda
132 
133  void defineExtra (const BoxArray& fba, const DistributionMapping& fdm);
134 };
135 
136 }
137 
138 #endif
void amrex_eb_disable_reredistribution()
amrex::Real amrex_eb_get_reredistribution_threshold()
A collection of Boxes stored in an Array.
Definition: AMReX_BoxArray.H:530
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
Definition: AMReX_EBFluxRegister.H:58
void Reflux(MultiFab &crse_state, const amrex::MultiFab &crse_vfrac, MultiFab &fine_state, const amrex::MultiFab &fine_vfrac)
Definition: AMReX_EBFluxRegister.cpp:289
void defineExtra(const BoxArray &fba, const DistributionMapping &fdm)
Definition: AMReX_EBFluxRegister.cpp:44
const IArrayBox * getCrseFlag(const MFIter &mfi) const
Definition: AMReX_EBFluxRegister.H:123
void FineAdd(const MFIter &mfi, const std::array< FAB const *, AMREX_SPACEDIM > &flux, const Real *dx, Real dt, RunOn runon) noexcept
Definition: AMReX_YAFluxRegister.H:401
void CrseAdd(const MFIter &mfi, const std::array< FAB const *, AMREX_SPACEDIM > &flux, const Real *dx, Real dt, RunOn runon) noexcept
Definition: AMReX_YAFluxRegister.H:343
iMultiFab m_cfp_inside_mask
Definition: AMReX_EBFluxRegister.H:129
FArrayBox * getCrseData(const MFIter &mfi)
Definition: AMReX_EBFluxRegister.H:119
void define(const BoxArray &fba, const BoxArray &cba, const DistributionMapping &fdm, const DistributionMapping &cdm, const Geometry &fgeom, const Geometry &cgeom, const IntVect &ref_ratio, int fine_lev, int nvar)
Definition: AMReX_EBFluxRegister.cpp:33
A Fortran Array of REALs.
Definition: AMReX_FArrayBox.H:229
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
A Fortran Array of ints.
Definition: AMReX_IArrayBox.H:47
Definition: AMReX_MFIter.H:57
A collection (stored as an array) of FArrayBox objects.
Definition: AMReX_MultiFab.H:38
Definition: AMReX_YAFluxRegister.H:28
iMultiFab m_crse_flag
Definition: AMReX_YAFluxRegister.H:96
MF m_crse_data
Definition: AMReX_YAFluxRegister.H:95
void FineAdd(const MFIter &mfi, const std::array< FAB const *, AMREX_SPACEDIM > &flux, const Real *dx, Real dt, RunOn runon) noexcept
Definition: AMReX_YAFluxRegister.H:401
void CrseAdd(const MFIter &mfi, const std::array< FAB const *, AMREX_SPACEDIM > &flux, const Real *dx, Real dt, RunOn runon) noexcept
Definition: AMReX_YAFluxRegister.H:343
Definition: AMReX_iMultiFab.H:32
Definition: AMReX_Amr.cpp:49
RunOn
Definition: AMReX_GpuControl.H:69