REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_FillPatcher.H
Go to the documentation of this file.
1#ifndef REMORA_FILLPATCHER_H_
2#define REMORA_FILLPATCHER_H_
3
4#include <AMReX_FillPatchUtil.H>
5#include <AMReX_Interp_C.H>
6#include <AMReX_MFInterp_C.H>
7
9{
10public:
11
12 /** \brief Fill valid and ghost data with the "state data" at the given time */
13 REMORAFillPatcher (amrex::BoxArray const& fba, amrex::DistributionMapping const& fdm,
14 amrex::Geometry const& fgeom,
15 amrex::BoxArray const& cba, amrex::DistributionMapping const& cdm,
16 amrex::Geometry const& cgeom,
17 int nghost, int nghost_set, int ncomp, amrex::InterpBase* interp);
18
19 /** \brief Redefine the coarse and fine patch MultiFabs. */
20 void Define (amrex::BoxArray const& fba, amrex::DistributionMapping const& fdm,
21 amrex::Geometry const& fgeom,
22 amrex::BoxArray const& cba, amrex::DistributionMapping const& cdm,
23 amrex::Geometry const& cgeom,
24 int nghost, int nghost_set, int ncomp,
25 amrex::InterpBase* interp);
26
27 /** \brief Generate masking array */
28 void BuildMask (amrex::BoxArray const& fba, int nghost, int nghost_set);
29
30 /** \brief Register the coarse data to be used by the REMORAFillPatcher */
31 void RegisterCoarseData (amrex::Vector<amrex::MultiFab const*> const& crse_data,
32 amrex::Vector<amrex::Real> const& crse_time);
33
34 /** \brief Interpolate to cell faces */
35 void InterpFace (amrex::MultiFab& fine,
36 amrex::MultiFab const& crse,
37 int mask_val);
38
39 /** \brief Interpolate to cell centers */
40 void InterpCell (amrex::MultiFab& fine,
41 amrex::MultiFab const& crse,
42 amrex::Vector<amrex::BCRec> const& bcr,
43 int mask_val);
44
45 int GetSetMaskVal () { return m_set_mask; }
46
47 int GetRelaxMaskVal () { return m_relax_mask; }
48
49 amrex::iMultiFab* GetMask () { return m_cf_mask.get(); }
50
51 /** \brief Fill fine data in the set region */
52 template <typename BC>
53 void FillSet (amrex::MultiFab& mf, amrex::Real time,
54 BC& cbc, amrex::Vector<amrex::BCRec> const& bcs);
55
56 /** \brief Fill fine data in the relax region */
57 template <typename BC>
58 void FillRelax (amrex::MultiFab& mf, amrex::Real time,
59 BC& cbc, amrex::Vector<amrex::BCRec> const& bcs);
60 /** \brief Fill fine data in the relax region */
61 template <typename BC>
62 void Fill (amrex::MultiFab& mf, amrex::Real time,
63 BC& cbc, amrex::Vector<amrex::BCRec> const& bcs, int mask_val);
64
65private:
66
67 amrex::BoxArray m_fba;
68 amrex::BoxArray m_cba;
69 amrex::DistributionMapping m_fdm;
70 amrex::DistributionMapping m_cdm;
71 amrex::Geometry m_fgeom;
72 amrex::Geometry m_cgeom;
76 amrex::InterpBase* m_interp;
77 amrex::IntVect m_ratio;
78 std::unique_ptr<amrex::MultiFab> m_cf_crse_data_old;
79 std::unique_ptr<amrex::MultiFab> m_cf_crse_data_new;
80 std::unique_ptr<amrex::iMultiFab> m_cf_mask;
81 amrex::Vector<amrex::Real> m_crse_times;
82 amrex::Real m_dt_crse;
83 int m_set_mask{2};
85};
86
87/**
88 * @param[out] mf MultiFab to be filled
89 * @param[in] time Time at which to fill data
90 * @param[in] cbc Coarse boundary condition
91 * @param[in] bcs Vector of boundary conditions
92 */
93template <typename BC>
94void
95REMORAFillPatcher::FillSet (amrex::MultiFab& mf, amrex::Real time,
96 BC& cbc, amrex::Vector<amrex::BCRec> const& bcs)
97{
98 Fill(mf,time,cbc,bcs,m_set_mask);
99}
100
101/**
102 * @param[out] mf MultiFab to be filled
103 * @param[in] time Time at which to fill data
104 * @param[in] cbc Coarse boundary condition
105 * @param[in] bcs Vector of boundary conditions
106 */
107template <typename BC>
108void
109REMORAFillPatcher::FillRelax (amrex::MultiFab& mf, amrex::Real time,
110 BC& cbc, amrex::Vector<amrex::BCRec> const& bcs)
111{
112 Fill(mf,time,cbc,bcs,m_relax_mask);
113}
114
115/**
116 * @param[out] mf MultiFab to be filled
117 * @param[in] time Time at which to fill data
118 * @param[in] cbc Coarse boundary condition
119 * @param[in] bcs Vector of boundary conditions
120 * @param[in] mask_val Value to assign mask array
121 */
122template <typename BC>
123void
124REMORAFillPatcher::Fill (amrex::MultiFab& mf, amrex::Real time,
125 BC& cbc, amrex::Vector<amrex::BCRec> const& bcs, int mask_val)
126{
127 constexpr amrex::Real eps = std::numeric_limits<float>::epsilon();
128
129 AMREX_ALWAYS_ASSERT((time >= m_crse_times[0]-eps) && (time <= m_crse_times[1]+eps));
130
131 // Time interpolation factors
132 amrex::Real fac_new = (time - m_crse_times[0]) / m_dt_crse;
133 amrex::Real fac_old = 1.0 - fac_new;
134
135 // Boundary condition operator
136 cbc(*(m_cf_crse_data_old), 0, m_ncomp, amrex::IntVect(0), time, 0);
137
138 // Coarse MF to hold time interpolated data
139 amrex::MultiFab crse_data_time_interp(m_cf_crse_data_old->boxArray(),
140 m_cf_crse_data_old->DistributionMap(),
141 m_ncomp, amrex::IntVect{0});
142
143 // Time interpolate the coarse data
144 amrex::MultiFab::LinComb(crse_data_time_interp,
145 fac_old, *(m_cf_crse_data_old), 0,
146 fac_new, *(m_cf_crse_data_new), 0,
147 0, m_ncomp, amrex::IntVect{0});
148
149 // Call correct spatial interpolation type
150 amrex::IndexType m_ixt = mf.boxArray().ixType();
151 int ixt_sum = m_ixt[0]+m_ixt[1]+m_ixt[2];
152 if (ixt_sum == 0) {
153 InterpCell(mf,crse_data_time_interp,bcs,mask_val);
154 } else if (ixt_sum == 1) {
155 InterpFace(mf,crse_data_time_interp,mask_val);
156 } else {
157 amrex::Abort("REMORA_FillPatcher only supports face linear and cell cons linear interp!");
158 }
159}
160#endif
amrex::InterpBase * m_interp
void FillRelax(amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs)
Fill fine data in the relax region.
amrex::DistributionMapping m_fdm
amrex::iMultiFab * GetMask()
void Define(amrex::BoxArray const &fba, amrex::DistributionMapping const &fdm, amrex::Geometry const &fgeom, amrex::BoxArray const &cba, amrex::DistributionMapping const &cdm, amrex::Geometry const &cgeom, int nghost, int nghost_set, int ncomp, amrex::InterpBase *interp)
Redefine the coarse and fine patch MultiFabs.
void RegisterCoarseData(amrex::Vector< amrex::MultiFab const * > const &crse_data, amrex::Vector< amrex::Real > const &crse_time)
Register the coarse data to be used by the REMORAFillPatcher.
amrex::IntVect m_ratio
void InterpCell(amrex::MultiFab &fine, amrex::MultiFab const &crse, amrex::Vector< amrex::BCRec > const &bcr, int mask_val)
Interpolate to cell centers.
amrex::Geometry m_fgeom
amrex::BoxArray m_cba
amrex::BoxArray m_fba
void Fill(amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs, int mask_val)
Fill fine data in the relax region.
amrex::Vector< amrex::Real > m_crse_times
void InterpFace(amrex::MultiFab &fine, amrex::MultiFab const &crse, int mask_val)
Interpolate to cell faces.
std::unique_ptr< amrex::MultiFab > m_cf_crse_data_new
void BuildMask(amrex::BoxArray const &fba, int nghost, int nghost_set)
Generate masking array.
amrex::DistributionMapping m_cdm
amrex::Geometry m_cgeom
std::unique_ptr< amrex::iMultiFab > m_cf_mask
std::unique_ptr< amrex::MultiFab > m_cf_crse_data_old
void FillSet(amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs)
Fill fine data in the set region.