IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
Static Public Member Functions | Static Public Attributes | List of all members
IBTK::LDataManager Class Reference

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>

Inheritance diagram for IBTK::LDataManager:
Inheritance graph
[legend]

Static Public Member Functions

static LDataManagergetManager (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

levels associated with this manager class.

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::stringgetDefaultInterpKernelFunction () const
 Return the default kernel function associated with the Eulerian-to-Lagrangian interpolation scheme.
 
const std::stringgetDefaultSpreadKernelFunction () 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< LMeshgetLMesh (int level_number) const
 Get the Lagrangian mesh associated with the given patch hierarchy level.
 
SAMRAI::tbox::Pointer< LDatagetLData (const std::string &quantity_name, int level_number) const
 Get the specified Lagrangian quantity data on the given patch hierarchy level.
 
SAMRAI::tbox::Pointer< LDatacreateLData (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::stringgetLagrangianStructureNames (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 &parallel_vec, Vec &sequential_vec) const
 Scatter data from a distributed PETSc vector to all processors. More...
 
void scatterToZero (Vec &parallel_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 > >(NULL), 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)
 

Detailed Description

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.

void
update_target_points(Pointer<PatchHierarchy<NDIM>> hierarchy,
LDataManager* const l_data_manager,
const double current_time,
const double dt)
{
const int finest_ln = hierarchy->getFinestLevelNumber();
// Note we assume the circle is the 0th structure and on the finest level.
const std::pair<int, int>& circle_lag_idxs = l_data_manager->getLagrangianStructureIndexRange(0, finest_ln);
Pointer<LMesh> mesh = l_data_manager->getLMesh(finest_ln);
const std::vector<LNode*>& local_nodes = mesh->getLocalNodes();
const std::vector<LNode*>& ghost_nodes = mesh->getGhostNodes();
std::vector<LNode*> nodes = local_nodes;
nodes.insert(nodes.end(), ghost_nodes.begin(), ghost_nodes.end());
for (auto& node : nodes)
{
auto force_spec = node->getNodeDataItem<IBTargetPointForceSpec>();
if (force_spec == nullptr) continue;
const int lag_idx = node->getLagrangianIndex();
Point& X_target = force_spec->getTargetPointPosition();
if (circle_lag_idxs.first <= lag_idx && lag_idx < circle_lag_idxs.second)
{
X_target[0] += 0.1 * dt; // Move point to the right with speed 0.1
}
}
return;
}
constexpr const_iterator begin() const noexcept
constexpr const_iterator end() const noexcept
constexpr iterator insert(const_iterator __position, _InputIterator __first, _InputIterator __last)

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
print_eul_data(Pointer<PatchHierarchy<NDIM>> hierarchy,
LDataManager* const l_data_manager)
{
const int finest_ln = hierarchy->getFinestLevelNumber();
Pointer<LData> X_data = l_data_manager->getLData("X", finest_ln);
Vec X_vec = X_data->getVec();
double* x_vals;
int ierr = VecGetArray(X_vec, &x_vals);
IBTK_CHKERRQ(ierr);
Pointer<LMesh> mesh = l_data_manager->getLMesh(finest_ln);
const std::vector<LNode*>& nodes = mesh->getLocalNodes();
for (const auto& node : nodes)
{
const int lag_idx = node->getLagrangianIndex();
// Note we use the local index instead of the global index.
// The global index is the local index plust the total indices on all processors of lower rank.
const int petsc_idx = node->getLocalPETScIndex();
Eigen::Map<VectorNd> X(&x_vals[petsc_idx * NDIM]);
pout << "Euerian location of node " << lag_idx << ":\n"
<< X << "\n";
}
return;
}
std::ostream pout
Note
Multiple LDataManager objects may be instantiated simultaneously.
See also
LMesh, LNode, LData

Member Function Documentation

◆ activateLagrangianStructures()

void IBTK::LDataManager::activateLagrangianStructures ( const std::vector< int > &  structure_ids,
int  level_number 
)

Activate the Lagrangian structures with the specified ID numbers.

Note
This method is collective (i.e., must be called by all MPI processes); however, each MPI process may provide a different collection of structures to activate.

◆ addWorkloadEstimate()

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.

◆ applyGradientDetector()

void IBTK::LDataManager::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 
)
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 >.

◆ beginDataRedistribution()

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.

Note
This routine assumes that the time interval between node redistribution satisfies a timestep restriction of the form dt <= C*dx*|U| with C <= 1. This restriction prevents nodes from moving more than one cell width per timestep.
See also
endDataRedistribution

◆ computeLagrangianStructureBoundingBox()

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.

Note
Returns the entire range of double precision values in the case that the Lagrangian structure ID is not associated with any Lagrangian structure.

◆ computeLagrangianStructureCenterOfMass()

Point IBTK::LDataManager::computeLagrangianStructureCenterOfMass ( int  structure_id,
int  level_number 
)

Get the center of mass of the Lagrangian structure with the specified ID.

Note
The center of mass X of a particular structure is computed as

X = (1/N) Sum_{k in structure} X_k

in which N is the number of nodes associated with that structure.

Note
Returns Point::Zero() in the case that the Lagrangian structure ID is not associated with any Lagrangian structure.

◆ createLData()

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.

Note
Quantities maintained by the LDataManager must have unique names. The name "X" is reserved for the nodal coordinates.

◆ displaceLagrangianStructure()

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.

Note
This operation must be performed immediately before a regridding operation, otherwise the results are undefined.
Warning
All displacements must involve shifts that do not cross periodic boundaries.

◆ endDataRedistribution()

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.

Note
Since this routine potentially results in a large amount of inter-processor communication, it may be worth putting it off for as long as possible. If the timestep dt satisfies a condition of the form dt <= C*dx*|U| with C << 1, it may be possible to redistribute the Lagrangian data less frequently than every timestep.
See also
beginDataRedistribution

◆ freeAllManagers()

void IBTK::LDataManager::freeAllManagers ( )
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.

◆ freeLInitStrategy()

void IBTK::LDataManager::freeLInitStrategy ( )

Free the concrete initialization strategy object.

Note
Be sure to call this method only once the initialization object is no longer needed.

◆ getGlobalNodeOffset()

unsigned int IBTK::LDataManager::getGlobalNodeOffset ( int  level_number) const
inline
Returns
The number of nodes on all processors with MPI rank less than the current process on the specified level of the patch hierarchy.
Note
This count does not include nodes that only lie in ghost cells for the current process.

◆ getLagrangianStructureID() [1/2]

int IBTK::LDataManager::getLagrangianStructureID ( const std::string structure_name,
int  level_number 
) const
inline

Get the ID of the Lagrangian structure with the specified name.

Note
Returns -1 in the case that the Lagrangian structure name is not associated with any Lagrangian structure.

◆ getLagrangianStructureID() [2/2]

int IBTK::LDataManager::getLagrangianStructureID ( int  lagrangian_index,
int  level_number 
) const
inline

Get the ID of the Lagrangian structure associated with the specified Lagrangian index.

Note
Returns -1 in the case that the Lagrangian index is not associated with any Lagrangian structure.

◆ getLagrangianStructureIndexRange()

std::pair< int, int > IBTK::LDataManager::getLagrangianStructureIndexRange ( int  structure_id,
int  level_number 
) const
inline

Get the range of Lagrangian indices for the Lagrangian structure with the specified ID.

Returns
A pair of indices such that if pair.first <= lag_idx < pair.second, then lag_idx is associated with the specified structure; otherwise, lag_idx is not associated with the specified structure.
Note
Returns std::make_pair(-1,-1) in the case that the Lagrangian structure ID is not associated with any Lagrangian structure.

◆ getLagrangianStructureName()

std::string IBTK::LDataManager::getLagrangianStructureName ( int  structure_id,
int  level_number 
) const
inline

Get the name of the Lagrangian structure with the specified ID.

Note
Returns "UNKNOWN" in the case that the Lagrangian structure ID is not associated with any Lagrangian structure.

◆ getManager()

LDataManager * IBTK::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

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.

Returns
A pointer to the data manager instance.
Note
By default, the ghost cell width is set according to the interpolation and spreading kernel functions.

◆ getNumberOfGhostNodes()

unsigned int IBTK::LDataManager::getNumberOfGhostNodes ( int  level_number) const
inline
Returns
The number of ghost (i.e., off-processor) nodes of the Lagrangian data for the specified level of the patch hierarchy.
See also
getNumberOfNodes
getNumberOfLocalNodes

◆ getNumberOfLocalNodes()

unsigned int IBTK::LDataManager::getNumberOfLocalNodes ( int  level_number) const
inline
Returns
The number of local (i.e., on-processor) nodes of the Lagrangian data for the specified level of the patch hierarchy.
Note
This count does not include nodes that only lie in ghost cells for the current process.
See also
getNumberOfNodes
getNumberOfGhostNodes

◆ getNumberOfNodes()

unsigned int IBTK::LDataManager::getNumberOfNodes ( int  level_number) const
inline
Returns
The number of total nodes of the Lagrangian data for the specified level of the patch hierarchy.

◆ getPatchLevels()

std::pair< int, int > IBTK::LDataManager::getPatchLevels ( ) const

Get the range of patch levels used by this object.

Note
Returns [coarsest_ln,finest_ln+1).

◆ getWorkloadPatchDescriptorIndex()

int IBTK::LDataManager::getWorkloadPatchDescriptorIndex ( ) const
inline

Get the patch data descriptor index for the workload cell data.

Deprecated:
This method is deprecated since, in future versions of IBAMR, this value will no longer be stored and will only be available via the parent hierarchy integrator.

◆ inactivateLagrangianStructures()

void IBTK::LDataManager::inactivateLagrangianStructures ( const std::vector< int > &  structure_ids,
int  level_number 
)

Inactivate the Lagrangian structures with the specified ID numbers.

Note
This method is collective (i.e., must be called by all MPI processes); however, each MPI process may provide a different collection of structures to inactivate.

◆ initializeLevelData()

void IBTK::LDataManager::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> >(NULL),
bool  allocate_data = true 
)
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 >.

◆ putToDatabase()

void IBTK::LDataManager::putToDatabase ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db)
overridevirtual

Write out object state to the given database.

When assertion checking is active, database pointer must be non-null.

Implements SAMRAI::tbox::Serializable.

◆ registerLInitStrategy()

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.

Note
This function calls LInitStrategy::init(). All preprocessing should be completed before registering a LInitStrategy object.

◆ registerLoadBalancer()

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.

Deprecated:
This method is deprecated since the current strategy for handling non-uniform load balancing does not require that this object store a pointer to the load balancer.

◆ registerUserDefinedLData()

void IBTK::LDataManager::registerUserDefinedLData ( const std::string data_name,
int  depth 
)

Register user defined Lagrangian data to be maintained

◆ reinitLagrangianStructure()

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.

Note
This operation must be performed immediately before a regridding operation, otherwise the results are undefined.

◆ resetHierarchyConfiguration()

void IBTK::LDataManager::resetHierarchyConfiguration ( SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > >  hierarchy,
int  coarsest_ln,
int  finest_ln 
)
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 >.

◆ scatterLagrangianToPETSc()

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.

Todo:
Optimize the implementation of this method.

◆ scatterPETScToLagrangian()

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.

Todo:
Optimize the implementation of this method.

◆ scatterToAll()

void IBTK::LDataManager::scatterToAll ( Vec &  parallel_vec,
Vec &  sequential_vec 
) const

Scatter data from a distributed PETSc vector to all processors.

Todo:
Optimize the implementation of this method.

◆ scatterToZero()

void IBTK::LDataManager::scatterToZero ( Vec &  parallel_vec,
Vec &  sequential_vec 
) const

Scatter data from a distributed PETSc vector to processor zero.

Todo:
Optimize the implementation of this method.

◆ setPatchLevels()

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.

◆ spread() [1/8]

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.

Note
This spreading operation does NOT include the scale factor corresponding to the curvilinear volume element (dq dr ds). The spreading formula is
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.

◆ spread() [2/8]

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.

Note
This spreading operation does NOT include the scale factor corresponding to the curvilinear volume element (dq dr ds). The spreading formula is
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.

◆ spread() [3/8]

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.

Note
This spreading operation does include the scale factor corresponding to the curvilinear volume element (dq dr ds). The spreading formula is
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.

◆ spread() [4/8]

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.

Note
This spreading operation does include the scale factor corresponding to the curvilinear volume element (dq dr ds). The spreading formula is
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.

◆ spread() [5/8]

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.

Note
This spreading operation does NOT include the scale factor corresponding to the curvilinear volume element (dq dr ds). The spreading formula is
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.

◆ spread() [6/8]

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.

Note
This spreading operation does NOT include the scale factor corresponding to the curvilinear volume element (dq dr ds). The spreading formula is
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.

◆ spread() [7/8]

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.

Note
This spreading operation does include the scale factor corresponding to the curvilinear volume element (dq dr ds). The spreading formula is
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.

◆ spread() [8/8]

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.

Note
This spreading operation does include the scale factor corresponding to the curvilinear volume element (dq dr ds). The spreading formula is
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.

◆ updateNodeCountData()

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.

Member Data Documentation

◆ INIT_POSN_DATA_NAME

const std::string IBTK::LDataManager::INIT_POSN_DATA_NAME = "X0"
static

The name of the LData that specifies the initial positions of the curvilinear mesh nodes.

◆ POSN_DATA_NAME

const std::string IBTK::LDataManager::POSN_DATA_NAME = "X"
static

The name of the LData that specifies the current positions of the curvilinear mesh nodes.

◆ VEL_DATA_NAME

const std::string IBTK::LDataManager::VEL_DATA_NAME = "U"
static

The name of the LData that specifies the velocities of the curvilinear mesh nodes.


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