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