Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
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
14namespace amrex {
15
59{
60public:
61
63 EdgeFluxRegister () = default;
64
76 EdgeFluxRegister (const BoxArray& fba, const BoxArray& cba,
77 const DistributionMapping& fdm, const DistributionMapping& cdm,
78 const Geometry& fgeom, const Geometry& cgeom,
79 int nvar = 1);
80
82 void define (const BoxArray& fba, const BoxArray& cba,
83 const DistributionMapping& fdm, const DistributionMapping& cdm,
84 const Geometry& fgeom, const Geometry& cgeom,
85 int nvar = 1);
86
88 void reset ();
89
92
93#if (AMREX_SPACEDIM == 3)
94
102 void CrseAdd (MFIter const& mfi, const Array<FArrayBox const*,3>& E_crse, Real dt_crse);
103
111 void FineAdd (MFIter const& mfi, const Array<FArrayBox const*,3>& E_fine, Real dt_fine);
112
113#else /* 2D */
114
116 void CrseAdd (MFIter const& mfi, FArrayBox const& E_crse, Real dt_crse);
117
119 void FineAdd (MFIter const& mfi, FArrayBox const& E_fine, Real dt_fine);
120
121#endif
122
128 void Reflux (Array<MultiFab*,AMREX_SPACEDIM> const& B_crse) const;
129
130private:
131
132 Geometry m_fine_geom;
133 Geometry m_crse_geom;
134
135 IntVect m_ratio;
136 int m_ncomp;
137
138#if (AMREX_SPACEDIM == 3)
139
140 Array<MultiFab,AMREX_SPACEDIM> m_E_crse; // on original grids
141
142 // There are AMREX_SPACEDIM*2 faces. For each face, we need to store two
143 // component. For example, at the x-faces, we need to store Ey and Ez.
144 Array<Array<MultiFab,2>,AMREX_SPACEDIM*2> m_E_fine;
145
146 // Mask on the coarse level indicating overlap with m_E_fine
148
149#else
150
151 MultiFab m_E_crse;
153 iMultiFab m_fine_mask;
154
155#endif
156
157 LayoutData<int> m_has_cf; // Flag on the coarse level indicating c/f interface
158};
159
160}
161
162#endif
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
Definition AMReX_EdgeFluxRegister.H:59
void reset()
Reset accumulated state at the start of a coarse step.
Definition AMReX_EdgeFluxRegister.cpp:118
void Reflux(Array< MultiFab *, 3 > const &B_crse) const
Apply the accumulated edge fluxes to coarse face-centered data.
Definition AMReX_EdgeFluxRegister.cpp:311
void FineAdd(MFIter const &mfi, const Array< FArrayBox const *, 3 > &E_fine, Real dt_fine)
Accumulate fine edge-centered fluxes for tile mfi.
Definition AMReX_EdgeFluxRegister.cpp:242
EdgeFluxRegister()=default
Construct an empty register; call define() before using it.
void define(const BoxArray &fba, const BoxArray &cba, const DistributionMapping &fdm, const DistributionMapping &cdm, const Geometry &fgeom, const Geometry &cgeom, int nvar=1)
Define the register using fba, cba, fdm, cdm, fgeom, cgeom, and component count nvar.
Definition AMReX_EdgeFluxRegister.cpp:19
void CrseAdd(MFIter const &mfi, const Array< FArrayBox const *, 3 > &E_crse, Real dt_crse)
Accumulate coarse edge-centered fluxes for tile mfi.
Definition AMReX_EdgeFluxRegister.cpp:212
EdgeFluxRegister & plus(EdgeFluxRegister const &rhs)
Accumulate another register with identical layout into this one.
Definition AMReX_EdgeFluxRegister.cpp:186
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
a one-thingy-per-box distributed object
Definition AMReX_LayoutData.H:13
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
A Collection of IArrayBoxes.
Definition AMReX_iMultiFab.H:34
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
std::array< T, N > Array
Definition AMReX_Array.H:26
Definition AMReX_Amr.cpp:49