Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
amrex::StateData Class Reference

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.
 
StateDataoperator= (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.
 
StateDataoperator= (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 StateDescriptordescriptor () const noexcept
 Returns the StateDescriptor.
 
const BoxgetDomain () const noexcept
 Returns the valid domain.
 
const BoxArrayboxArray () const noexcept
 Returns the BoxArray.
 
void setBoxArray (BoxArray const &a_grids) noexcept
 Replace the stored BoxArray layout with a_grids (no redistribution occurs).
 
const DistributionMappingDistributionMap () 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.
 
MultiFabnewData () noexcept
 Returns the new data.
 
const MultiFabnewData () const noexcept
 Returns the new data.
 
MultiFaboldData () noexcept
 Returns the old data.
 
const MultiFaboldData () const noexcept
 Returns the old data.
 
FArrayBoxnewGrid (int i) noexcept
 Returns the FAB of new data at grid index i.
 
FArrayBoxoldGrid (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.
 
ArenagetArena () 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
 

Detailed Description

Current and previous level-time data.

StateData holds state data on a level for the current and previous time step.

Constructor & Destructor Documentation

◆ StateData() [1/4]

amrex::StateData::StateData ( )

The default constructor.

◆ StateData() [2/4]

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.

Parameters
p_domainDomain box of the AMR problem (used for BCs).
grdsCell-centered BoxArray for this level.
dmDistribution mapping describing ownership of grids.
dState descriptor providing metadata (components, BCs, etc.).
cur_timeCurrent simulation time.
dtInitial time step (dt_old = dt_new = dt).
factoryFab factory controlling allocation (EB-aware, etc.).

◆ ~StateData()

amrex::StateData::~StateData ( )

The destructor.

◆ StateData() [3/4]

amrex::StateData::StateData ( StateData &&  rhs)
noexcept

Move constructor; transfers ownership of stored MultiFabs and metadata.

◆ StateData() [4/4]

amrex::StateData::StateData ( StateData const &  rhs)
delete

Copy construction is disabled to avoid ambiguous ownership of state buffers.

Member Function Documentation

◆ allocOldData()

void amrex::StateData::allocOldData ( )

Allocates space for old timestep data.

◆ boxArray()

const BoxArray & amrex::StateData::boxArray ( ) const
inlinenoexcept

Returns the BoxArray.

◆ checkPoint()

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.

Parameters
nameState name (prefix for component folders).
fullpathnameDirectory that will hold the data.
osOutput stream for metadata.
howVisMF write mode.
dump_oldWhether to include the old-time slot.

◆ ClearFabArrayHeaderNames()

static void amrex::StateData::ClearFabArrayHeaderNames ( )
inlinestatic

Clear the cached list of FabArray header names.

◆ copyNew()

void amrex::StateData::copyNew ( const StateData state)

Copy the "new" time slice from state, allocating storage if needed.

Parameters
stateSource StateData whose new data/time level are cloned.

◆ copyOld()

void amrex::StateData::copyOld ( const StateData state)

Copy the "old" time slice from state, allocating storage if needed.

Parameters
stateSource StateData whose old data/time level are cloned.

◆ curTime()

Real amrex::StateData::curTime ( ) const
inlinenoexcept

Returns the current time.

◆ define()

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.

Parameters
p_domainDomain box of the AMR problem (used for BCs).
grdsCell-centered BoxArray for this level.
dmDistribution mapping for that grid layout.
dState descriptor describing component layout.
cur_timeCurrent simulation time.
dtInitial dt (applied to both old/new slots).
factoryFab factory used for allocation.

◆ descriptor()

const StateDescriptor * amrex::StateData::descriptor ( ) const
inlinenoexcept

Returns the StateDescriptor.

◆ DistributionMap()

const DistributionMapping & amrex::StateData::DistributionMap ( ) const
inlinenoexcept

◆ FabArrayHeaderNames()

static const Vector< std::string > & amrex::StateData::FabArrayHeaderNames ( )
inlinestatic

These facilitate prereading FabArray headers to avoid synchronization when reading multiple FabArrays.

View of FabArray header filenames recorded during the last checkpoint write.

◆ Factory()

const FabFactory< FArrayBox > & amrex::StateData::Factory ( ) const
inlinenoexcept

◆ FillBoundary() [1/2]

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.

◆ FillBoundary() [2/2]

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.

Parameters
destDestination FAB to modify.
timePhysical time at which BCs are evaluated.
dxCell spacing array.
prob_domainReal-space domain extents.
dest_compDestination component offset in dest.
src_compSource component offset in the state data.
num_compNumber of components to update (default 1).

◆ getArena()

Arena * amrex::StateData::getArena ( ) const
inlinenoexcept

Get the Arena used.

◆ getBC()

BCRec amrex::StateData::getBC ( int  comp,
int  i 
) const
noexcept

Returns boundary conditions of specified component on the specified grid.

Parameters
compComponent index.
iGrid index.

◆ getData()

void amrex::StateData::getData ( Vector< MultiFab * > &  data,
Vector< Real > &  datatime,
Real  time 
) const

Gather data pointers and their time stamps corresponding to time.

Parameters
dataOutput list receiving pointers to MultiFabs.
datatimeOutput list receiving the associated times.
timeRequested interpolation time.

◆ getDomain()

const Box & amrex::StateData::getDomain ( ) const
inlinenoexcept

Returns the valid domain.

◆ hasNewData()

bool amrex::StateData::hasNewData ( ) const
inlinenoexcept

True if there is any new data available.

◆ hasOldData()

bool amrex::StateData::hasOldData ( ) const
inlinenoexcept

True if there is any old data available.

◆ InterpAddBox()

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.

Parameters
multiFabCopyDescDescriptor returned by RegisterData().
mfidIdentifiers previously produced by RegisterData().
returnedUnfillableBoxesOptional output list of boxes that could not be filled.
returnedFillBoxIdsOutput fill-box identifiers referencing queued requests.
subboxRegion that needs to be populated.
timePhysical time at which data should be sampled.
src_compSource component index in the state data.
dest_compDestination component index inside the caller's buffers.
num_compNumber of components requested.
extrapIf true, allow extrapolation beyond old/new time intervals.

◆ InterpFillFab()

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.

Parameters
fabCopyDescDescriptor produced by RegisterData().
mfidIdentifiers previously produced by RegisterData().
fillBoxIdsIdentifiers returned by InterpAddBox().
destDestination FAB.
timePhysical time at which data should be sampled.
src_compSource component offset in the state data.
dest_compDestination component offset inside dest.
num_compNumber of components to fill.
extrapIf true, allow extrapolation beyond old/new time intervals.

◆ newData() [1/2]

const MultiFab & amrex::StateData::newData ( ) const
inlinenoexcept

Returns the new data.

◆ newData() [2/2]

MultiFab & amrex::StateData::newData ( )
inlinenoexcept

Returns the new data.

◆ newGrid()

FArrayBox & amrex::StateData::newGrid ( int  i)
inlinenoexcept

Returns the FAB of new data at grid index i.

◆ oldData() [1/2]

const MultiFab & amrex::StateData::oldData ( ) const
inlinenoexcept

Returns the old data.

◆ oldData() [2/2]

MultiFab & amrex::StateData::oldData ( )
inlinenoexcept

Returns the old data.

◆ oldGrid()

FArrayBox & amrex::StateData::oldGrid ( int  i)
inlinenoexcept

Returns the FAB of old data at grid index i.

◆ operator=() [1/2]

StateData & amrex::StateData::operator= ( StateData &&  rhs)
delete

Move assignment is disabled (use define()/swap helpers instead).

◆ operator=() [2/2]

StateData & amrex::StateData::operator= ( StateData const &  rhs)

Copy-assignment; deep-copies metadata and allocates new FAB storage.

◆ prevTime()

Real amrex::StateData::prevTime ( ) const
inlinenoexcept

Returns the previous time.

◆ printTimeInterval()

void amrex::StateData::printTimeInterval ( std::ostream &  os) const

Prints out the time interval.

Parameters
osOutput stream receiving a textual summary of old/new times.

◆ RegisterData()

void amrex::StateData::RegisterData ( MultiFabCopyDescriptor multiFabCopyDesc,
Vector< MultiFabId > &  mfid 
)

Register the state data with a MultiFabCopyDescriptor for later interpolation/fills.

Parameters
multiFabCopyDescDescriptor that tracks MultiFabs participating in copy/interp operations.
mfidOutput vector receiving identifiers for the registered MultiFabs.

◆ removeOldData()

void amrex::StateData::removeOldData ( )
inline

Deletes the space used by the old timestep data.

◆ replaceNewData() [1/2]

void amrex::StateData::replaceNewData ( MultiFab &&  mf)

Swaps new data with a new MultiFab. Deletes the previous new data.

Parameters
mfMultiFab that will replace the stored new data (ownership transferred).

◆ replaceNewData() [2/2]

void amrex::StateData::replaceNewData ( StateData s)

Swaps new data with another StateData. Does not delete the previous new data.

Parameters
sPeer whose new data will be exchanged with this instance.

◆ replaceOldData() [1/2]

void amrex::StateData::replaceOldData ( MultiFab &&  mf)

Swaps old data with a new MultiFab. Deletes the previous old data.

Parameters
mfMultiFab that will replace the stored old data (ownership transferred).

◆ replaceOldData() [2/2]

void amrex::StateData::replaceOldData ( StateData s)

Swaps old data with another StateData. Does not delete the previous old data.

Parameters
sPeer whose old data will be exchanged with this instance.

◆ reset()

void amrex::StateData::reset ( )

Reverts back to initial state.

◆ restart() [1/2]

void amrex::StateData::restart ( const StateDescriptor d,
const StateData rhs 
)

Initialize from another StateData sharing descriptor d.

Parameters
dState descriptor describing the components.
rhsSource data that will be cloned.

◆ restart() [2/2]

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.

Parameters
isInput stream tied to the checkpoint data.
p_domainDomain box for the level.
grdsCell-centered grid layout.
dmDistribution mapping describing ownership of grids.
factoryFab factory controlling allocation strategy.
dState descriptor describing components/BCs.
chkfileName of the checkpoint file (for logging).

◆ setArena()

void amrex::StateData::setArena ( Arena ar)
inlinenoexcept

Set the Arena used for subsequent allocations to ar.

◆ setBoxArray()

void amrex::StateData::setBoxArray ( BoxArray const &  a_grids)
inlinenoexcept

Replace the stored BoxArray layout with a_grids (no redistribution occurs).

◆ setDistributionMap()

void amrex::StateData::setDistributionMap ( DistributionMapping new_dmap)
inlinenoexcept

Replace the stored DistributionMapping with new_dmap (no redistribution occurs).

◆ setFactory()

void amrex::StateData::setFactory ( FabFactory< FArrayBox > const &  a_factory)
inline

Replace the internal factory with a clone of a_factory.

◆ SetFAHeaderMapPtr()

static void amrex::StateData::SetFAHeaderMapPtr ( std::map< std::string, Vector< char > > *  fahmp)
inlinestatic

Install a pointer fahmp to the global map storing pre-read FabArray headers.

◆ setNewTimeLevel()

void amrex::StateData::setNewTimeLevel ( Real  time)

Sets time of new data.

Parameters
timeNew-time value to record.

◆ setOldTimeLevel()

void amrex::StateData::setOldTimeLevel ( Real  time)

Sets time of old data.

Parameters
timeOld-time value to record.

◆ setTimeLevel()

void amrex::StateData::setTimeLevel ( Real  time,
Real  dt_old,
Real  dt_new 
)

Sets time of old and new data.

Parameters
timeNew-time value.
dt_oldDt associated with the old time interval.
dt_newDt associated with the new time interval.

◆ swapTimeLevels()

void amrex::StateData::swapTimeLevels ( Real  dt)

Old data becomes new data and new time is incremented by dt.

Parameters
dtTime step separating the old and new slots.

◆ syncNewTimeLevel()

void amrex::StateData::syncNewTimeLevel ( Real  time)

Set both start/stop endpoints of the new-time interval to time.

Friends And Related Symbol Documentation

◆ StateDataPhysBCFunct

friend class StateDataPhysBCFunct
friend

The documentation for this class was generated from the following files: