Block-Structured AMR Software Framework
AMReX_StateData.H
Go to the documentation of this file.
1 
2 #ifndef AMREX_StateData_H_
3 #define AMREX_StateData_H_
4 #include <AMReX_Config.H>
5 
6 #include <AMReX_Box.H>
7 #include <AMReX_BoxArray.H>
8 #include <AMReX_MultiFab.H>
10 #include <AMReX_BCRec.H>
11 #include <AMReX_Array.H>
12 #include <AMReX_Vector.H>
13 #include <AMReX_VisMF.H>
15 #include <AMReX_PhysBCFunct.H>
16 #include <AMReX_Geometry.H>
17 #include <AMReX_RealBox.H>
18 #include <AMReX_StateDescriptor.H>
19 
20 #include <memory>
21 
22 namespace amrex {
23 
24 class StateDataPhysBCFunct;
25 
26 
32 class StateData
33 {
34  friend class StateDataPhysBCFunct;
35 
36 public:
37 
41  StateData ();
42 
54  StateData (const Box& p_domain,
55  const BoxArray& grds,
56  const DistributionMapping& dm,
57  const StateDescriptor* d,
58  Real cur_time,
59  Real dt,
60  const FabFactory<FArrayBox>& factory);
61 
65  ~StateData ();
66 
67  StateData (StateData&& rhs) noexcept;
68  StateData& operator= (StateData const& rhs);
69 
70  StateData (StateData const& rhs) = delete;
71  StateData& operator= (StateData && rhs) = delete;
72 
84  void define (const Box& p_domain,
85  const BoxArray& grds,
86  const DistributionMapping& dm,
87  const StateDescriptor& d,
88  Real cur_time,
89  Real dt,
90  const FabFactory<FArrayBox>& factory);
91 
98  void copyOld (const StateData& state);
99 
106  void copyNew (const StateData& state);
107 
111  void allocOldData ();
112 
116  void removeOldData () { old_data.reset(); }
117 
121  void reset ();
122 
128  void swapTimeLevels (Real dt);
129 
136  void replaceOldData ( MultiFab&& mf );
137 
144  void replaceOldData ( StateData& s );
145 
152  void replaceNewData ( MultiFab&& mf );
153 
160  void replaceNewData ( StateData& s );
161 
162 
170  void setTimeLevel (Real time,
171  Real dt_old,
172  Real dt_new);
173 
179  void setOldTimeLevel (Real time);
180 
186  void setNewTimeLevel (Real time);
187 
188  void syncNewTimeLevel (Real time);
189 
190  void RegisterData (MultiFabCopyDescriptor& multiFabCopyDesc,
191  Vector<MultiFabId>& mfid);
192 
193  void InterpAddBox (MultiFabCopyDescriptor& multiFabCopyDesc,
194  Vector<MultiFabId>& mfid,
195  BoxList* returnedUnfillableBoxes,
196  Vector<FillBoxId>& returnedFillBoxIds,
197  const Box& subbox,
198  Real time,
199  int src_comp,
200  int dest_comp,
201  int num_comp,
202  bool extrap = false);
203 
204  void InterpFillFab (MultiFabCopyDescriptor& fabCopyDesc,
205  const Vector<MultiFabId>& mfid,
206  const Vector<FillBoxId>& fillBoxIds,
207  FArrayBox& dest,
208  Real time,
209  int src_comp,
210  int dest_comp,
211  int num_comp,
212  bool extrap = false);
213 
214 
226  void FillBoundary (FArrayBox& dest,
227  Real time,
228  const Real* dx,
229  const RealBox& prob_domain,
230  int dest_comp,
231  int src_comp,
232  int num_comp = 1);
233 
234  void FillBoundary (Box const& bx,
235  FArrayBox& dest,
236  Real time,
237  Geometry const& geom,
238  int dest_comp,
239  int src_comp,
240  int num_comp);
241 
251  void checkPoint (const std::string& name,
252  const std::string& fullpathname,
253  std::ostream& os,
254  VisMF::How how,
255  bool dump_old = true);
256 
268  void restart (std::istream& is,
269  const Box& p_domain,
270  const BoxArray& grds,
271  const DistributionMapping& dm,
272  const FabFactory<FArrayBox>& factory,
273  const StateDescriptor& d,
274  const std::string& chkfile);
275 
282  void restart (const StateDescriptor& d,
283  const StateData& rhs);
284 
288  const StateDescriptor* descriptor () const noexcept { return desc; }
289 
293  const Box& getDomain () const noexcept { return domain; }
294 
298  const BoxArray& boxArray () const noexcept { return grids; }
299 
300  const DistributionMapping& DistributionMap () const noexcept { return dmap; }
301 
306  void setDistributionMap ( DistributionMapping& new_dmap ) noexcept { dmap = new_dmap; }
307 
308  const FabFactory<FArrayBox>& Factory () const noexcept { return *m_factory; }
309 
313  Real curTime () const noexcept {
314  return (desc->timeType() == StateDescriptor::Point) ?
315  new_time.stop : 0.5_rt*(new_time.start + new_time.stop);
316  }
317 
321  Real prevTime () const noexcept {
322  return (desc->timeType() == StateDescriptor::Point) ?
323  old_time.stop : 0.5_rt*(old_time.start + old_time.stop);
324  }
325 
329  MultiFab& newData () noexcept { BL_ASSERT(new_data != nullptr); return *new_data; }
330 
334  const MultiFab& newData () const noexcept { BL_ASSERT(new_data != nullptr); return *new_data; }
335 
339  MultiFab& oldData () noexcept { BL_ASSERT(old_data != nullptr); return *old_data; }
340 
344  const MultiFab& oldData () const noexcept { BL_ASSERT(old_data != nullptr); return *old_data; }
345 
351  FArrayBox& newGrid (int i) noexcept { BL_ASSERT(new_data != nullptr); return (*new_data)[i]; }
352 
358  FArrayBox& oldGrid (int i) noexcept { BL_ASSERT(old_data != nullptr); return (*old_data)[i]; }
359 
366  BCRec getBC (int comp, int i) const noexcept;
367 
373  void printTimeInterval (std::ostream& os) const;
374 
378  bool hasOldData () const noexcept { return old_data != nullptr; }
379 
383  bool hasNewData () const noexcept { return new_data != nullptr; }
384 
385  void getData (Vector<MultiFab*>& data,
386  Vector<Real>& datatime,
387  Real time) const;
388 
392  Arena* getArena () const noexcept { return arena; }
393 
397  void setArena (Arena* ar) noexcept { arena = ar; }
398 
405 
406  static void SetFAHeaderMapPtr(std::map<std::string, Vector<char> > *fahmp) { faHeaderMap = fahmp; }
407 
408 
409 private:
410 
411  std::unique_ptr<FabFactory<FArrayBox> > m_factory;
412 
414  {
415  Real start;
416  Real stop;
417  };
418 
420  const StateDescriptor* desc{nullptr};
421 
424 
427  //
429 
432 
435 
437  std::unique_ptr<MultiFab> new_data;
438 
440  std::unique_ptr<MultiFab> old_data;
441 
443  Arena* arena{nullptr};
444 
450 
452  static std::map<std::string, Vector<char> > *faHeaderMap; // ---- [faheader name, the header]
453 
454  void restartDoit (std::istream& is, const std::string& chkfile);
455 };
456 
458 {
459 public:
460 
461  StateDataPhysBCFunct (StateData& sd, int sc, const Geometry& geom_);
462 
463  void operator() (MultiFab& mf, int dest_comp, int num_comp, IntVect const& nghost,
464  Real time, int bccomp);
465 
466  void FillBoundary (MultiFab& mf, int dest_comp, int num_comp, IntVect const& nghost,
467  Real time, int bccomp) {
468  this->operator()(mf,dest_comp,num_comp,nghost,time,bccomp); // NOLINT(readability-suspicious-call-argument)
469  }
470 private:
472  int src_comp;
473  const Geometry* geom;
474 
475 };
476 
477 }
478 
479 #endif /*_StateData_H_*/
#define BL_ASSERT(EX)
Definition: AMReX_BLassert.H:39
A virtual base class for objects that manage their own dynamic memory allocation.
Definition: AMReX_Arena.H:100
Boundary Condition Records. Necessary information and functions for computing boundary conditions.
Definition: AMReX_BCRec.H:17
A collection of Boxes stored in an Array.
Definition: AMReX_BoxArray.H:550
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition: AMReX_BoxList.H:52
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
A Fortran Array of REALs.
Definition: AMReX_FArrayBox.H:229
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
Definition: AMReX_MFCopyDescriptor.H:46
A collection (stored as an array) of FArrayBox objects.
Definition: AMReX_MultiFab.H:38
A Box with real dimensions. A RealBox is OK iff volume >= 0.
Definition: AMReX_RealBox.H:21
Definition: AMReX_StateData.H:458
void FillBoundary(MultiFab &mf, int dest_comp, int num_comp, IntVect const &nghost, Real time, int bccomp)
Definition: AMReX_StateData.H:466
StateData * statedata
Definition: AMReX_StateData.H:471
const Geometry * geom
Definition: AMReX_StateData.H:473
void operator()(MultiFab &mf, int dest_comp, int num_comp, IntVect const &nghost, Real time, int bccomp)
Definition: AMReX_StateData.cpp:870
int src_comp
Definition: AMReX_StateData.H:472
StateDataPhysBCFunct(StateData &sd, int sc, const Geometry &geom_)
Definition: AMReX_StateData.cpp:863
Current and previous level-time data.
Definition: AMReX_StateData.H:33
static const Vector< std::string > & FabArrayHeaderNames()
These facilitate prereading FabArray headers to avoid synchronization when reading multiple FabArrays...
Definition: AMReX_StateData.H:403
void removeOldData()
Deletes the space used by the old timestep data.
Definition: AMReX_StateData.H:116
std::unique_ptr< MultiFab > new_data
Pointer to new-time data.
Definition: AMReX_StateData.H:437
void replaceNewData(MultiFab &&mf)
Swaps new data with a new MultiFab. Deletes the previous new data.
Definition: AMReX_StateData.cpp:420
static void SetFAHeaderMapPtr(std::map< std::string, Vector< char > > *fahmp)
Definition: AMReX_StateData.H:406
void setOldTimeLevel(Real time)
Sets time of old data.
Definition: AMReX_StateData.cpp:327
void replaceOldData(MultiFab &&mf)
Swaps old data with a new MultiFab. Deletes the previous old data.
Definition: AMReX_StateData.cpp:406
void copyOld(const StateData &state)
Copies old data from another StateData object and sets the same time level. If old data is uninitiali...
Definition: AMReX_StateData.cpp:138
const DistributionMapping & DistributionMap() const noexcept
Definition: AMReX_StateData.H:300
MultiFab & newData() noexcept
Returns the new data.
Definition: AMReX_StateData.H:329
StateData & operator=(StateData const &rhs)
Definition: AMReX_StateData.cpp:65
FArrayBox & oldGrid(int i) noexcept
Returns the FAB of old data at grid index ā€˜iā€™.
Definition: AMReX_StateData.H:358
void restart(std::istream &is, const Box &p_domain, const BoxArray &grds, const DistributionMapping &dm, const FabFactory< FArrayBox > &factory, const StateDescriptor &d, const std::string &chkfile)
Restart with domain box, grids, and dmap provided.
Definition: AMReX_StateData.cpp:178
const MultiFab & oldData() const noexcept
Returns the old data.
Definition: AMReX_StateData.H:344
~StateData()
The destructor.
Definition: AMReX_StateData.cpp:302
Real curTime() const noexcept
Returns the current time.
Definition: AMReX_StateData.H:313
TimeInterval old_time
Time variable assoc with old data.
Definition: AMReX_StateData.H:434
void define(const Box &p_domain, const BoxArray &grds, const DistributionMapping &dm, const StateDescriptor &d, Real cur_time, Real dt, const FabFactory< FArrayBox > &factory)
Initializes data members if you used default constructor.
Definition: AMReX_StateData.cpp:92
Box domain
Problem domain.
Definition: AMReX_StateData.H:423
const StateDescriptor * desc
Pointer to data descriptor.
Definition: AMReX_StateData.H:420
Arena * getArena() const noexcept
Get the Arena used.
Definition: AMReX_StateData.H:392
const MultiFab & newData() const noexcept
Returns the new data.
Definition: AMReX_StateData.H:334
static std::map< std::string, Vector< char > > * faHeaderMap
This is used to store preread FabArray headers.
Definition: AMReX_StateData.H:452
void InterpFillFab(MultiFabCopyDescriptor &fabCopyDesc, const Vector< MultiFabId > &mfid, const Vector< FillBoxId > &fillBoxIds, FArrayBox &dest, Real time, int src_comp, int dest_comp, int num_comp, bool extrap=false)
Definition: AMReX_StateData.cpp:662
void syncNewTimeLevel(Real time)
Definition: AMReX_StateData.cpp:353
void setDistributionMap(DistributionMapping &new_dmap) noexcept
Definition: AMReX_StateData.H:306
void getData(Vector< MultiFab * > &data, Vector< Real > &datatime, Real time) const
Definition: AMReX_StateData.cpp:717
Arena * arena
Arena we should use for allocating the data.
Definition: AMReX_StateData.H:443
void FillBoundary(FArrayBox &dest, Real time, const Real *dx, const RealBox &prob_domain, int dest_comp, int src_comp, int num_comp=1)
Set physical bndry values.
Definition: AMReX_StateData.cpp:434
const BoxArray & boxArray() const noexcept
Returns the BoxArray.
Definition: AMReX_StateData.H:298
const FabFactory< FArrayBox > & Factory() const noexcept
Definition: AMReX_StateData.H:308
void checkPoint(const std::string &name, const std::string &fullpathname, std::ostream &os, VisMF::How how, bool dump_old=true)
Write the state data to a checkpoint file.
Definition: AMReX_StateData.cpp:773
BCRec getBC(int comp, int i) const noexcept
Returns boundary conditions of specified component on the specified grid.
Definition: AMReX_StateData.cpp:319
FArrayBox & newGrid(int i) noexcept
Returns the FAB of new data at grid index ā€˜iā€™.
Definition: AMReX_StateData.H:351
void setTimeLevel(Real time, Real dt_old, Real dt_new)
Sets time of old and new data.
Definition: AMReX_StateData.cpp:370
StateData()
The default constructor.
Definition: AMReX_StateData.cpp:31
std::unique_ptr< MultiFab > old_data
Pointer to previous time data.
Definition: AMReX_StateData.H:440
void setNewTimeLevel(Real time)
Sets time of new data.
Definition: AMReX_StateData.cpp:340
void RegisterData(MultiFabCopyDescriptor &multiFabCopyDesc, Vector< MultiFabId > &mfid)
Definition: AMReX_StateData.cpp:579
StateData(StateData const &rhs)=delete
DistributionMapping dmap
Definition: AMReX_StateData.H:428
void printTimeInterval(std::ostream &os) const
Prints out the time interval.
Definition: AMReX_StateData.cpp:849
Real prevTime() const noexcept
Returns the previous time.
Definition: AMReX_StateData.H:321
MultiFab & oldData() noexcept
Returns the old data.
Definition: AMReX_StateData.H:339
std::unique_ptr< FabFactory< FArrayBox > > m_factory
Definition: AMReX_StateData.H:411
void allocOldData()
Allocates space for old timestep data.
Definition: AMReX_StateData.cpp:308
void swapTimeLevels(Real dt)
Old data becomes new data and new time is incremented by dt.
Definition: AMReX_StateData.cpp:389
void reset()
Reverts back to initial state.
Definition: AMReX_StateData.cpp:170
static void ClearFabArrayHeaderNames()
Definition: AMReX_StateData.H:404
void copyNew(const StateData &state)
Copies new data from another StateData object and sets the same time level. If new data is uninitiali...
Definition: AMReX_StateData.cpp:154
TimeInterval new_time
Time variable assoc with new data.
Definition: AMReX_StateData.H:431
void InterpAddBox(MultiFabCopyDescriptor &multiFabCopyDesc, Vector< MultiFabId > &mfid, BoxList *returnedUnfillableBoxes, Vector< FillBoxId > &returnedFillBoxIds, const Box &subbox, Real time, int src_comp, int dest_comp, int num_comp, bool extrap=false)
Definition: AMReX_StateData.cpp:588
void restartDoit(std::istream &is, const std::string &chkfile)
Definition: AMReX_StateData.cpp:213
const Box & getDomain() const noexcept
Returns the valid domain.
Definition: AMReX_StateData.H:293
void setArena(Arena *ar) noexcept
Set the Arena used.
Definition: AMReX_StateData.H:397
static Vector< std::string > fabArrayHeaderNames
This is used as a temporary collection of FabArray header names written during a checkpoint.
Definition: AMReX_StateData.H:449
BoxArray grids
Grids defined at this level.
Definition: AMReX_StateData.H:426
const StateDescriptor * descriptor() const noexcept
Returns the StateDescriptor.
Definition: AMReX_StateData.H:288
bool hasOldData() const noexcept
True if there is any old data available.
Definition: AMReX_StateData.H:378
bool hasNewData() const noexcept
True if there is any new data available.
Definition: AMReX_StateData.H:383
Attributes of StateData.
Definition: AMReX_StateDescriptor.H:33
StateDescriptor::TimeCenter timeType() const noexcept
Returns StateDescriptor::TimeCenter.
Definition: AMReX_StateDescriptor.cpp:233
@ Point
Definition: AMReX_StateDescriptor.H:39
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition: AMReX_Vector.H:27
How
How we write out FabArray<FArrayBox>s. These are deprecated, we always use NFiles....
Definition: AMReX_VisMF.H:41
Definition: AMReX_Amr.cpp:49
Definition: AMReX_StateData.H:414
Real start
Definition: AMReX_StateData.H:415
Real stop
Definition: AMReX_StateData.H:416