Alamo
BaseField.H
Go to the documentation of this file.
1 #ifndef INTEGRATOR_BASEFIELD_H
2 #define INTEGRATOR_BASEFIELD_H
3 
4 #include "AMReX_FillPatchUtil.H"
5 
6 #include "Util/Util.H"
7 #include "Set/Set.H"
9 
10 namespace Integrator
11 {
12 
13 class BaseField
14 {
15 public:
16  virtual void RemakeLevel (int lev,
17  amrex::Real time,
18  const amrex::BoxArray& cgrids,
19  const amrex::DistributionMapping& dm) = 0;
20  virtual void MakeNewLevelFromCoarse (int lev,
21  amrex::Real time,
22  const amrex::BoxArray& cgrids,
23  const amrex::DistributionMapping& dm) = 0;
24  virtual void MakeNewLevelFromScratch (int lev,
25  amrex::Real t,
26  const amrex::BoxArray& cgrids,
27  const amrex::DistributionMapping& dm) = 0;
28  virtual void SetFinestLevel(const int a_finestlevel) = 0;
29  virtual void FillPatch(const int lev, const Set::Scalar time)=0;
30  virtual void FillBoundary(const int lev, Set::Scalar time)=0;
31  virtual void AverageDown(const int lev, amrex::IntVect refRatio) = 0;
32 
33  virtual int NComp() = 0;
34  virtual void Copy(int /*a_lev*/, amrex::MultiFab &/*a_dst*/, int /*a_dstcomp*/, int /*a_nghost*/) = 0;
35  bool writeout = false;
36  virtual std::string Name(int) = 0;
37  virtual void setName(std::string a_name) = 0;
38  bool evolving = true;
39  virtual void setBC(void * a_bc) = 0;
40  virtual void * getBC() = 0;
42 };
43 
44 template<class T>
45 class Field : public BaseField
46 {
47 public:
48  class EmptyBC
49  {
50  public:
51  virtual void operator () (amrex::FabArray<amrex::BaseFab<T>> & /*mf*/,
52  int /*dcomp*/, int /*ncomp*/, amrex::IntVect const& /*nghost*/,
53  amrex::Real /*time*/, int /*bccomp*/)
54  { /* Do nothing - this is a shell to satisfy the FillPatch functions.*/ }
55  amrex::BCRec GetBCRec() {return amrex::BCRec();}
56  };
57 
58  Field(Set::Field<T> & a_field,
59  const amrex::Vector<amrex::Geometry> &a_geom,
60  const amrex::Vector<amrex::IntVect> &a_refRatio,
61  int a_ncomp, int a_nghost) :
62  m_field(a_field), m_geom(a_geom), m_refRatio(a_refRatio),
63  m_ncomp(a_ncomp), m_nghost(a_nghost)
64  {}
65 
66  void
67  FillPatch (int lev, amrex::Real time,
68  amrex::Vector<std::unique_ptr<amrex::FabArray<amrex::BaseFab<T>>>> &source_mf,
69  amrex::FabArray<amrex::BaseFab<T>> &destination_mf,
70  int icomp)
71  {
72  EmptyBC physbc;
73  if (lev == 0)
74  {
75  amrex::Vector<amrex::FabArray<amrex::BaseFab<T>>*> smf;
76  smf.push_back(source_mf[lev].get());
77  amrex::Vector<amrex::Real> stime;
78  stime.push_back(time);
79 
81  {
82  m_bc->define(m_geom[lev]);
83  amrex::FillPatchSingleLevel(destination_mf, time, smf, stime,
84  0, icomp, destination_mf.nComp(), m_geom[lev],
85  *m_bc, 0);
86  }
87  else if (m_gridtype == Set::Hypercube::Node)
88  {
89  amrex::FillPatchSingleLevel(destination_mf, time, smf, stime,
90  0, icomp, destination_mf.nComp(), m_geom[lev],
91  physbc, 0);
92  }
93  }
94  else
95  {
96  amrex::Vector<amrex::FabArray<amrex::BaseFab<T>>*> cmf, fmf;
97  cmf.push_back(source_mf[lev-1].get());
98  fmf.push_back(source_mf[lev].get());
99  amrex::Vector<amrex::Real> ctime, ftime;
100  ctime.push_back(time);
101  ftime.push_back(time);
102 
103  amrex::Vector<amrex::BCRec> bcs(destination_mf.nComp(), physbc.GetBCRec());
104 
105  if (destination_mf.boxArray().ixType() == amrex::IndexType::TheNodeType())
106  {
109 
110  amrex::FillPatchTwoLevels(destination_mf, time, cmf, ctime, fmf, ftime,
111  0, icomp, destination_mf.nComp(), m_geom[lev-1], m_geom[lev],
112  physbc, 0,
113  physbc, 0,
114  m_refRatio[lev-1],
115  mapper, bcs, 0);
116 
117  }
118  else
119  {
120  Util::Abort(INFO,"This implementation is not yet complete.");
121  /*
122  Numeric::Interpolator::CellBilinear<T> *mapper
123  = new Numeric::Interpolator::CellBilinear<T>();
124  m_bc->define(m_geom[lev]);
125  amrex::FillPatchTwoLevels(destination_mf, time, cmf, ctime, fmf, ftime,
126  0, icomp, destination_mf.nComp(), m_geom[lev-1], m_geom[lev],
127  *m_bc, 0,
128  *m_bc, 0,
129  m_refRatio[lev-1],
130  mapper, bcs, 0);
131  */
132  }
133  }
134 
135  }
136  void
137  FillPatch (int lev, amrex::Real time) override
138  {
139  this->FillPatch(lev,time,m_field,*(m_field[lev]),0);
140  }
141  void
142  FillCoarsePatch (int lev,
143  amrex::Real time,
144  int icomp,
145  int ncomp)
146  {
147  BL_PROFILE("Integrator::FillCoarsePatch");
148  AMREX_ASSERT(lev > 0);
149 
151  {
152  Util::Assert(INFO, TEST(m_field[lev]->boxArray().ixType() == amrex::IndexType::TheNodeType()));
153  EmptyBC physbc;
154  amrex::Vector<amrex::FabArray<amrex::BaseFab<T>> *> cmf;
155  cmf.push_back(m_field[lev-1].get());
156  amrex::Vector<amrex::Real> ctime;
157  ctime.push_back(time);
159  amrex::Vector<amrex::BCRec> bcs(ncomp, physbc.GetBCRec());
160  amrex::InterpFromCoarseLevel(*m_field[lev], time, *cmf[0], 0, icomp, ncomp,
161  m_geom[lev-1], m_geom[lev],
162  physbc, 0,
163  physbc, 0,
164  m_refRatio[lev-1],
165  mapper, bcs, 0);
166  }
167  else if (m_gridtype == Set::Hypercube::Cell)
168  {
169  Util::Abort(INFO,"This implementation is not yet complete");
170  /*
171  amrex::Vector<amrex::FabArray<amrex::BaseFab<T>> *> cmf;
172  cmf.push_back(m_field[lev-1].get());
173  amrex::Vector<amrex::Real> ctime;
174  ctime.push_back(time);
175  m_bc->define(m_geom[lev]);
176  Numeric::Interpolator::CellBilinear<T> *mapper = new Numeric::Interpolator::CellBilinear<T>();
177  amrex::Vector<amrex::BCRec> bcs(ncomp, m_bc->GetBCRec());
178  amrex::InterpFromCoarseLevel(*m_field[lev], time, *cmf[0], 0, icomp, ncomp,
179  m_geom[lev-1], m_geom[lev],
180  *m_bc, 0,
181  *m_bc, 0,
182  m_refRatio[lev-1],
183  mapper, bcs, 0);
184  */
185  }
186  else
187  Util::Abort(INFO,"Invalid grid");
188  }
189 
190 
191  virtual void RemakeLevel (int lev,
192  amrex::Real time,
193  const amrex::BoxArray& cgrids,
194  const amrex::DistributionMapping& dm) override
195  {
197  {
198  Util::Assert(INFO, TEST(m_field[lev]->boxArray().ixType() == amrex::IndexType::TheNodeType()));
199  amrex::BoxArray ngrids = cgrids;
200  ngrids.convert(amrex::IntVect::TheNodeVector());
201  amrex::FabArray<amrex::BaseFab<T>> new_state(ngrids, dm, m_ncomp, m_nghost);
202  this->FillPatch(lev, time, m_field, new_state, 0);
203  std::swap(new_state, *m_field[lev]);
204  }
205  else if (m_gridtype == Set::Hypercube::Cell)
206  {
207  amrex::FabArray<amrex::BaseFab<T>> new_state(cgrids, dm, m_ncomp, m_nghost);
208  this->FillPatch(lev, time, m_field, new_state, 0);
209  std::swap(new_state, *m_field[lev]);
210  }
211  else
212  Util::Abort(INFO,"Invalid grid");
213  }
214 
215  virtual void MakeNewLevelFromCoarse (int lev,
216  amrex::Real time,
217  const amrex::BoxArray& cgrids,
218  const amrex::DistributionMapping& dm) override
219  {
221  {
222  amrex::BoxArray ngrids = cgrids;
223  ngrids.convert(amrex::IntVect::TheNodeVector());
224  m_field[lev].reset(new amrex::FabArray<amrex::BaseFab<T>>(ngrids,dm,m_ncomp,m_nghost));
225  FillCoarsePatch(lev,time,0,m_ncomp);
226  }
227  else if (m_gridtype == Set::Hypercube::Cell)
228  {
229  m_field[lev].reset(new amrex::FabArray<amrex::BaseFab<T>>(cgrids,dm,m_ncomp,m_nghost));
230  FillCoarsePatch(lev,time,0,m_ncomp);
231  }
232  else
233  Util::Abort(INFO,"Invalid grid");
234  }
235 
236  virtual void MakeNewLevelFromScratch (int lev,
237  amrex::Real /*time*/,
238  const amrex::BoxArray& cgrids,
239  const amrex::DistributionMapping& dm) override
240  {
242  {
243  amrex::BoxArray ngrids = cgrids;
244  ngrids.convert(amrex::IntVect::TheNodeVector());
245  m_field[lev].reset(new amrex::FabArray<amrex::BaseFab<T>>(ngrids,dm,m_ncomp,m_nghost));
246  m_field[lev]->setVal(T::Zero());
247  }
248  else if (m_gridtype == Set::Hypercube::Cell)
249  {
250  m_field[lev].reset(new amrex::FabArray<amrex::BaseFab<T>>(cgrids,dm,m_ncomp,m_nghost));
251  m_field[lev]->setVal(T::Zero());
252  }
253  else
254  Util::Abort(INFO,"Invalid grid");
255  }
256 
257  virtual void SetFinestLevel(const int a_finestlevel) override
258  {
259  m_field.finest_level = a_finestlevel;
260  }
261 
262  virtual void FillBoundary(const int lev, Set::Scalar time) override
263  {
266  {
267  m_bc->define(m_geom[lev]);
268  for (amrex::MFIter mfi(*m_field[lev],true); mfi.isValid(); ++mfi)
269  {
270  const amrex::Box& box = mfi.tilebox();
271  amrex::BaseFab<T> &patch = (*(m_field[lev]))[mfi];
272  m_bc->FillBoundary(patch,box,m_nghost,0,m_ncomp,time);
273  }
274  //
275  }
276  }
277 
278  virtual void AverageDown(const int lev, amrex::IntVect refRatio) override
279  {
281  amrex::average_down_nodal(*m_field[lev+1],*m_field[lev],refRatio);
282  else if (m_gridtype == Set::Hypercube::Cell)
283  Util::Abort(INFO,"This is not yet working");
284  //amrex::average_down(*m_field[lev+1],*m_field[lev],0,m_ncomp,refRatio);
285  //amrex::average_down(*m_field[lev+1],*m_field[lev],0,m_ncomp,refRatio);
286  else
287  Util::Abort(INFO,"Incorrect grid");
288  }
289 
290 private:
291 
293  const amrex::Vector<amrex::Geometry> &m_geom;
294  const amrex::Vector<amrex::IntVect> &m_refRatio;
295  const int m_ncomp, m_nghost;
297 
298 public:
299  virtual int NComp() override {
300  return m_field.NComp();
301  }
302  virtual void Copy(int a_lev, amrex::MultiFab &a_dst, int a_dstcomp, int a_nghost) override
303  {
304  m_field.Copy(a_lev,a_dst,a_dstcomp,a_nghost);
305  }
306  virtual void setName(std::string a_name) override {
307  m_field.name = a_name;
308  }
309  virtual std::string Name(int i) override {
310  return m_field.Name(i);
311  }
312  virtual void setBC(void * a_bc) override {
313  m_bc = (BC::BC<T> *)(a_bc);
314  }
315  virtual void * getBC() override {
316  return (void*)m_bc;
317  }
318 
319 };
320 
321 }
322 
323 #endif
Integrator::Field::FillPatch
void FillPatch(int lev, amrex::Real time) override
Definition: BaseField.H:137
Integrator::BaseField::RemakeLevel
virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &cgrids, const amrex::DistributionMapping &dm)=0
Node
@ Node
Definition: Operator.H:13
Cell
@ Cell
Definition: Operator.H:13
Integrator::Field::AverageDown
virtual void AverageDown(const int lev, amrex::IntVect refRatio) override
Definition: BaseField.H:278
Util::RealFillBoundary
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T >> &a_mf, const amrex::Geometry &)
Definition: Util.H:307
Integrator::BaseField::getBC
virtual void * getBC()=0
BC::BC
Definition: BC.H:43
Util.H
Integrator::BaseField::NComp
virtual int NComp()=0
Integrator::Field::getBC
virtual void * getBC() override
Definition: BaseField.H:315
TEST
#define TEST(x)
Definition: Util.H:21
Integrator::Field::m_ncomp
const int m_ncomp
Definition: BaseField.H:295
Numeric::Interpolator::NodeBilinear
Definition: NodeBilinear.H:25
Integrator::Field::setName
virtual void setName(std::string a_name) override
Definition: BaseField.H:306
Integrator::BaseField::FillPatch
virtual void FillPatch(const int lev, const Set::Scalar time)=0
t
std::time_t t
Definition: FileNameParse.cpp:12
Integrator::Field::EmptyBC::operator()
virtual void operator()(amrex::FabArray< amrex::BaseFab< T >> &, int, int, amrex::IntVect const &, amrex::Real, int)
Definition: BaseField.H:51
Integrator::Field::EmptyBC::GetBCRec
amrex::BCRec GetBCRec()
Definition: BaseField.H:55
Integrator::BaseField::SetFinestLevel
virtual void SetFinestLevel(const int a_finestlevel)=0
Integrator::Field::MakeNewLevelFromCoarse
virtual void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &cgrids, const amrex::DistributionMapping &dm) override
Definition: BaseField.H:215
Integrator::BaseField
Definition: BaseField.H:13
Set::Hypercube
Hypercube
Definition: Set.H:31
Set::Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition: Base.H:20
Integrator::Field::m_refRatio
const amrex::Vector< amrex::IntVect > & m_refRatio
Definition: BaseField.H:294
Integrator::Field::EmptyBC
Definition: BaseField.H:48
Integrator::Field
Definition: BaseField.H:45
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
Integrator::Field::FillCoarsePatch
void FillCoarsePatch(int lev, amrex::Real time, int icomp, int ncomp)
Definition: BaseField.H:142
Integrator::BaseField::FillBoundary
virtual void FillBoundary(const int lev, Set::Scalar time)=0
Integrator::BaseField::evolving
bool evolving
Definition: BaseField.H:38
Integrator::BaseField::setName
virtual void setName(std::string a_name)=0
Integrator::BaseField::setBC
virtual void setBC(void *a_bc)=0
Integrator::Field::m_field
Set::Field< T > & m_field
Definition: BaseField.H:292
Integrator::Field::FillPatch
void FillPatch(int lev, amrex::Real time, amrex::Vector< std::unique_ptr< amrex::FabArray< amrex::BaseFab< T >>>> &source_mf, amrex::FabArray< amrex::BaseFab< T >> &destination_mf, int icomp)
Definition: BaseField.H:67
Integrator::BaseField::Copy
virtual void Copy(int, amrex::MultiFab &, int, int)=0
Integrator::Field::m_geom
const amrex::Vector< amrex::Geometry > & m_geom
Definition: BaseField.H:293
Integrator::Field::FillBoundary
virtual void FillBoundary(const int lev, Set::Scalar time) override
Definition: BaseField.H:262
Util::Assert
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition: Util.H:65
Integrator::BaseField::m_gridtype
Set::Hypercube m_gridtype
Definition: BaseField.H:41
Integrator::BaseField::Name
virtual std::string Name(int)=0
Integrator::Field::RemakeLevel
virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &cgrids, const amrex::DistributionMapping &dm) override
Definition: BaseField.H:191
NodeBilinear.H
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:165
Integrator
Collection of numerical integrator objects.
Definition: AllenCahn.H:39
Integrator::Field::m_nghost
const int m_nghost
Definition: BaseField.H:295
Integrator::Field::NComp
virtual int NComp() override
Definition: BaseField.H:299
Integrator::BaseField::MakeNewLevelFromScratch
virtual void MakeNewLevelFromScratch(int lev, amrex::Real t, const amrex::BoxArray &cgrids, const amrex::DistributionMapping &dm)=0
Integrator::Field::setBC
virtual void setBC(void *a_bc) override
Definition: BaseField.H:312
Set.H
Integrator::Field::Name
virtual std::string Name(int i) override
Definition: BaseField.H:309
Integrator::Field::Field
Field(Set::Field< T > &a_field, const amrex::Vector< amrex::Geometry > &a_geom, const amrex::Vector< amrex::IntVect > &a_refRatio, int a_ncomp, int a_nghost)
Definition: BaseField.H:58
Integrator::Field::SetFinestLevel
virtual void SetFinestLevel(const int a_finestlevel) override
Definition: BaseField.H:257
INFO
#define INFO
Definition: Util.H:20
Integrator::BaseField::AverageDown
virtual void AverageDown(const int lev, amrex::IntVect refRatio)=0
Set::Field
Definition: Set.H:42
Integrator::Field::MakeNewLevelFromScratch
virtual void MakeNewLevelFromScratch(int lev, amrex::Real, const amrex::BoxArray &cgrids, const amrex::DistributionMapping &dm) override
Definition: BaseField.H:236
Integrator::Field::m_bc
BC::BC< T > * m_bc
Definition: BaseField.H:296
Integrator::BaseField::MakeNewLevelFromCoarse
virtual void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &cgrids, const amrex::DistributionMapping &dm)=0
Integrator::Field::Copy
virtual void Copy(int a_lev, amrex::MultiFab &a_dst, int a_dstcomp, int a_nghost) override
Definition: BaseField.H:302
Integrator::BaseField::writeout
bool writeout
Definition: BaseField.H:35