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

Class AdvDiffConservativeMassScalarTransportRKIntegrator is a concrete class that integrates the collocated density field. More...

#include </home/runner/work/IBAMR/IBAMR/include/ibamr/AdvDiffConservativeMassScalarTransportRKIntegrator.h>

Inheritance diagram for IBAMR::AdvDiffConservativeMassScalarTransportRKIntegrator:
Inheritance graph
[legend]

Public Member Functions

 AdvDiffConservativeMassScalarTransportRKIntegrator (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db)
 Class constructor.
 
 ~AdvDiffConservativeMassScalarTransportRKIntegrator ()
 Destructor.
 
void integrate (double dt) override
 Integrate density and momentum field.
 
- Public Member Functions inherited from IBAMR::STSMassFluxIntegrator
 STSMassFluxIntegrator (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db)
 Class constructor.
 
virtual ~STSMassFluxIntegrator ()=default
 Destructor.
 
void setDensityPatchDataIndex (int rho_idx)
 Set the current cell-centered density patch data index.
 
void setConvectiveDerivativePatchDataIndex (int N_idx)
 Set the new convective derivative patch data index.
 
void setDensityBoundaryConditions (const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &rho_sc_bc_coefs)
 Set the boundary condition object for the density.
 
int getUpdatedDensityPatchDataIndex ()
 Get the newly constructed density patch data index.
 
void setFluidVelocityPatchDataIndices (int V_old_idx, int V_current_idx, int V_new_idx)
 Set the patch data indices corresponding to the velocity at the previous time step to be used when computing the density update. More...
 
void setCycleNumber (int cycle_num)
 Set the cycle number currently being executed by the INS integrator. This will determine the rho advection velocity.
 
void setSolutionTime (double solution_time)
 Set solution time.
 
void setTimeInterval (double current_time, double new_time)
 Set the current time interval.
 
std::pair< double, double > getTimeInterval () const
 Get the current time interval.
 
double getTimeStepSize () const
 Get the current time step size.
 
void setPreviousTimeStepSize (double dt_prev)
 Set the previous time step value between times n - 1 (old) and n (current). This is used during velocity extrapolation.
 
void setHierarchyMathOps (SAMRAI::tbox::Pointer< IBTK::HierarchyMathOps > hier_math_ops)
 Set the HierarchyMathOps object used by the operator.
 
SAMRAI::tbox::Pointer< IBTK::HierarchyMathOpsgetHierarchyMathOps () const
 Get the HierarchyMathOps object used by the operator.
 

General operator functionality.

void initializeSTSIntegrator (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > base_hierarchy) override
 Compute hierarchy dependent data required for time integrating variables.
 
void deallocateSTSIntegrator () override
 Remove all hierarchy dependent data allocated by initializeTimeIntegrator(). More...
 
void setCellCenteredMaterialPropertyPatchDataIndex (int gamma_cc_idx)
 Set the current cell-centered material property patch data index.
 
void setCellCenteredTransportQuantityPatchDataIndex (int Q_cc_idx)
 Set the current cell-centered transport quantity patch data index.
 
void setCellCenteredDensityBoundaryConditions (SAMRAI::solv::RobinBcCoefStrategy< NDIM > *&rho_cc_bc_coefs)
 Set the boundary condition object for the cell-centered density.
 
void setCellCenteredMaterialPropertyBoundaryConditions (SAMRAI::solv::RobinBcCoefStrategy< NDIM > *&gamma_cc_bc_coefs)
 Set the boundary condition object for the cell-centered material property.
 
void setCellCenteredTransportQuantityBoundaryConditions (SAMRAI::solv::RobinBcCoefStrategy< NDIM > *&Q_cc_bc_coefs)
 Set the boundary condition object for the cell-centered transport quantity.
 
int getUpdatedCellCenteredDensityPatchDataIndex ()
 Get the newly constructed cell-centered density patch data index.
 
void setMaterialPropertyPatchDataIndices (int gamma_current_idx, int gamma_new_idx)
 Set the patch data indices corresponding to the material property to be used when computing the convective derivative. More...
 
void setTransportQuantityPatchDataIndices (int Q_current_idx, int Q_new_idx)
 Set the patch data indices corresponding to the transport quantity to be used when computing the convective derivative. More...
 
void setMaterialPropertyVariable (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > gamma_var)
 Set the material property variable.
 

Additional Inherited Members

- Protected Attributes inherited from IBAMR::STSMassFluxIntegrator
std::string d_object_name
 
SAMRAI::tbox::Pointer< StaggeredStokesPhysicalBoundaryHelperd_bc_helper
 
std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > d_bc_coefs
 
std::string d_density_bdry_extrap_type = "CONSTANT"
 
std::vector< IBTK::HierarchyGhostCellInterpolation::InterpolationTransactionComponentd_rho_transaction_comps
 
SAMRAI::tbox::Pointer< IBTK::HierarchyGhostCellInterpolationd_hier_rho_bdry_fill
 
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > d_hierarchy
 
int d_coarsest_ln = IBTK::invalid_level_number
 
int d_finest_ln = IBTK::invalid_level_number
 
std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > d_u_bc_coefs
 
std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > d_rho_bc_coefs
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyFaceDataOpsReal< NDIM, double > > d_hier_fc_data_ops
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchySideDataOpsReal< NDIM, double > > d_hier_sc_data_ops
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyCellDataOpsReal< NDIM, double > > d_hier_cc_data_ops
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > d_V_var
 
int d_V_scratch_idx = IBTK::invalid_index
 
int d_V_old_idx = IBTK::invalid_index
 
int d_V_current_idx = IBTK::invalid_index
 
int d_V_new_idx = IBTK::invalid_index
 
int d_V_composite_idx = IBTK::invalid_index
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > d_rho_var
 
int d_rho_current_idx = IBTK::invalid_index
 
int d_rho_scratch_idx = IBTK::invalid_index
 
int d_rho_new_idx = IBTK::invalid_index
 
int d_rho_composite_idx = IBTK::invalid_index
 
int d_N_idx = IBTK::invalid_index
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > d_S_var
 
int d_S_scratch_idx = IBTK::invalid_index
 
SAMRAI::tbox::Pointer< IBTK::CartGridFunctiond_S_fcn
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > d_E_var
 
int d_E_scratch_idx = IBTK::invalid_index
 
SAMRAI::tbox::Pointer< IBTK::HierarchyMathOpsd_hier_math_ops
 
bool d_hier_math_ops_external = false
 
bool d_is_initialized = false
 
bool d_enable_logging = false
 
LimiterType d_density_convective_limiter = CUI
 
int d_density_limiter_gcw = 1
 
TimeSteppingType d_density_time_stepping_type = FORWARD_EULER
 
int d_cycle_num = -1
 
double d_current_time = std::numeric_limits<double>::quiet_NaN()
 
double d_new_time = std::numeric_limits<double>::quiet_NaN()
 
double d_solution_time = std::numeric_limits<double>::quiet_NaN()
 
double d_dt = std::numeric_limits<double>::quiet_NaN()
 
double d_dt_prev = -1.0
 
std::vector< SAMRAI::hier::CoarseFineBoundary< NDIM > > d_cf_boundary
 

Detailed Description

Class AdvDiffConservativeMassScalarTransportRKIntegrator is a concrete class that integrates the collocated density field.

$ \frac{\partial \rho}{\partial t} + \nabla \cdot (\rho u) = S(x,t) $

and computes the conservative form of the convective operator $ \nabla \cdot (\rho \gamma u Q)$ where $\gamma$ is a material property and $Q$ is a transport quantity that is being transported. If $\gamma$ is not registered to this class, then this class computes $ \nabla \cdot (\rho u Q)$.

This class implements the Forward Euler (RK-1) for single cycle and midpoint rule (RK-2) for multiple cycles as the time-stepping schemes.

This class computes the convective derivative of a cell-centered transport quantity using various bounded-limiters described by Patel and Natarajan.

References Patel, JK. and Natarajan, G., A generic framework for design of interface capturing schemes for multi-fluid flows

A cell-centered density update is provided by this class, which is used in the conservative discretization of the energy equation.

Note
This class is specialized in that it computes a conservative discretization of the form $N = \nabla \cdot (u \rho \gamma Q)$, where the density $\rho$ can vary in space and time. This operator is to be used in conjuction with the conservative form of the variable coefficient energy equations.

Member Function Documentation

◆ deallocateSTSIntegrator()

void IBAMR::AdvDiffConservativeMassScalarTransportRKIntegrator::deallocateSTSIntegrator ( )
overridevirtual

Remove all hierarchy dependent data allocated by initializeTimeIntegrator().

Note
It is safe to call deallocateTimeIntegrator() when the time integrator is already deallocated.
See also
initializeTimeIntegrator

Implements IBAMR::STSMassFluxIntegrator.

◆ setMaterialPropertyPatchDataIndices()

void IBAMR::AdvDiffConservativeMassScalarTransportRKIntegrator::setMaterialPropertyPatchDataIndices ( int  gamma_current_idx,
int  gamma_new_idx 
)

Set the patch data indices corresponding to the material property to be used when computing the convective derivative.

Note
This material property will be used to compute an approximation to material property required for computing convective derivative. gamma_current_idx = n, gamma_new_idx = n+1,k (after a cycle). If gamma_new_idx is not set, then it will degenerate to gamma_current_idx automatically, for the very first simulation time step and cases where an AdvDiff cycle has not been executed, respectively.

◆ setTransportQuantityPatchDataIndices()

void IBAMR::AdvDiffConservativeMassScalarTransportRKIntegrator::setTransportQuantityPatchDataIndices ( int  Q_current_idx,
int  Q_new_idx 
)

Set the patch data indices corresponding to the transport quantity to be used when computing the convective derivative.

Note
This transport quantity will be used to compute an approximation to transport quantity required for computing convective derivative. Q_current_idx = n, Q_new_idx = n+1,k (after a cycle). If Q_new_idx is not set, then it will degenerate to Q_current_idx automatically, for the very first simulation time step and cases where an AdvDiff cycle has not been executed, respectively.

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