IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
Class LDataManager coordinates the irregular distribution of LNode and LData on the patch hierarchy. More...
#include </home/runner/work/IBAMR/IBAMR/ibtk/include/ibtk/LDataManager.h>
Static Public Member Functions | |
static LDataManager * | getManager (const std::string &name, const std::string &default_interp_kernel_fcn, const std::string &default_spread_kernel_fcn, bool error_if_points_leave_domain=false, const SAMRAI::hier::IntVector< NDIM > &min_ghost_width=SAMRAI::hier::IntVector< NDIM >(0), bool register_for_restart=true) |
static void | freeAllManagers () |
Static Public Attributes | |
static const std::string | POSN_DATA_NAME = "X" |
static const std::string | INIT_POSN_DATA_NAME = "X0" |
static const std::string | VEL_DATA_NAME = "U" |
Methods to set and get the patch hierarchy and range of patch | |
void | setPatchHierarchy (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy) |
Reset patch hierarchy over which operations occur. | |
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > | getPatchHierarchy () const |
Get the patch hierarchy used by this object. | |
void | setPatchLevels (int coarsest_ln, int finest_ln) |
Reset range of patch levels over which operations occur. More... | |
std::pair< int, int > | getPatchLevels () const |
Get the range of patch levels used by this object. More... | |
const SAMRAI::hier::IntVector< NDIM > & | getGhostCellWidth () const |
Return the ghost cell width associated with the interaction scheme. | |
const std::string & | getDefaultInterpKernelFunction () const |
Return the default kernel function associated with the Eulerian-to-Lagrangian interpolation scheme. | |
const std::string & | getDefaultSpreadKernelFunction () const |
Return the default kernel function associated with the Lagrangian-to-Eulerian spreading scheme. | |
void | spread (int f_data_idx, SAMRAI::tbox::Pointer< LData > F_data, SAMRAI::tbox::Pointer< LData > X_data, SAMRAI::tbox::Pointer< LData > ds_data, RobinPhysBdryPatchStrategy *f_phys_bdry_op, int level_num, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0, bool F_data_ghost_node_update=true, bool X_data_ghost_node_update=true, bool ds_data_ghost_node_update=true) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the default spreading kernel function. More... | |
void | spread (int f_data_idx, SAMRAI::tbox::Pointer< LData > F_data, SAMRAI::tbox::Pointer< LData > X_data, SAMRAI::tbox::Pointer< LData > ds_data, const std::string &spread_kernel_fcn, RobinPhysBdryPatchStrategy *f_phys_bdry_op, int level_num, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0, bool F_data_ghost_node_update=true, bool X_data_ghost_node_update=true, bool ds_data_ghost_node_update=true) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using a specified spreading kernel function. More... | |
void | spread (int f_data_idx, std::vector< SAMRAI::tbox::Pointer< LData > > &F_data, std::vector< SAMRAI::tbox::Pointer< LData > > &X_data, std::vector< SAMRAI::tbox::Pointer< LData > > &ds_data, RobinPhysBdryPatchStrategy *f_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0, bool F_data_ghost_node_update=true, bool X_data_ghost_node_update=true, bool ds_data_ghost_node_update=true, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the default spreading kernel function. More... | |
void | spread (int f_data_idx, std::vector< SAMRAI::tbox::Pointer< LData > > &F_data, std::vector< SAMRAI::tbox::Pointer< LData > > &X_data, std::vector< SAMRAI::tbox::Pointer< LData > > &ds_data, const std::string &spread_kernel_fcn, RobinPhysBdryPatchStrategy *f_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0, bool F_data_ghost_node_update=true, bool X_data_ghost_node_update=true, bool ds_data_ghost_node_update=true, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using a specified spreading kernel function. More... | |
void | spread (int f_data_idx, SAMRAI::tbox::Pointer< LData > F_data, SAMRAI::tbox::Pointer< LData > X_data, RobinPhysBdryPatchStrategy *f_phys_bdry_op, int level_num, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0, bool F_data_ghost_node_update=true, bool X_data_ghost_node_update=true) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the default spreading kernel function. More... | |
void | spread (int f_data_idx, SAMRAI::tbox::Pointer< LData > F_data, SAMRAI::tbox::Pointer< LData > X_data, const std::string &spread_kernel_fcn, RobinPhysBdryPatchStrategy *f_phys_bdry_op, int level_num, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0, bool F_data_ghost_node_update=true, bool X_data_ghost_node_update=true) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the specified spreading kernel function. More... | |
void | spread (int f_data_idx, std::vector< SAMRAI::tbox::Pointer< LData > > &F_data, std::vector< SAMRAI::tbox::Pointer< LData > > &X_data, RobinPhysBdryPatchStrategy *f_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0, bool F_data_ghost_node_update=true, bool X_data_ghost_node_update=true, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the default spreading kernel function. More... | |
void | spread (int f_data_idx, std::vector< SAMRAI::tbox::Pointer< LData > > &F_data, std::vector< SAMRAI::tbox::Pointer< LData > > &X_data, const std::string &spread_kernel_fcn, RobinPhysBdryPatchStrategy *f_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0, bool F_data_ghost_node_update=true, bool X_data_ghost_node_update=true, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the specified spreading kernel function. More... | |
void | interp (int f_data_idx, SAMRAI::tbox::Pointer< LData > F_data, SAMRAI::tbox::Pointer< LData > X_data, int level_num, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &f_synch_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > >(), const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_ghost_fill_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0) |
Interpolate a quantity from the Eulerian grid to the Lagrangian mesh using the default interpolation kernel function. | |
void | interp (int f_data_idx, SAMRAI::tbox::Pointer< LData > F_data, SAMRAI::tbox::Pointer< LData > X_data, const std::string &interp_kernel_fcn, int level_num, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &f_synch_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > >(), const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_ghost_fill_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0) |
Interpolate a quantity from the Eulerian grid to the Lagrangian mesh using the specified interpolation kernel function. | |
void | interp (int f_data_idx, std::vector< SAMRAI::tbox::Pointer< LData > > &F_data, std::vector< SAMRAI::tbox::Pointer< LData > > &X_data, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &f_synch_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > >(), const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_ghost_fill_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
Interpolate a quantity from the Eulerian grid to the Lagrangian mesh using the default interpolation kernel function. | |
void | interp (int f_data_idx, std::vector< SAMRAI::tbox::Pointer< LData > > &F_data, std::vector< SAMRAI::tbox::Pointer< LData > > &X_data, const std::string &interp_kernel_fcn, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &f_synch_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > >(), const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_ghost_fill_scheds=std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > >(), double fill_data_time=0.0, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
Interpolate a quantity from the Eulerian grid to the Lagrangian mesh using the specified interpolation kernel function. | |
void | registerLInitStrategy (SAMRAI::tbox::Pointer< LInitStrategy > lag_init) |
void | freeLInitStrategy () |
void | registerVisItDataWriter (SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > visit_writer) |
Register a VisIt data writer with the manager. | |
void | registerLSiloDataWriter (SAMRAI::tbox::Pointer< LSiloDataWriter > silo_writer) |
Register a Silo data writer with the manager. | |
void | registerLoadBalancer (SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > load_balancer, int workload_data_idx) |
Register a load balancer for non-uniform load balancing. More... | |
bool | levelContainsLagrangianData (int level_number) const |
Indicates whether there is Lagrangian data on the given patch hierarchy level. | |
unsigned int | getNumberOfNodes (int level_number) const |
unsigned int | getNumberOfLocalNodes (int level_number) const |
unsigned int | getNumberOfGhostNodes (int level_number) const |
unsigned int | getGlobalNodeOffset (int level_number) const |
SAMRAI::tbox::Pointer< LMesh > | getLMesh (int level_number) const |
Get the Lagrangian mesh associated with the given patch hierarchy level. | |
SAMRAI::tbox::Pointer< LData > | getLData (const std::string &quantity_name, int level_number) const |
Get the specified Lagrangian quantity data on the given patch hierarchy level. | |
SAMRAI::tbox::Pointer< LData > | createLData (const std::string &quantity_name, int level_number, unsigned int depth=1, bool maintain_data=false) |
Allocate new Lagrangian level data with the specified name and depth. If specified, the quantity is maintained as the patch hierarchy evolves. More... | |
int | getLNodePatchDescriptorIndex () const |
Get the patch data descriptor index for the Lagrangian index data. | |
int | getWorkloadPatchDescriptorIndex () const |
Get the patch data descriptor index for the workload cell data. More... | |
int | getNodeCountPatchDescriptorIndex () const |
Get the patch data descriptor index for the Lagrangian node count cell data. | |
std::vector< std::string > | getLagrangianStructureNames (int level_number) const |
Get a list of Lagrangian structure names for the specified level of the patch hierarchy. | |
std::vector< int > | getLagrangianStructureIDs (int level_number) const |
Get a list of Lagrangian structure IDs for the specified level of the patch hierarchy. | |
int | getLagrangianStructureID (int lagrangian_index, int level_number) const |
Get the ID of the Lagrangian structure associated with the specified Lagrangian index. More... | |
int | getLagrangianStructureID (const std::string &structure_name, int level_number) const |
Get the ID of the Lagrangian structure with the specified name. More... | |
std::string | getLagrangianStructureName (int structure_id, int level_number) const |
Get the name of the Lagrangian structure with the specified ID. More... | |
std::pair< int, int > | getLagrangianStructureIndexRange (int structure_id, int level_number) const |
Get the range of Lagrangian indices for the Lagrangian structure with the specified ID. More... | |
Point | computeLagrangianStructureCenterOfMass (int structure_id, int level_number) |
Get the center of mass of the Lagrangian structure with the specified ID. More... | |
std::pair< Point, Point > | computeLagrangianStructureBoundingBox (int structure_id, int level_number) |
Get the bounding box of the Lagrangian structure with the specified ID. More... | |
void | reinitLagrangianStructure (const Point &X_center, int structure_id, int level_number) |
Reset the positions of the nodes of the Lagrangian structure with the specified ID to be equal to the initial positions but shifted so that the bounding box of the structure is centered about X_center. More... | |
void | displaceLagrangianStructure (const Vector &dX, int structure_id, int level_number) |
Shift the positions of the nodes of the Lagrangian structure with the specified ID by a displacement dX. More... | |
void | activateLagrangianStructures (const std::vector< int > &structure_ids, int level_number) |
Activate the Lagrangian structures with the specified ID numbers. More... | |
void | inactivateLagrangianStructures (const std::vector< int > &structure_ids, int level_number) |
Inactivate the Lagrangian structures with the specified ID numbers. More... | |
bool | getLagrangianStructureIsActivated (int structure_id, int level_number) const |
Determine whether the Lagrangian structure with the specified ID number is activated. | |
void | zeroInactivatedComponents (SAMRAI::tbox::Pointer< LData > lag_data, int level_number) const |
Set the components of the supplied LData object to zero for those entries that correspond to inactivated structures. | |
void | mapLagrangianToPETSc (std::vector< int > &inds, int level_number) const |
Map the collection of Lagrangian indices to the corresponding global PETSc indices. | |
void | mapPETScToLagrangian (std::vector< int > &inds, int level_number) const |
Map the collection of global PETSc indices to the corresponding Lagrangian indices. | |
AO & | getAO (int level_number) |
Return the AO object that maps Lagrangian indices to PETSc indices. | |
void | scatterLagrangianToPETSc (Vec &lagrangian_vec, Vec &petsc_vec, int level_number) const |
Scatter data from the Lagrangian ordering to the global PETSc ordering. More... | |
void | scatterPETScToLagrangian (Vec &petsc_vec, Vec &lagrangian_vec, int level_number) const |
Scatter data from the global PETSc ordering to the Lagrangian ordering. More... | |
void | scatterToAll (Vec ¶llel_vec, Vec &sequential_vec) const |
Scatter data from a distributed PETSc vector to all processors. More... | |
void | scatterToZero (Vec ¶llel_vec, Vec &sequential_vec) const |
Scatter data from a distributed PETSc vector to processor zero. More... | |
void | beginDataRedistribution (int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
Start the process of redistributing the Lagrangian data. More... | |
void | endDataRedistribution (int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
Finish the process of redistributing the Lagrangian data. More... | |
void | addWorkloadEstimate (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, const int workload_data_idx, const int coarsest_ln=invalid_level_number, const int finest_ln=invalid_level_number) |
Update the workload and count of nodes per cell. More... | |
void | updateNodeCountData (int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
Update the count of nodes per cell. More... | |
void | initializeLevelData (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double init_data_time, bool can_be_refined, bool initial_time, SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchLevel< NDIM > > old_level=SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchLevel< NDIM > >(nullptr), bool allocate_data=true) override |
void | resetHierarchyConfiguration (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int coarsest_ln, int finest_ln) override |
void | applyGradientDetector (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double error_data_time, int tag_index, bool initial_time, bool uses_richardson_extrapolation_too) override |
void | putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override |
void | registerUserDefinedLData (const std::string &data_name, int depth) |
LDataManager (std::string object_name, std::string default_interp_kernel_fcn, std::string default_spread_kernel_fcn, bool error_if_points_leave_domain, SAMRAI::hier::IntVector< NDIM > ghost_width, bool register_for_restart=true) | |
Constructor. | |
~LDataManager () | |
The LDataManager destructor cleans up any remaining PETSc AO objects. | |
Additional Inherited Members | |
Public Member Functions inherited from SAMRAI::mesh::StandardTagAndInitStrategy< NDIM > | |
virtual double | getLevelDt (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const double dt_time, const bool initial_time) |
virtual double | advanceLevel (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const tbox::Pointer< hier::BasePatchHierarchy< NDIM > > hierarchy, const double current_time, const double new_time, const bool first_step, const bool last_step, const bool regrid_advance=false) |
virtual void | resetTimeDependentData (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const double new_time, const bool can_be_refined) |
virtual void | resetDataToPreadvanceState (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level) |
virtual void | applyRichardsonExtrapolation (const tbox::Pointer< hier::PatchLevel< NDIM > > level, const double error_data_time, const int tag_index, const double deltat, const int error_coarsen_ratio, const bool initial_time, const bool uses_gradient_detector_too) |
virtual void | coarsenDataForRichardsonExtrapolation (const tbox::Pointer< hier::PatchHierarchy< NDIM > > hierarchy, const int level_number, const tbox::Pointer< hier::PatchLevel< NDIM > > coarser_level, const double coarsen_data_time, const bool before_advance) |
Class LDataManager coordinates the irregular distribution of LNode and LData on the patch hierarchy.
The manager class is responsible for maintaining this data distribution and for all inter-processor communications. All access to instantiated LDataManager objects is via the static method getManager().
Each Lagrangian point is associated with an LNode object which provides an interface between the Lagrangian and PETSc ordering as well as the data storage for force specification objects and other objects associated with the point. From each LNode object, we can retrieve the Lagrangian index as well as any associated node data. The LNode objects for each processor are contained in an LMesh object. The LMesh object contains two vectors of LNodes: one containing just the local LNode objects, and another containing 'ghost' LNodes. The 'ghost' LNodes are Lagrangian points assigned to other processors but have data needed for calculations. Interacting with Lagrangian data is mediated through the static LDataManager object. From this class, we can get the local LMesh object for each level of the Cartesian grid, and then loop through all associated Lagrangian nodes.
As an example, suppose we have a circle of Lagrangian points that are tethered to target points, and we wish to specify the motion of the target points, and hence the circle, to move in a straight line. We can do this by looping through all the nodes, pull out the IBAMR::IBTargetPointForceSpec, and then update the target point location. In parallel, we also need to update the target point locations of the ghost data.
For some applications, it might be useful to get the Eulerian data of the Lagrangian points. To do this, it is important to understand the Lagrangian and PETSc index ordering. The Lagrangian ordering is fixed at the beginning of the simulation and is set by the reading of the vertices, whether through an input file or programmatically. The PETSc ordering, however, will change over the course of the simulation. Whenever a regridding operation is triggered, the PETSc indexing is changed to ensure efficient memory usage. The Eulerian data is stored in PETSc ordering, so in order to access that data, we need to map Lagrangian indices to their corresponding PETSc indices. We can then access the corresponding components of the PETSc data vector which is wrapped in an LData object. As an example, here we loop over and print out the physical coordinates of all the Lagrangian points on the finest level.
void IBTK::LDataManager::activateLagrangianStructures | ( | const std::vector< int > & | structure_ids, |
int | level_number | ||
) |
Activate the Lagrangian structures with the specified ID numbers.
void IBTK::LDataManager::addWorkloadEstimate | ( | SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > | hierarchy, |
const int | workload_data_idx, | ||
const int | coarsest_ln = invalid_level_number , |
||
const int | finest_ln = invalid_level_number |
||
) |
Update the workload and count of nodes per cell.
This routine updates cell data that is maintained on the patch hierarchy to track the number of nodes in each cell of the AMR index space. The node count data is used to tag cells for refinement, and to specify non-uniform load balancing. The workload per cell is defined by
workload(i) = 1 + beta_work*node_count(i)
in which alpha and beta are parameters that each default to the value 1.
|
overridevirtual |
Set integer tags to "one" in cells where refinement of the given level should occur due to the presence of Lagrangian data. The double time argument is the regrid time. The integer "tag_index" argument is the patch descriptor index of the cell centered integer tag array on each patch in the hierarchy. The boolean argument initial_time indicates whether the level is being subject to refinement at the initial simulation time. If it is false, then the error estimation process is being invoked at some later time after the AMR hierarchy was initially constructed. The boolean argument uses_richardson_extrapolation_too is true when Richardson extrapolation error estimation is used in addition to the gradient detector, and false otherwise. This argument helps the user to manage multiple regridding criteria.
When assertion checking is active, an unrecoverable exception will result if the hierarchy pointer is null or the level number does not match any existing level in the hierarchy.
Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.
void IBTK::LDataManager::beginDataRedistribution | ( | int | coarsest_ln = invalid_level_number , |
int | finest_ln = invalid_level_number |
||
) |
Start the process of redistributing the Lagrangian data.
This method uses the present location of each Lagrangian mesh node to redistribute the LNodeData managed by this object.
std::pair< Point, Point > IBTK::LDataManager::computeLagrangianStructureBoundingBox | ( | int | structure_id, |
int | level_number | ||
) |
Get the bounding box of the Lagrangian structure with the specified ID.
Point IBTK::LDataManager::computeLagrangianStructureCenterOfMass | ( | int | structure_id, |
int | level_number | ||
) |
Get the center of mass of the Lagrangian structure with the specified ID.
X = (1/N) Sum_{k in structure} X_k
in which N is the number of nodes associated with that structure.
Pointer< LData > IBTK::LDataManager::createLData | ( | const std::string & | quantity_name, |
int | level_number, | ||
unsigned int | depth = 1 , |
||
bool | maintain_data = false |
||
) |
Allocate new Lagrangian level data with the specified name and depth. If specified, the quantity is maintained as the patch hierarchy evolves.
void IBTK::LDataManager::displaceLagrangianStructure | ( | const Vector & | dX, |
int | structure_id, | ||
int | level_number | ||
) |
Shift the positions of the nodes of the Lagrangian structure with the specified ID by a displacement dX.
void IBTK::LDataManager::endDataRedistribution | ( | int | coarsest_ln = invalid_level_number , |
int | finest_ln = invalid_level_number |
||
) |
Finish the process of redistributing the Lagrangian data.
This method redistributes the quantities associated with each node in the Lagrangian mesh according to the data distribution defined by the LNodeData managed by this object. This routine potentially involves SUBSTANTIAL inter-processor communication.
|
static |
Deallocate all of the LDataManager instances.
It is not necessary to call this function at program termination since it is automatically called by the ShutdownRegistry class.
void IBTK::LDataManager::freeLInitStrategy | ( | ) |
Free the concrete initialization strategy object.
|
inline |
|
inline |
Get the ID of the Lagrangian structure with the specified name.
|
inline |
Get the ID of the Lagrangian structure associated with the specified Lagrangian index.
|
inline |
Get the range of Lagrangian indices for the Lagrangian structure with the specified ID.
|
inline |
Get the name of the Lagrangian structure with the specified ID.
|
static |
Return a pointer to the instance of the Lagrangian data manager corresponding to the specified name. Access to LDataManager objects is mediated by the getManager() function.
Note that when a manager is accessed for the first time, the freeAllManagers static method is registered with the ShutdownRegistry class. Consequently, all allocated managers are freed at program completion. Thus, users of this class do not explicitly allocate or deallocate the LDataManager instances.
|
inline |
|
inline |
|
inline |
std::pair< int, int > IBTK::LDataManager::getPatchLevels | ( | ) | const |
Get the range of patch levels used by this object.
|
inline |
Get the patch data descriptor index for the workload cell data.
void IBTK::LDataManager::inactivateLagrangianStructures | ( | const std::vector< int > & | structure_ids, |
int | level_number | ||
) |
Inactivate the Lagrangian structures with the specified ID numbers.
|
overridevirtual |
Initialize data on a new level after it is inserted into an AMR patch hierarchy by the gridding algorithm. The level number indicates that of the new level. The old_level pointer corresponds to the level that resided in the hierarchy before the level with the specified number was introduced. If the pointer is null, there was no level in the hierarchy prior to the call and the level data is set based on the user routines and the simulation time. Otherwise, the specified level replaces the old level and the new level receives data from the old level appropriately before it is destroyed.
The boolean argument initial_time indicates whether the level is being introduced for the first time (i.e., at initialization time) or after some regrid process during the calculation beyond the initial hierarchy construction. This information is provided since the initialization of the data on a patch may be different in each of those circumstances. The can_be_refined boolean argument indicates whether the level is the finest level allowed in the hierarchy. This may or may not affect the data initialization process depending on the problem.
When assertion checking is active, an unrecoverable exception will result if the hierarchy pointer is null, the level number does not match any level in the hierarchy, or the old level number does not match the level number (if the old level pointer is non-null).
Implements SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.
|
overridevirtual |
Write out object state to the given database.
When assertion checking is active, database pointer must be non-null.
Implements SAMRAI::tbox::Serializable.
void IBTK::LDataManager::registerLInitStrategy | ( | SAMRAI::tbox::Pointer< LInitStrategy > | lag_init | ) |
Register a concrete strategy object with the integrator that specifies the initial configuration of the curvilinear mesh nodes.
void IBTK::LDataManager::registerLoadBalancer | ( | SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > | load_balancer, |
int | workload_data_idx | ||
) |
Register a load balancer for non-uniform load balancing.
void IBTK::LDataManager::registerUserDefinedLData | ( | const std::string & | data_name, |
int | depth | ||
) |
Register user defined Lagrangian data to be maintained
void IBTK::LDataManager::reinitLagrangianStructure | ( | const Point & | X_center, |
int | structure_id, | ||
int | level_number | ||
) |
Reset the positions of the nodes of the Lagrangian structure with the specified ID to be equal to the initial positions but shifted so that the bounding box of the structure is centered about X_center.
|
overridevirtual |
Reset cached communication schedules after the hierarchy has changed (for example, due to regridding) and the data has been initialized on the new levels. The intent is that the cost of data movement on the hierarchy will be amortized across multiple communication cycles, if possible. The level numbers indicate the range of levels in the hierarchy that have changed. However, this routine updates communication schedules every level finer than and including that indexed by the coarsest level number given.
When assertion checking is active, an unrecoverable exception will result if the hierarchy pointer is null, any pointer to a level in the hierarchy that is coarser than the finest level is null, or the given level numbers not specified properly; e.g., coarsest_ln > finest_ln.
Implements SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.
void IBTK::LDataManager::scatterLagrangianToPETSc | ( | Vec & | lagrangian_vec, |
Vec & | petsc_vec, | ||
int | level_number | ||
) | const |
Scatter data from the Lagrangian ordering to the global PETSc ordering.
void IBTK::LDataManager::scatterPETScToLagrangian | ( | Vec & | petsc_vec, |
Vec & | lagrangian_vec, | ||
int | level_number | ||
) | const |
Scatter data from the global PETSc ordering to the Lagrangian ordering.
void IBTK::LDataManager::scatterToAll | ( | Vec & | parallel_vec, |
Vec & | sequential_vec | ||
) | const |
Scatter data from a distributed PETSc vector to all processors.
void IBTK::LDataManager::scatterToZero | ( | Vec & | parallel_vec, |
Vec & | sequential_vec | ||
) | const |
Scatter data from a distributed PETSc vector to processor zero.
void IBTK::LDataManager::setPatchLevels | ( | int | coarsest_ln, |
int | finest_ln | ||
) |
Reset range of patch levels over which operations occur.
The levels must exist in the hierarchy or an assertion failure will result.
void IBTK::LDataManager::spread | ( | int | f_data_idx, |
SAMRAI::tbox::Pointer< LData > | F_data, | ||
SAMRAI::tbox::Pointer< LData > | X_data, | ||
const std::string & | spread_kernel_fcn, | ||
RobinPhysBdryPatchStrategy * | f_phys_bdry_op, | ||
int | level_num, | ||
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > & | f_prolongation_scheds = std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >() , |
||
double | fill_data_time = 0.0 , |
||
bool | F_data_ghost_node_update = true , |
||
bool | X_data_ghost_node_update = true |
||
) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the specified spreading kernel function.
f(i,j,k) = f(i,j,k) + Sum_{q,r,s} F(q,r,s) delta_h(x(i,j,k) - X(q,r,s))Unlike the standard regularized delta function spreading operation, the implemented operation spreads values, NOT densities.
void IBTK::LDataManager::spread | ( | int | f_data_idx, |
SAMRAI::tbox::Pointer< LData > | F_data, | ||
SAMRAI::tbox::Pointer< LData > | X_data, | ||
RobinPhysBdryPatchStrategy * | f_phys_bdry_op, | ||
int | level_num, | ||
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > & | f_prolongation_scheds = std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >() , |
||
double | fill_data_time = 0.0 , |
||
bool | F_data_ghost_node_update = true , |
||
bool | X_data_ghost_node_update = true |
||
) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the default spreading kernel function.
f(i,j,k) = f(i,j,k) + Sum_{q,r,s} F(q,r,s) delta_h(x(i,j,k) - X(q,r,s))Unlike the standard regularized delta function spreading operation, the implemented operation spreads values, NOT densities.
void IBTK::LDataManager::spread | ( | int | f_data_idx, |
SAMRAI::tbox::Pointer< LData > | F_data, | ||
SAMRAI::tbox::Pointer< LData > | X_data, | ||
SAMRAI::tbox::Pointer< LData > | ds_data, | ||
const std::string & | spread_kernel_fcn, | ||
RobinPhysBdryPatchStrategy * | f_phys_bdry_op, | ||
int | level_num, | ||
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > & | f_prolongation_scheds = std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >() , |
||
double | fill_data_time = 0.0 , |
||
bool | F_data_ghost_node_update = true , |
||
bool | X_data_ghost_node_update = true , |
||
bool | ds_data_ghost_node_update = true |
||
) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using a specified spreading kernel function.
f(i,j,k) = f(i,j,k) + Sum_{q,r,s} F(q,r,s) delta_h(x(i,j,k) - X(q,r,s)) ds(q,r,s)This is the standard regularized delta function spreading operation, which spreads densities, NOT values.
void IBTK::LDataManager::spread | ( | int | f_data_idx, |
SAMRAI::tbox::Pointer< LData > | F_data, | ||
SAMRAI::tbox::Pointer< LData > | X_data, | ||
SAMRAI::tbox::Pointer< LData > | ds_data, | ||
RobinPhysBdryPatchStrategy * | f_phys_bdry_op, | ||
int | level_num, | ||
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > & | f_prolongation_scheds = std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >() , |
||
double | fill_data_time = 0.0 , |
||
bool | F_data_ghost_node_update = true , |
||
bool | X_data_ghost_node_update = true , |
||
bool | ds_data_ghost_node_update = true |
||
) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the default spreading kernel function.
f(i,j,k) = f(i,j,k) + Sum_{q,r,s} F(q,r,s) delta_h(x(i,j,k) - X(q,r,s)) ds(q,r,s)This is the standard regularized delta function spreading operation, which spreads densities, NOT values.
void IBTK::LDataManager::spread | ( | int | f_data_idx, |
std::vector< SAMRAI::tbox::Pointer< LData > > & | F_data, | ||
std::vector< SAMRAI::tbox::Pointer< LData > > & | X_data, | ||
const std::string & | spread_kernel_fcn, | ||
RobinPhysBdryPatchStrategy * | f_phys_bdry_op, | ||
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > & | f_prolongation_scheds = std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >() , |
||
double | fill_data_time = 0.0 , |
||
bool | F_data_ghost_node_update = true , |
||
bool | X_data_ghost_node_update = true , |
||
int | coarsest_ln = invalid_level_number , |
||
int | finest_ln = invalid_level_number |
||
) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the specified spreading kernel function.
f(i,j,k) = f(i,j,k) + Sum_{q,r,s} F(q,r,s) delta_h(x(i,j,k) - X(q,r,s))Unlike the standard regularized delta function spreading operation, the implemented operation spreads values, NOT densities.
void IBTK::LDataManager::spread | ( | int | f_data_idx, |
std::vector< SAMRAI::tbox::Pointer< LData > > & | F_data, | ||
std::vector< SAMRAI::tbox::Pointer< LData > > & | X_data, | ||
RobinPhysBdryPatchStrategy * | f_phys_bdry_op, | ||
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > & | f_prolongation_scheds = std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >() , |
||
double | fill_data_time = 0.0 , |
||
bool | F_data_ghost_node_update = true , |
||
bool | X_data_ghost_node_update = true , |
||
int | coarsest_ln = invalid_level_number , |
||
int | finest_ln = invalid_level_number |
||
) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the default spreading kernel function.
f(i,j,k) = f(i,j,k) + Sum_{q,r,s} F(q,r,s) delta_h(x(i,j,k) - X(q,r,s))Unlike the standard regularized delta function spreading operation, the implemented operation spreads values, NOT densities.
void IBTK::LDataManager::spread | ( | int | f_data_idx, |
std::vector< SAMRAI::tbox::Pointer< LData > > & | F_data, | ||
std::vector< SAMRAI::tbox::Pointer< LData > > & | X_data, | ||
std::vector< SAMRAI::tbox::Pointer< LData > > & | ds_data, | ||
const std::string & | spread_kernel_fcn, | ||
RobinPhysBdryPatchStrategy * | f_phys_bdry_op, | ||
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > & | f_prolongation_scheds = std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >() , |
||
double | fill_data_time = 0.0 , |
||
bool | F_data_ghost_node_update = true , |
||
bool | X_data_ghost_node_update = true , |
||
bool | ds_data_ghost_node_update = true , |
||
int | coarsest_ln = invalid_level_number , |
||
int | finest_ln = invalid_level_number |
||
) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using a specified spreading kernel function.
f(i,j,k) = f(i,j,k) + Sum_{q,r,s} F(q,r,s) delta_h(x(i,j,k) - X(q,r,s)) ds(q,r,s)This is the standard regularized delta function spreading operation, which spreads densities, NOT values.
void IBTK::LDataManager::spread | ( | int | f_data_idx, |
std::vector< SAMRAI::tbox::Pointer< LData > > & | F_data, | ||
std::vector< SAMRAI::tbox::Pointer< LData > > & | X_data, | ||
std::vector< SAMRAI::tbox::Pointer< LData > > & | ds_data, | ||
RobinPhysBdryPatchStrategy * | f_phys_bdry_op, | ||
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > & | f_prolongation_scheds = std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >() , |
||
double | fill_data_time = 0.0 , |
||
bool | F_data_ghost_node_update = true , |
||
bool | X_data_ghost_node_update = true , |
||
bool | ds_data_ghost_node_update = true , |
||
int | coarsest_ln = invalid_level_number , |
||
int | finest_ln = invalid_level_number |
||
) |
Spread a quantity from the Lagrangian mesh to the Eulerian grid using the default spreading kernel function.
f(i,j,k) = f(i,j,k) + Sum_{q,r,s} F(q,r,s) delta_h(x(i,j,k) - X(q,r,s)) ds(q,r,s)This is the standard regularized delta function spreading operation, which spreads densities, NOT values.
void IBTK::LDataManager::updateNodeCountData | ( | int | coarsest_ln = invalid_level_number , |
int | finest_ln = invalid_level_number |
||
) |
Update the count of nodes per cell.
This routine updates cell data that is maintained on the patch hierarchy to track the number of nodes in each cell of the AMR index space. The node count data is used to tag cells for refinement, and to specify non-uniform load balancing.
|
static |
The name of the LData that specifies the initial positions of the curvilinear mesh nodes.
|
static |
The name of the LData that specifies the current positions of the curvilinear mesh nodes.
|
static |
The name of the LData that specifies the velocities of the curvilinear mesh nodes.