|
IBAMR
IBAMR version 0.19.
|
Class LDataManager coordinates the irregular distribution of LNode and LData on the patch hierarchy. More...
#include <ibtk/LDataManager.h>

Public Member Functions | |
| virtual double | getLevelDt (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const double dt_time, const bool initial_time) |
| virtual double | advanceLevel (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const tbox::Pointer< hier::BasePatchHierarchy< DIM > > 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< DIM > > level, const double new_time, const bool can_be_refined) |
| virtual void | resetDataToPreadvanceState (const tbox::Pointer< hier::BasePatchLevel< DIM > > level) |
| virtual void | initializeLevelData (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int level_number, const double init_data_time, const bool can_be_refined, const bool initial_time, const tbox::Pointer< hier::BasePatchLevel< DIM > > old_level=tbox::Pointer< hier::BasePatchLevel< DIM > >(NULL), const bool allocate_data=true)=0 |
| virtual void | resetHierarchyConfiguration (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int coarsest_level, const int finest_level)=0 |
| virtual void | applyGradientDetector (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int level_number, const double error_data_time, const int tag_index, const bool initial_time, const bool uses_richardson_extrapolation_too) |
| virtual void | applyRichardsonExtrapolation (const tbox::Pointer< hier::PatchLevel< DIM > > 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< DIM > > hierarchy, const int level_number, const tbox::Pointer< hier::PatchLevel< DIM > > coarser_level, const double coarsen_data_time, const bool before_advance) |
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 |
| static const std::string | INIT_POSN_DATA_NAME |
| static const std::string | VEL_DATA_NAME |
Methods to set and get the patch hierarchy and range of patch | |
| static std::map< std::string, LDataManager * > | s_data_manager_instances |
| static bool | s_registered_callback |
| static unsigned char | s_shutdown_priority |
| std::string | d_object_name |
| bool | d_registered_for_restart |
| SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > | d_hierarchy |
| SAMRAI::tbox::Pointer< SAMRAI::geom::CartesianGridGeometry< NDIM > > | d_grid_geom |
| int | d_coarsest_ln = IBTK::invalid_level_number |
| int | d_finest_ln = IBTK::invalid_level_number |
| SAMRAIDataCache | d_cached_eulerian_data |
| SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > | d_visit_writer |
| SAMRAI::tbox::Pointer< LSiloDataWriter > | d_silo_writer |
| SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > | d_load_balancer |
| SAMRAI::tbox::Pointer< LInitStrategy > | d_lag_init |
| std::vector< bool > | d_level_contains_lag_data |
| SAMRAI::tbox::Pointer< LNodeSetVariable > | d_lag_node_index_var |
| int | d_lag_node_index_current_idx = IBTK::invalid_index |
| int | d_lag_node_index_scratch_idx = IBTK::invalid_index |
| std::vector< SAMRAI::tbox::Pointer< std::vector< LNode > > > | d_local_and_ghost_nodes |
| double | d_beta_work = 1.0 |
| SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > | d_workload_var |
| int | d_workload_idx = IBTK::invalid_index |
| bool | d_output_workload = false |
| SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > | d_node_count_var |
| int | d_node_count_idx = IBTK::invalid_index |
| bool | d_output_node_count = false |
| const std::string | d_default_interp_kernel_fcn |
| const std::string | d_default_spread_kernel_fcn |
| bool | d_error_if_points_leave_domain |
| const SAMRAI::hier::IntVector< NDIM > | d_ghost_width |
| SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > > | d_lag_node_index_bdry_fill_alg |
| std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > | d_lag_node_index_bdry_fill_scheds |
| SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenAlgorithm< NDIM > > | d_node_count_coarsen_alg |
| std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > | d_node_count_coarsen_scheds |
| SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > | d_current_context |
| SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > | d_scratch_context |
| SAMRAI::hier::ComponentSelector | d_current_data |
| SAMRAI::hier::ComponentSelector | d_scratch_data |
| void | setPatchHierarchy (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy) |
| Reset patch hierarchy over which operations occur. More... | |
| SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > | getPatchHierarchy () const |
| Get the patch hierarchy used by this object. More... | |
| 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. More... | |
| const std::string & | getDefaultInterpKernelFunction () const |
| Return the default kernel function associated with the Eulerian-to-Lagrangian interpolation scheme. More... | |
| const std::string & | getDefaultSpreadKernelFunction () const |
| Return the default kernel function associated with the Lagrangian-to-Eulerian spreading scheme. 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, 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. More... | |
| 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. More... | |
| 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. More... | |
| 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. More... | |
| 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. More... | |
| void | registerLSiloDataWriter (SAMRAI::tbox::Pointer< LSiloDataWriter > silo_writer) |
| Register a Silo data writer with the manager. More... | |
| 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. More... | |
| 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. More... | |
| 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. More... | |
| 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. More... | |
| 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. More... | |
| std::vector< std::string > | getLagrangianStructureNames (int level_number) const |
| Get a list of Lagrangian structure names for the specified level of the patch hierarchy. More... | |
| std::vector< int > | getLagrangianStructureIDs (int level_number) const |
| Get a list of Lagrangian structure IDs for the specified level of the patch hierarchy. More... | |
| 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. More... | |
| 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. More... | |
| void | mapLagrangianToPETSc (std::vector< int > &inds, int level_number) const |
| Map the collection of Lagrangian indices to the corresponding global PETSc indices. More... | |
| void | mapPETScToLagrangian (std::vector< int > &inds, int level_number) const |
| Map the collection of global PETSc indices to the corresponding Lagrangian indices. More... | |
| AO & | getAO (int level_number) |
| Return the AO object that maps Lagrangian indices to PETSc indices. More... | |
| 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. More... | |
| ~LDataManager () | |
| The LDataManager destructor cleans up any remaining PETSc AO objects. More... | |
| LDataManager () | |
| Default constructor. More... | |
| LDataManager (const LDataManager &from)=delete | |
| Copy constructor. More... | |
| LDataManager & | operator= (const LDataManager &that)=delete |
| Assignment operator. More... | |
| void | scatterData (Vec &lagrangian_vec, Vec &petsc_vec, int level_number, ScatterMode mode) const |
| Common implementation of scatterPETScToLagrangian() and scatterLagrangianToPETSc(). More... | |
| void | beginNonlocalDataFill (int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
| Begin the process of refilling nonlocal Lagrangian quantities over the specified range of levels in the patch hierarchy. More... | |
| void | endNonlocalDataFill (int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) |
| End the process of refilling nonlocal Lagrangian quantities over the specified range of levels in the patch hierarchy. More... | |
| void | computeNodeDistribution (AO &ao, std::vector< int > &local_lag_indices, std::vector< int > &nonlocal_lag_indices, std::vector< int > &local_petsc_indices, std::vector< int > &nonlocal_petsc_indices, unsigned int &num_nodes, unsigned int &node_offset, int level_number) |
| void | getFromRestart () |
| static void | computeNodeOffsets (unsigned int &num_nodes, unsigned int &node_offset, unsigned int num_local_nodes) |
Data that is separately maintained for each level of the patch | |
| static std::vector< int > | s_ao_dummy |
| std::vector< std::map< std::string, int > > | d_strct_name_to_strct_id_map |
| std::vector< std::map< int, std::string > > | d_strct_id_to_strct_name_map |
| std::vector< std::map< int, std::pair< int, int > > > | d_strct_id_to_lag_idx_range_map |
| std::vector< std::map< int, int > > | d_last_lag_idx_to_strct_id_map |
| std::vector< ParallelSet > | d_inactive_strcts |
| std::vector< std::vector< int > > | d_displaced_strct_ids |
| std::vector< std::vector< std::pair< Point, Point > > > | d_displaced_strct_bounding_boxes |
| std::vector< std::vector< LNodeSet::value_type > > | d_displaced_strct_lnode_idxs |
| std::vector< std::vector< Point > > | d_displaced_strct_lnode_posns |
| std::vector< SAMRAI::tbox::Pointer< LMesh > > | d_lag_mesh |
| std::vector< std::map< std::string, SAMRAI::tbox::Pointer< LData > > > | d_lag_mesh_data |
| std::vector< bool > | d_needs_synch |
| std::vector< AO > | d_ao |
| std::vector< unsigned int > | d_num_nodes |
| std::vector< unsigned int > | d_node_offset |
| std::vector< std::vector< int > > | d_local_lag_indices |
| std::vector< std::vector< int > > | d_nonlocal_lag_indices |
| std::vector< std::vector< int > > | d_local_petsc_indices |
| std::vector< std::vector< int > > | d_nonlocal_petsc_indices |
| std::map< std::string, int > | d_user_defined_ldata |
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.
|
protected |
|
protected |
|
private |
|
privatedelete |
| from | The value to copy to this object. |
|
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.
|
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::setPatchHierarchy | ( | SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > | hierarchy | ) |
| SAMRAI::tbox::Pointer<SAMRAI::hier::PatchHierarchy<NDIM> > IBTK::LDataManager::getPatchHierarchy | ( | ) | const |
The levels must exist in the hierarchy or an assertion failure will result.
|
inline |
|
inline |
|
inline |
| 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 |
||
| ) |
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, | ||
| 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 |
||
| ) |
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 |
||
| ) |
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, | ||
| 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 |
||
| ) |
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, | ||
| 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 |
||
| ) |
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, | ||
| 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 |
||
| ) |
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 |
||
| ) |
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, | ||
| 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 |
||
| ) |
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::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 |
||
| ) |
| void IBTK::LDataManager::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 |
||
| ) |
| void IBTK::LDataManager::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 |
||
| ) |
| void IBTK::LDataManager::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 |
||
| ) |
| 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::freeLInitStrategy | ( | ) |
Free the concrete initialization strategy object.
| void IBTK::LDataManager::registerVisItDataWriter | ( | SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > | visit_writer | ) |
| void IBTK::LDataManager::registerLSiloDataWriter | ( | SAMRAI::tbox::Pointer< LSiloDataWriter > | silo_writer | ) |
| void IBTK::LDataManager::registerLoadBalancer | ( | SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > | load_balancer, |
| int | workload_data_idx | ||
| ) |
|
inline |
|
inline |
| SAMRAI::tbox::Pointer<LData> IBTK::LDataManager::createLData | ( | const std::string & | quantity_name, |
| int | level_number, | ||
| unsigned int | depth = 1, |
||
| bool | maintain_data = false |
||
| ) |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
| Point IBTK::LDataManager::computeLagrangianStructureCenterOfMass | ( | int | structure_id, |
| int | level_number | ||
| ) |
X = (1/N) Sum_{k in structure} X_k
in which N is the number of nodes associated with that structure.
| std::pair<Point, Point> IBTK::LDataManager::computeLagrangianStructureBoundingBox | ( | int | structure_id, |
| int | level_number | ||
| ) |
| void IBTK::LDataManager::reinitLagrangianStructure | ( | const Point & | X_center, |
| int | structure_id, | ||
| int | level_number | ||
| ) |
| void IBTK::LDataManager::displaceLagrangianStructure | ( | const Vector & | dX, |
| int | structure_id, | ||
| int | level_number | ||
| ) |
| void IBTK::LDataManager::activateLagrangianStructures | ( | const std::vector< int > & | structure_ids, |
| int | level_number | ||
| ) |
| void IBTK::LDataManager::inactivateLagrangianStructures | ( | const std::vector< int > & | structure_ids, |
| int | level_number | ||
| ) |
|
inline |
| void IBTK::LDataManager::zeroInactivatedComponents | ( | SAMRAI::tbox::Pointer< LData > | lag_data, |
| int | level_number | ||
| ) | const |
| AO& IBTK::LDataManager::getAO | ( | int | level_number | ) |
| void IBTK::LDataManager::scatterLagrangianToPETSc | ( | Vec & | lagrangian_vec, |
| Vec & | petsc_vec, | ||
| int | level_number | ||
| ) | const |
| void IBTK::LDataManager::scatterPETScToLagrangian | ( | Vec & | petsc_vec, |
| Vec & | lagrangian_vec, | ||
| int | level_number | ||
| ) | const |
| void IBTK::LDataManager::scatterToAll | ( | Vec & | parallel_vec, |
| Vec & | sequential_vec | ||
| ) | const |
| void IBTK::LDataManager::scatterToZero | ( | Vec & | parallel_vec, |
| Vec & | sequential_vec | ||
| ) | const |
| void IBTK::LDataManager::beginDataRedistribution | ( | int | coarsest_ln = invalid_level_number, |
| int | finest_ln = invalid_level_number |
||
| ) |
This method uses the present location of each Lagrangian mesh node to redistribute the LNodeData managed by this object.
| void IBTK::LDataManager::endDataRedistribution | ( | int | coarsest_ln = invalid_level_number, |
| int | finest_ln = invalid_level_number |
||
| ) |
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.
| 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 |
||
| ) |
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.
| void IBTK::LDataManager::updateNodeCountData | ( | int | coarsest_ln = invalid_level_number, |
| int | finest_ln = invalid_level_number |
||
| ) |
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.
|
override |
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).
|
override |
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.
|
override |
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.
|
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::registerUserDefinedLData | ( | const std::string & | data_name, |
| int | depth | ||
| ) |
Register user defined Lagrangian data to be maintained
|
privatedelete |
| that | The value to assign to this object. |
|
private |
|
private |
The operation is essentially equivalent to refilling ghost cells for structured (SAMRAI native) data.
|
private |
The operation is essentially equivalent to refilling ghost cells for structured (SAMRAI native) data.
|
private |
Determines the global Lagrangian and PETSc indices of the local and nonlocal nodes associated with the processor as well as the local PETSc indices of the interior and ghost nodes in each patch of the specified level.
Since each processor may own multiple patches in a given level, nodes appearing in the ghost cell region of a patch may or may not be owned by this processor.
|
staticprivate |
Determine the number of local Lagrangian nodes on all MPI processes with rank less than the rank of the current MPI process.
|
private |
Read object state from the restart file and initialize class data members. The database from which the restart data is read is determined by the object_name specified in the constructor.
Unrecoverable Errors:
|
virtualinherited |
Determine time increment to advance data on level. The recompute_dt option specifies whether to compute the timestep using the current level data or to return the value stored by the time integrator. The default true setting means the timestep will be computed if no value is supplied.
This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
virtualinherited |
Advance data on all patches on specified patch level from current time (current_time) to new time (new_time). This routine is called only during time-dependent regridding procedures, such as Richardson extrapolation. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed. The boolean arguments are used to determine the state of the algorithm and the data when the advance routine is called. Note that this advance function is also used during normal time integration steps.
When this function is called, the level data required to begin the advance must be allocated and be defined appropriately. Typically, this is equivalent to what is needed to initialize a new level after regridding. Upon exiting this routine, both current and new data may exist on the level. This data is needed until level synchronization occurs, in general. Current and new data may be reset by calling the member function resetTimeDependentData().
This routine is called from two different points within the Richardson exptrapolation process: to advance a temporary level that is coarser than the hierarchy level on which error estimation is performed, and to advance the hierarchy level itself. In the first case, the values of the boolean flags are:
In the second case, the values of the boolean flags are:
|
virtualinherited |
Reset time-dependent data storage for the specified patch level.
This routine only applies when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
virtualinherited |
Reset data on the patch level by destroying all patch data other than that which is needed to initialize the solution on that level. In other words, this is the data needed to begin a time integration step on the level.
This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
pure virtualinherited |
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.
Generally, when data is set, it is interpolated from coarser levels in the hierarchy. If the old level pointer in the argument list is non-null, then data is copied from the old level to the new level on regions of intersection between those levels before interpolation occurs. In this case, the level number must match that of the old level. The specific operations that occur when initializing level data are determined by the particular solution methods in use; i.e., in the subclass of this abstract base class.
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 may be different in each of those circumstances. The can_be_refined boolean argument indicates whether the level is the finest allowable level in the hierarchy.
|
pure virtualinherited |
After hierarchy levels have changed and data has been initialized on the new levels, this routine can be used to reset any information needed by the solution method that is particular to the hierarchy configuration. For example, the solution procedure may cache communication schedules to amortize the cost of data movement on the AMR patch hierarchy. This function will be called by the gridding algorithm after the initialization occurs so that the algorithm-specific subclass can reset such things. Also, if the solution method must make the solution consistent across multiple levels after the hierarchy is changed, this process may be invoked by this routine. Of course the details of these processes are determined by the particular solution methods in use.
The level number arguments indicate the coarsest and finest levels in the current hierarchy configuration that have changed. It should be assumed that all intermediate levels have changed as well.
|
virtualinherited |
Set integer tags to "one" in cells where refinement of the given level should occur according to some user-supplied gradient criteria. 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. Typically, this information is passed to the user's patch tagging routines since the error estimator or gradient detector may be different in each case.
The boolean 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.
This routine is only when gradient detector is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
virtualinherited |
Set integer tags to "one" in cells where refinement of the given level should occur according to some user-supplied Richardson extrapolation criteria. The "error_data_time" argument is the regrid time. The "deltat" argument is the time increment to advance the solution on the level to be refined. Note that that level is finer than the level in the argument list, in general. The ratio between the argument level and the actual hierarchy level is given by the integer "coarsen ratio".
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. Typically, this information is passed to the user's patch tagging routines since the application of the Richardson extrapolation process may be different in each case.
The boolean uses_gradient_detector_too is true when a gradient detector procedure is used in addition to Richardson extrapolation, and false otherwise. This argument helps the user to manage multiple regridding criteria.
This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
virtualinherited |
Coarsen solution data from level to coarse_level for Richardson extrapolation. Note that this routine will be called twice during the Richardson extrapolation error estimation process, once to set data on the coarser level and once to coarsen data from after advancing the fine level. The init_coarse_level boolean argument indicates whether data is set on the coarse level by coarsening the "old" time level solution or by coarsening the "new" solution on the fine level (i.e., after it has been advanced).
This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
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 initial positions of the curvilinear mesh nodes.
|
static |
The name of the LData that specifies the velocities of the curvilinear mesh nodes.
|
staticprivate |
Static data members used to control access to and destruction of singleton data manager instance.
|
staticprivate |
|
staticprivate |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
Information about the names and IDs of the various Lagrangian structures.
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
Lagrangian mesh data.
|
private |
The Lagrangian mesh data owned by the manager object.
|
private |
Indicates whether the LData is in synch with the LNodeData.
|
private |
PETSc AO objects provide mappings between the fixed global Lagrangian node IDs and the ever-changing global PETSc ordering.
|
staticprivate |
|
private |
The total number of nodes for all processors.
|
private |
The total number of local nodes for all processors with rank less than the rank of the current processor.
|
private |
The Lagrangian node indices of all local and nonlocal nodes on each level of the patch hierarchy.
A local node is one that is owned by a patch on this processor, while a nonlocal node is one that is owned by a patch on another processor, but found in the ghost region of some patch owned by this processor.
Note that these sets of indices provide the information necessary to determine the local PETSc index for all nodes. Local node d_local_lag_indices[ln][j] has local PETSc index j, while nonlocal node d_nonlocal_lag_indices[ln][k] has local PETSc index d_local_lag_indices.size()+j.
It is possible to determine the global PETSc index of a local node by making use of d_node_offset. Local node d_local_lag_indices[ln][j] has global PETSc index j+d_node_offset[ln]. A similar mapping for nonlocal nodes is not well defined.
|
private |
|
private |
The node indices of all local nodes (i.e. the nodes owned by this processor) on each level of the hierarchy. The indices are in the global PETSc ordering corresponding to a depth of 1.
|
private |
The node indices of all nonlocal nodes (i.e. the nodes owned by another processor that appear in the ghost region of some patch owned by this processor) on each level of the hierarchy. The indices are in the global PETSc ordering corresponding to a depth of 1.
|
private |
Container for additional user defined Lagrangian data
1.8.17