Block-Structured AMR Software Framework
AMReX_EdgeFluxRegister.H
Go to the documentation of this file.
1 #ifndef AMREX_EDGE_FLUX_REGISTER_H_
2 #define AMREX_EDGE_FLUX_REGISTER_H_
3 #include <AMReX_Config.H>
4 
5 #include <AMReX_iMultiFab.H>
6 #include <AMReX_LayoutData.H>
7 #include <AMReX_MultiFab.H>
8 
9 namespace amrex {
10 
54 {
55 public:
56 
57  EdgeFluxRegister () = default;
58 
59  EdgeFluxRegister (const BoxArray& fba, const BoxArray& cba,
60  const DistributionMapping& fdm, const DistributionMapping& cdm,
61  const Geometry& fgeom, const Geometry& cgeom,
62  int nvar = 1);
63 
64  void define (const BoxArray& fba, const BoxArray& cba,
65  const DistributionMapping& fdm, const DistributionMapping& cdm,
66  const Geometry& fgeom, const Geometry& cgeom,
67  int nvar = 1);
68 
69  void reset ();
70 
71 #if (AMREX_SPACEDIM == 3)
72 
73  void CrseAdd (MFIter const& mfi, const Array<FArrayBox const*,3>& E_crse, Real dt_crse);
74  void FineAdd (MFIter const& mfi, const Array<FArrayBox const*,3>& E_fine, Real dt_fine);
75 
76 #else /* 2D */
77 
78  void CrseAdd (MFIter const& mfi, FArrayBox const& E_crse, Real dt_crse);
79  void FineAdd (MFIter const& mfi, FArrayBox const& E_fine, Real dt_fine);
80 
81 #endif
82 
83  void Reflux (Array<MultiFab*,AMREX_SPACEDIM> const& B_crse) const;
84 
85 private:
86 
89 
91  int m_ncomp;
92 
93 #if (AMREX_SPACEDIM == 3)
94 
95  Array<MultiFab,AMREX_SPACEDIM> m_E_crse; // on original grids
96 
97  // There are AMREX_SPACEDIM*2 faces. For each face, we need to store two
98  // component. For example, at the x-faces, we need to store Ey and Ez.
99  Array<Array<MultiFab,2>,AMREX_SPACEDIM*2> m_E_fine;
100 
101  // Mask on the coarse level indicating overlap with m_E_fine
103 
104 #else
105 
109 
110 #endif
111 
112  LayoutData<int> m_has_cf; // Flag on the coarse level indicating c/f interface
113 };
114 
115 }
116 
117 #endif
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_EdgeFluxRegister.H:54
LayoutData< int > m_has_cf
Definition: AMReX_EdgeFluxRegister.H:112
void reset()
Definition: AMReX_EdgeFluxRegister.cpp:115
void Reflux(Array< MultiFab *, AMREX_SPACEDIM > const &B_crse) const
Definition: AMReX_EdgeFluxRegister.cpp:283
MultiFab m_E_crse
Definition: AMReX_EdgeFluxRegister.H:106
iMultiFab m_fine_mask
Definition: AMReX_EdgeFluxRegister.H:108
Geometry m_crse_geom
Definition: AMReX_EdgeFluxRegister.H:88
IntVect m_ratio
Definition: AMReX_EdgeFluxRegister.H:90
int m_ncomp
Definition: AMReX_EdgeFluxRegister.H:91
void CrseAdd(MFIter const &mfi, FArrayBox const &E_crse, Real dt_crse)
Definition: AMReX_EdgeFluxRegister.cpp:249
void define(const BoxArray &fba, const BoxArray &cba, const DistributionMapping &fdm, const DistributionMapping &cdm, const Geometry &fgeom, const Geometry &cgeom, int nvar=1)
Definition: AMReX_EdgeFluxRegister.cpp:18
Array< MultiFab, AMREX_SPACEDIM *2 > m_E_fine
Definition: AMReX_EdgeFluxRegister.H:107
Geometry m_fine_geom
Definition: AMReX_EdgeFluxRegister.H:87
void FineAdd(MFIter const &mfi, FArrayBox const &E_fine, Real dt_fine)
Definition: AMReX_EdgeFluxRegister.cpp:263
A Fortran Array of REALs.
Definition: AMReX_FArrayBox.H:229
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
Definition: AMReX_MFIter.H:57
A collection (stored as an array) of FArrayBox objects.
Definition: AMReX_MultiFab.H:38
Definition: AMReX_iMultiFab.H:32
Definition: AMReX_Amr.cpp:49
std::array< T, N > Array
Definition: AMReX_Array.H:23