Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
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>
19
25#include <memory>
26
27namespace amrex {
28
29class StateDataPhysBCFunct;
30
31
38{
40
41public:
42
46 StateData ();
47
59 StateData (const Box& p_domain,
60 const BoxArray& grds,
61 const DistributionMapping& dm,
62 const StateDescriptor* d,
63 Real cur_time,
64 Real dt,
65 const FabFactory<FArrayBox>& factory);
66
70 ~StateData ();
71
73 StateData (StateData&& rhs) noexcept;
75 StateData& operator= (StateData const& rhs);
76
78 StateData (StateData const& rhs) = delete;
80 StateData& operator= (StateData && rhs) = delete;
81
93 void define (const Box& p_domain,
94 const BoxArray& grds,
95 const DistributionMapping& dm,
96 const StateDescriptor& d,
97 Real cur_time,
98 Real dt,
99 const FabFactory<FArrayBox>& factory);
100
106 void copyOld (const StateData& state);
107
113 void copyNew (const StateData& state);
114
118 void allocOldData ();
119
123 void removeOldData () { old_data.reset(); }
124
128 void reset ();
129
135 void swapTimeLevels (Real dt);
136
143 void replaceOldData ( MultiFab&& mf );
144
151 void replaceOldData ( StateData& s );
152
159 void replaceNewData ( MultiFab&& mf );
160
167 void replaceNewData ( StateData& s );
168
169
177 void setTimeLevel (Real time,
178 Real dt_old,
179 Real dt_new);
180
186 void setOldTimeLevel (Real time);
187
193 void setNewTimeLevel (Real time);
194
196 void syncNewTimeLevel (Real time);
197
204 void RegisterData (MultiFabCopyDescriptor& multiFabCopyDesc,
205 Vector<MultiFabId>& mfid);
206
221 void InterpAddBox (MultiFabCopyDescriptor& multiFabCopyDesc,
222 Vector<MultiFabId>& mfid,
223 BoxList* returnedUnfillableBoxes,
224 Vector<FillBoxId>& returnedFillBoxIds,
225 const Box& subbox,
226 Real time,
227 int src_comp,
228 int dest_comp,
229 int num_comp,
230 bool extrap = false);
231
245 void InterpFillFab (MultiFabCopyDescriptor& fabCopyDesc,
246 const Vector<MultiFabId>& mfid,
247 const Vector<FillBoxId>& fillBoxIds,
248 FArrayBox& dest,
249 Real time,
250 int src_comp,
251 int dest_comp,
252 int num_comp,
253 bool extrap = false);
254
255
267 void FillBoundary (FArrayBox& dest,
268 Real time,
269 const Real* dx,
270 const RealBox& prob_domain,
271 int dest_comp,
272 int src_comp,
273 int num_comp = 1);
274
277 void FillBoundary (Box const& bx,
278 FArrayBox& dest,
279 Real time,
280 Geometry const& geom,
281 int dest_comp,
282 int src_comp,
283 int num_comp);
284
294 void checkPoint (const std::string& name,
295 const std::string& fullpathname,
296 std::ostream& os,
297 VisMF::How how,
298 bool dump_old = true);
299
311 void restart (std::istream& is,
312 const Box& p_domain,
313 const BoxArray& grds,
314 const DistributionMapping& dm,
315 const FabFactory<FArrayBox>& factory,
316 const StateDescriptor& d,
317 const std::string& chkfile);
318
325 void restart (const StateDescriptor& d,
326 const StateData& rhs);
327
331 const StateDescriptor* descriptor () const noexcept { return desc; }
332
336 const Box& getDomain () const noexcept { return domain; }
337
341 const BoxArray& boxArray () const noexcept { return grids; }
342
344 void setBoxArray (BoxArray const& a_grids) noexcept { grids = a_grids; }
345
346 const DistributionMapping& DistributionMap () const noexcept { return dmap; }
347
349 void setDistributionMap ( DistributionMapping& new_dmap ) noexcept { dmap = new_dmap; }
350
351 const FabFactory<FArrayBox>& Factory () const noexcept { return *m_factory; }
352
354 void setFactory (FabFactory<FArrayBox> const& a_factory) {
355 m_factory.reset(a_factory.clone());
356 }
357
361 Real curTime () const noexcept {
362 return (desc->timeType() == StateDescriptor::Point) ?
363 new_time.stop : 0.5_rt*(new_time.start + new_time.stop);
364 }
365
369 Real prevTime () const noexcept {
370 return (desc->timeType() == StateDescriptor::Point) ?
371 old_time.stop : 0.5_rt*(old_time.start + old_time.stop);
372 }
373
377 MultiFab& newData () noexcept { BL_ASSERT(new_data != nullptr); return *new_data; }
378
382 const MultiFab& newData () const noexcept { BL_ASSERT(new_data != nullptr); return *new_data; }
383
387 MultiFab& oldData () noexcept { BL_ASSERT(old_data != nullptr); return *old_data; }
388
392 const MultiFab& oldData () const noexcept { BL_ASSERT(old_data != nullptr); return *old_data; }
393
395 FArrayBox& newGrid (int i) noexcept { BL_ASSERT(new_data != nullptr); return (*new_data)[i]; }
396
398 FArrayBox& oldGrid (int i) noexcept { BL_ASSERT(old_data != nullptr); return (*old_data)[i]; }
399
406 BCRec getBC (int comp, int i) const noexcept;
407
413 void printTimeInterval (std::ostream& os) const;
414
418 bool hasOldData () const noexcept { return old_data != nullptr; }
419
423 bool hasNewData () const noexcept { return new_data != nullptr; }
424
432 void getData (Vector<MultiFab*>& data,
433 Vector<Real>& datatime,
434 Real time) const;
435
439 Arena* getArena () const noexcept { return arena; }
440
442 void setArena (Arena* ar) noexcept { arena = ar; }
443
449 static const Vector<std::string> &FabArrayHeaderNames() { return fabArrayHeaderNames; }
451 static void ClearFabArrayHeaderNames() { fabArrayHeaderNames.clear(); }
452
454 static void SetFAHeaderMapPtr(std::map<std::string, Vector<char> > *fahmp) { faHeaderMap = fahmp; }
455
456
457private:
458
459 std::unique_ptr<FabFactory<FArrayBox> > m_factory;
460
461 struct TimeInterval
462 {
463 Real start;
464 Real stop;
465 };
466
468 const StateDescriptor* desc{nullptr};
469
471 Box domain;
472
474 BoxArray grids;
475 //
476 DistributionMapping dmap;
477
479 TimeInterval new_time;
480
482 TimeInterval old_time;
483
485 std::unique_ptr<MultiFab> new_data;
486
488 std::unique_ptr<MultiFab> old_data;
489
491 Arena* arena{nullptr};
492
497 static Vector<std::string> fabArrayHeaderNames;
498
500 static std::map<std::string, Vector<char> > *faHeaderMap; // ---- [faheader name, the header]
501
503 void restartDoit (std::istream& is, const std::string& chkfile);
504};
505
507{
508public:
509
517 StateDataPhysBCFunct (StateData& sd, int sc, const Geometry& geom_);
518
529 void operator() (MultiFab& mf, int dest_comp, int num_comp, IntVect const& nghost,
530 Real time, int bccomp);
531
532 void FillBoundary (MultiFab& mf, int dest_comp, int num_comp, IntVect const& nghost,
533 Real time, int bccomp) {
534 this->operator()(mf,dest_comp,num_comp,nghost,time,bccomp); // NOLINT(readability-suspicious-call-argument)
535 }
536private:
537 StateData* statedata;
538 int src_comp;
539 const Geometry* geom;
540
541};
542
543}
544
545#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:127
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:568
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:43
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
Definition AMReX_FabFactory.H:50
virtual FabFactory< FAB > * clone() const =0
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
Definition AMReX_MFCopyDescriptor.H:46
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
A Box with real dimensions.
Definition AMReX_RealBox.H:28
Definition AMReX_StateData.H:507
void FillBoundary(MultiFab &mf, int dest_comp, int num_comp, IntVect const &nghost, Real time, int bccomp)
Definition AMReX_StateData.H:532
void operator()(MultiFab &mf, int dest_comp, int num_comp, IntVect const &nghost, Real time, int bccomp)
Apply BCs to mf for components [dest_comp, dest_comp+num_comp).
Definition AMReX_StateData.cpp:871
Current and previous level-time data.
Definition AMReX_StateData.H:38
void removeOldData()
Deletes the space used by the old timestep data.
Definition AMReX_StateData.H:123
MultiFab & oldData() noexcept
Returns the old data.
Definition AMReX_StateData.H:387
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)
Install a pointer fahmp to the global map storing pre-read FabArray headers.
Definition AMReX_StateData.H:454
MultiFab & newData() noexcept
Returns the new data.
Definition AMReX_StateData.H:377
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)
Copy the "old" time slice from state, allocating storage if needed.
Definition AMReX_StateData.cpp:138
const FabFactory< FArrayBox > & Factory() const noexcept
Definition AMReX_StateData.H:351
const Box & getDomain() const noexcept
Returns the valid domain.
Definition AMReX_StateData.H:336
StateData & operator=(StateData const &rhs)
Copy-assignment; deep-copies metadata and allocates new FAB storage.
Definition AMReX_StateData.cpp:65
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
~StateData()
The destructor.
Definition AMReX_StateData.cpp:302
Real curTime() const noexcept
Returns the current time.
Definition AMReX_StateData.H:361
void setFactory(FabFactory< FArrayBox > const &a_factory)
Replace the internal factory with a clone of a_factory.
Definition AMReX_StateData.H:354
const MultiFab & oldData() const noexcept
Returns the old data.
Definition AMReX_StateData.H:392
FArrayBox & newGrid(int i) noexcept
Returns the FAB of new data at grid index i.
Definition AMReX_StateData.H:395
const StateDescriptor * descriptor() const noexcept
Returns the StateDescriptor.
Definition AMReX_StateData.H:331
void define(const Box &p_domain, const BoxArray &grds, const DistributionMapping &dm, const StateDescriptor &d, Real cur_time, Real dt, const FabFactory< FArrayBox > &factory)
Initialize the object after default construction.
Definition AMReX_StateData.cpp:92
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)
Complete outstanding interpolation requests and fill dest.
Definition AMReX_StateData.cpp:663
void syncNewTimeLevel(Real time)
Set both start/stop endpoints of the new-time interval to time.
Definition AMReX_StateData.cpp:353
Arena * getArena() const noexcept
Get the Arena used.
Definition AMReX_StateData.H:439
const BoxArray & boxArray() const noexcept
Returns the BoxArray.
Definition AMReX_StateData.H:341
void setDistributionMap(DistributionMapping &new_dmap) noexcept
Replace the stored DistributionMapping with new_dmap (no redistribution occurs).
Definition AMReX_StateData.H:349
static const Vector< std::string > & FabArrayHeaderNames()
These facilitate prereading FabArray headers to avoid synchronization when reading multiple FabArrays...
Definition AMReX_StateData.H:449
void getData(Vector< MultiFab * > &data, Vector< Real > &datatime, Real time) const
Gather data pointers and their time stamps corresponding to time.
Definition AMReX_StateData.cpp:718
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:774
BCRec getBC(int comp, int i) const noexcept
Returns boundary conditions of specified component on the specified grid.
Definition AMReX_StateData.cpp:319
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
void setNewTimeLevel(Real time)
Sets time of new data.
Definition AMReX_StateData.cpp:340
void RegisterData(MultiFabCopyDescriptor &multiFabCopyDesc, Vector< MultiFabId > &mfid)
Register the state data with a MultiFabCopyDescriptor for later interpolation/fills.
Definition AMReX_StateData.cpp:580
StateData(StateData const &rhs)=delete
Copy construction is disabled to avoid ambiguous ownership of state buffers.
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_StateData.H:346
void printTimeInterval(std::ostream &os) const
Prints out the time interval.
Definition AMReX_StateData.cpp:850
Real prevTime() const noexcept
Returns the previous time.
Definition AMReX_StateData.H:369
const MultiFab & newData() const noexcept
Returns the new data.
Definition AMReX_StateData.H:382
void setBoxArray(BoxArray const &a_grids) noexcept
Replace the stored BoxArray layout with a_grids (no redistribution occurs).
Definition AMReX_StateData.H:344
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()
Clear the cached list of FabArray header names.
Definition AMReX_StateData.H:451
void copyNew(const StateData &state)
Copy the "new" time slice from state, allocating storage if needed.
Definition AMReX_StateData.cpp:154
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)
Request interpolation/copied data covering subbox.
Definition AMReX_StateData.cpp:589
void setArena(Arena *ar) noexcept
Set the Arena used for subsequent allocations to ar.
Definition AMReX_StateData.H:442
bool hasOldData() const noexcept
True if there is any old data available.
Definition AMReX_StateData.H:418
FArrayBox & oldGrid(int i) noexcept
Returns the FAB of old data at grid index i.
Definition AMReX_StateData.H:398
bool hasNewData() const noexcept
True if there is any new data available.
Definition AMReX_StateData.H:423
Attributes of StateData.
Definition AMReX_StateDescriptor.H:33
StateDescriptor::TimeCenter timeType() const noexcept
Returns StateDescriptor::TimeCenter.
Definition AMReX_StateDescriptor.cpp:235
@ 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:28
How
How we write out FabArray<FArrayBox>s. These are deprecated, we always use NFiles....
Definition AMReX_VisMF.H:42
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
Definition AMReX_Amr.cpp:49
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30