![]() |
Block-Structured AMR Software Framework
|
Current and previous level-time data. More...
#include <AMReX_StateData.H>
Public Member Functions | |
| StateData () | |
| The default constructor. | |
| StateData (const Box &p_domain, const BoxArray &grds, const DistributionMapping &dm, const StateDescriptor *d, Real cur_time, Real dt, const FabFactory< FArrayBox > &factory) | |
| Construct and allocate state storage for a level. | |
| ~StateData () | |
| The destructor. | |
| StateData (StateData &&rhs) noexcept | |
| Move constructor; transfers ownership of stored MultiFabs and metadata. | |
| StateData & | operator= (StateData const &rhs) |
| Copy-assignment; deep-copies metadata and allocates new FAB storage. | |
| StateData (StateData const &rhs)=delete | |
| Copy construction is disabled to avoid ambiguous ownership of state buffers. | |
| StateData & | operator= (StateData &&rhs)=delete |
| Move assignment is disabled (use define()/swap helpers instead). | |
| 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. | |
| void | copyOld (const StateData &state) |
Copy the "old" time slice from state, allocating storage if needed. | |
| void | copyNew (const StateData &state) |
Copy the "new" time slice from state, allocating storage if needed. | |
| void | allocOldData () |
| Allocates space for old timestep data. | |
| void | removeOldData () |
| Deletes the space used by the old timestep data. | |
| void | reset () |
| Reverts back to initial state. | |
| void | swapTimeLevels (Real dt) |
Old data becomes new data and new time is incremented by dt. | |
| void | replaceOldData (MultiFab &&mf) |
| Swaps old data with a new MultiFab. Deletes the previous old data. | |
| void | replaceOldData (StateData &s) |
| Swaps old data with another StateData. Does not delete the previous old data. | |
| void | replaceNewData (MultiFab &&mf) |
| Swaps new data with a new MultiFab. Deletes the previous new data. | |
| void | replaceNewData (StateData &s) |
| Swaps new data with another StateData. Does not delete the previous new data. | |
| void | setTimeLevel (Real time, Real dt_old, Real dt_new) |
| Sets time of old and new data. | |
| void | setOldTimeLevel (Real time) |
| Sets time of old data. | |
| void | setNewTimeLevel (Real time) |
| Sets time of new data. | |
| void | syncNewTimeLevel (Real time) |
Set both start/stop endpoints of the new-time interval to time. | |
| void | RegisterData (MultiFabCopyDescriptor &multiFabCopyDesc, Vector< MultiFabId > &mfid) |
| Register the state data with a MultiFabCopyDescriptor for later interpolation/fills. | |
| 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. | |
| 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. | |
| 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 boundary values for a FAB extracted from the MultiFab. | |
| void | FillBoundary (Box const &bx, FArrayBox &dest, Real time, Geometry const &geom, int dest_comp, int src_comp, int num_comp) |
| 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. | |
| 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. | |
| void | restart (const StateDescriptor &d, const StateData &rhs) |
Initialize from another StateData sharing descriptor d. | |
| const StateDescriptor * | descriptor () const noexcept |
| Returns the StateDescriptor. | |
| const Box & | getDomain () const noexcept |
| Returns the valid domain. | |
| const BoxArray & | boxArray () const noexcept |
| Returns the BoxArray. | |
| void | setBoxArray (BoxArray const &a_grids) noexcept |
Replace the stored BoxArray layout with a_grids (no redistribution occurs). | |
| const DistributionMapping & | DistributionMap () const noexcept |
| void | setDistributionMap (DistributionMapping &new_dmap) noexcept |
Replace the stored DistributionMapping with new_dmap (no redistribution occurs). | |
| const FabFactory< FArrayBox > & | Factory () const noexcept |
| void | setFactory (FabFactory< FArrayBox > const &a_factory) |
Replace the internal factory with a clone of a_factory. | |
| Real | curTime () const noexcept |
| Returns the current time. | |
| Real | prevTime () const noexcept |
| Returns the previous time. | |
| MultiFab & | newData () noexcept |
| Returns the new data. | |
| const MultiFab & | newData () const noexcept |
| Returns the new data. | |
| MultiFab & | oldData () noexcept |
| Returns the old data. | |
| const MultiFab & | oldData () const noexcept |
| Returns the old data. | |
| FArrayBox & | newGrid (int i) noexcept |
Returns the FAB of new data at grid index i. | |
| FArrayBox & | oldGrid (int i) noexcept |
Returns the FAB of old data at grid index i. | |
| BCRec | getBC (int comp, int i) const noexcept |
| Returns boundary conditions of specified component on the specified grid. | |
| void | printTimeInterval (std::ostream &os) const |
| Prints out the time interval. | |
| bool | hasOldData () const noexcept |
| True if there is any old data available. | |
| bool | hasNewData () const noexcept |
| True if there is any new data available. | |
| void | getData (Vector< MultiFab * > &data, Vector< Real > &datatime, Real time) const |
Gather data pointers and their time stamps corresponding to time. | |
| Arena * | getArena () const noexcept |
| Get the Arena used. | |
| void | setArena (Arena *ar) noexcept |
Set the Arena used for subsequent allocations to ar. | |
Static Public Member Functions | |
| static const Vector< std::string > & | FabArrayHeaderNames () |
| These facilitate prereading FabArray headers to avoid synchronization when reading multiple FabArrays. | |
| static void | ClearFabArrayHeaderNames () |
| Clear the cached list of FabArray header names. | |
| static void | SetFAHeaderMapPtr (std::map< std::string, Vector< char > > *fahmp) |
Install a pointer fahmp to the global map storing pre-read FabArray headers. | |
Friends | |
| class | StateDataPhysBCFunct |
Current and previous level-time data.
StateData holds state data on a level for the current and previous time step.
| amrex::StateData::StateData | ( | ) |
The default constructor.
| amrex::StateData::StateData | ( | const Box & | p_domain, |
| const BoxArray & | grds, | ||
| const DistributionMapping & | dm, | ||
| const StateDescriptor * | d, | ||
| Real | cur_time, | ||
| Real | dt, | ||
| const FabFactory< FArrayBox > & | factory | ||
| ) |
Construct and allocate state storage for a level.
| p_domain | Domain box of the AMR problem (used for BCs). |
| grds | Cell-centered BoxArray for this level. |
| dm | Distribution mapping describing ownership of grids. |
| d | State descriptor providing metadata (components, BCs, etc.). |
| cur_time | Current simulation time. |
| dt | Initial time step (dt_old = dt_new = dt). |
| factory | Fab factory controlling allocation (EB-aware, etc.). |
| amrex::StateData::~StateData | ( | ) |
The destructor.
|
noexcept |
Move constructor; transfers ownership of stored MultiFabs and metadata.
|
delete |
Copy construction is disabled to avoid ambiguous ownership of state buffers.
| void amrex::StateData::allocOldData | ( | ) |
Allocates space for old timestep data.
| void amrex::StateData::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.
| name | State name (prefix for component folders). |
| fullpathname | Directory that will hold the data. |
| os | Output stream for metadata. |
| how | VisMF write mode. |
| dump_old | Whether to include the old-time slot. |
|
inlinestatic |
Clear the cached list of FabArray header names.
| void amrex::StateData::copyNew | ( | const StateData & | state | ) |
Copy the "new" time slice from state, allocating storage if needed.
| state | Source StateData whose new data/time level are cloned. |
| void amrex::StateData::copyOld | ( | const StateData & | state | ) |
Copy the "old" time slice from state, allocating storage if needed.
| state | Source StateData whose old data/time level are cloned. |
|
inlinenoexcept |
Returns the current time.
| void amrex::StateData::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.
| p_domain | Domain box of the AMR problem (used for BCs). |
| grds | Cell-centered BoxArray for this level. |
| dm | Distribution mapping for that grid layout. |
| d | State descriptor describing component layout. |
| cur_time | Current simulation time. |
| dt | Initial dt (applied to both old/new slots). |
| factory | Fab factory used for allocation. |
|
inlinenoexcept |
Returns the StateDescriptor.
|
inlinenoexcept |
|
inlinestatic |
|
inlinenoexcept |
| void amrex::StateData::FillBoundary | ( | Box const & | bx, |
| FArrayBox & | dest, | ||
| Real | time, | ||
| Geometry const & | geom, | ||
| int | dest_comp, | ||
| int | src_comp, | ||
| int | num_comp | ||
| ) |
Variant that fills data on region bx into dest at components starting dest_comp, sampling state components beginning at src_comp for num_comp values using geometry geom at time.
| void amrex::StateData::FillBoundary | ( | FArrayBox & | dest, |
| Real | time, | ||
| const Real * | dx, | ||
| const RealBox & | prob_domain, | ||
| int | dest_comp, | ||
| int | src_comp, | ||
| int | num_comp = 1 |
||
| ) |
Set physical boundary values for a FAB extracted from the MultiFab.
| dest | Destination FAB to modify. |
| time | Physical time at which BCs are evaluated. |
| dx | Cell spacing array. |
| prob_domain | Real-space domain extents. |
| dest_comp | Destination component offset in dest. |
| src_comp | Source component offset in the state data. |
| num_comp | Number of components to update (default 1). |
Returns boundary conditions of specified component on the specified grid.
| comp | Component index. |
| i | Grid index. |
| void amrex::StateData::getData | ( | Vector< MultiFab * > & | data, |
| Vector< Real > & | datatime, | ||
| Real | time | ||
| ) | const |
Gather data pointers and their time stamps corresponding to time.
| data | Output list receiving pointers to MultiFabs. |
| datatime | Output list receiving the associated times. |
| time | Requested interpolation time. |
|
inlinenoexcept |
Returns the valid domain.
|
inlinenoexcept |
True if there is any new data available.
|
inlinenoexcept |
True if there is any old data available.
| void amrex::StateData::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.
| multiFabCopyDesc | Descriptor returned by RegisterData(). |
| mfid | Identifiers previously produced by RegisterData(). |
| returnedUnfillableBoxes | Optional output list of boxes that could not be filled. |
| returnedFillBoxIds | Output fill-box identifiers referencing queued requests. |
| subbox | Region that needs to be populated. |
| time | Physical time at which data should be sampled. |
| src_comp | Source component index in the state data. |
| dest_comp | Destination component index inside the caller's buffers. |
| num_comp | Number of components requested. |
| extrap | If true, allow extrapolation beyond old/new time intervals. |
| void amrex::StateData::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.
| fabCopyDesc | Descriptor produced by RegisterData(). |
| mfid | Identifiers previously produced by RegisterData(). |
| fillBoxIds | Identifiers returned by InterpAddBox(). |
| dest | Destination FAB. |
| time | Physical time at which data should be sampled. |
| src_comp | Source component offset in the state data. |
| dest_comp | Destination component offset inside dest. |
| num_comp | Number of components to fill. |
| extrap | If true, allow extrapolation beyond old/new time intervals. |
|
inlinenoexcept |
Returns the new data.
|
inlinenoexcept |
Returns the new data.
Returns the FAB of new data at grid index i.
|
inlinenoexcept |
Returns the old data.
|
inlinenoexcept |
Returns the old data.
Returns the FAB of old data at grid index i.
Move assignment is disabled (use define()/swap helpers instead).
Copy-assignment; deep-copies metadata and allocates new FAB storage.
|
inlinenoexcept |
Returns the previous time.
| void amrex::StateData::printTimeInterval | ( | std::ostream & | os | ) | const |
Prints out the time interval.
| os | Output stream receiving a textual summary of old/new times. |
| void amrex::StateData::RegisterData | ( | MultiFabCopyDescriptor & | multiFabCopyDesc, |
| Vector< MultiFabId > & | mfid | ||
| ) |
Register the state data with a MultiFabCopyDescriptor for later interpolation/fills.
| multiFabCopyDesc | Descriptor that tracks MultiFabs participating in copy/interp operations. |
| mfid | Output vector receiving identifiers for the registered MultiFabs. |
|
inline |
Deletes the space used by the old timestep data.
| void amrex::StateData::replaceNewData | ( | MultiFab && | mf | ) |
| void amrex::StateData::replaceNewData | ( | StateData & | s | ) |
Swaps new data with another StateData. Does not delete the previous new data.
| s | Peer whose new data will be exchanged with this instance. |
| void amrex::StateData::replaceOldData | ( | MultiFab && | mf | ) |
| void amrex::StateData::replaceOldData | ( | StateData & | s | ) |
Swaps old data with another StateData. Does not delete the previous old data.
| s | Peer whose old data will be exchanged with this instance. |
| void amrex::StateData::reset | ( | ) |
Reverts back to initial state.
| void amrex::StateData::restart | ( | const StateDescriptor & | d, |
| const StateData & | rhs | ||
| ) |
Initialize from another StateData sharing descriptor d.
| d | State descriptor describing the components. |
| rhs | Source data that will be cloned. |
| void amrex::StateData::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.
| is | Input stream tied to the checkpoint data. |
| p_domain | Domain box for the level. |
| grds | Cell-centered grid layout. |
| dm | Distribution mapping describing ownership of grids. |
| factory | Fab factory controlling allocation strategy. |
| d | State descriptor describing components/BCs. |
| chkfile | Name of the checkpoint file (for logging). |
|
inlinenoexcept |
Set the Arena used for subsequent allocations to ar.
|
inlinenoexcept |
Replace the stored BoxArray layout with a_grids (no redistribution occurs).
|
inlinenoexcept |
Replace the stored DistributionMapping with new_dmap (no redistribution occurs).
|
inline |
Replace the internal factory with a clone of a_factory.
|
inlinestatic |
Install a pointer fahmp to the global map storing pre-read FabArray headers.
| void amrex::StateData::setNewTimeLevel | ( | Real | time | ) |
Sets time of new data.
| time | New-time value to record. |
| void amrex::StateData::setOldTimeLevel | ( | Real | time | ) |
Sets time of old data.
| time | Old-time value to record. |
Sets time of old and new data.
| time | New-time value. |
| dt_old | Dt associated with the old time interval. |
| dt_new | Dt associated with the new time interval. |
| void amrex::StateData::swapTimeLevels | ( | Real | dt | ) |
Old data becomes new data and new time is incremented by dt.
| dt | Time step separating the old and new slots. |
| void amrex::StateData::syncNewTimeLevel | ( | Real | time | ) |
Set both start/stop endpoints of the new-time interval to time.
|
friend |