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

Class AdvectorPredictorCorrectorHyperbolicPatchOps is a concrete SAMRAI::algs::HyperbolicPatchStrategy that makes use of class AdvectorExplicitPredictorPatchOps to solve the linear advection equation. More...

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

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

Public Member Functions

 AdvectorPredictorCorrectorHyperbolicPatchOps (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, SAMRAI::tbox::Pointer< AdvectorExplicitPredictorPatchOps > explicit_predictor, SAMRAI::tbox::Pointer< SAMRAI::geom::CartesianGridGeometry< NDIM > > grid_geom, bool register_for_restart=true)
 
virtual ~AdvectorPredictorCorrectorHyperbolicPatchOps ()
 
const std::stringgetName () const
 
void registerVisItDataWriter (SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > visit_writer)
 
void registerAdvectionVelocity (SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > u_var)
 
void setAdvectionVelocityIsDivergenceFree (SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > u_var, bool is_div_free)
 
void setAdvectionVelocityFunction (SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > u_var, SAMRAI::tbox::Pointer< IBTK::CartGridFunction > u_fcn)
 
void registerSourceTerm (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > F_var)
 
void setSourceTermFunction (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > F_var, SAMRAI::tbox::Pointer< IBTK::CartGridFunction > F_fcn)
 
void registerTransportedQuantity (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var)
 
void setAdvectionVelocity (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > u_var)
 
void setSourceTerm (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > F_var)
 
void setConvectiveDifferencingType (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, ConvectiveDifferencingType difference_form)
 
void setInitialConditions (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, SAMRAI::tbox::Pointer< IBTK::CartGridFunction > Q_init)
 
void setPhysicalBcCoefs (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, SAMRAI::solv::RobinBcCoefStrategy< NDIM > *Q_bc_coef)
 
void setPhysicalBcCoefs (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > Q_bc_coef)
 
virtual void registerModelVariables (SAMRAI::algs::HyperbolicLevelIntegrator< NDIM > *integrator) override
 Register AdvectorPredictorCorrectorHyperbolicPatchOps model variables with the SAMRAI::algs::HyperbolicLevelIntegrator according to the variable registration function provided by the integrator. More...
 
virtual void initializeDataOnPatch (SAMRAI::hier::Patch< NDIM > &patch, double data_time, bool initial_time) override
 Set the data on the patch interior to some initial values via the concrete IBTK::CartGridFunction objects registered with the patch strategy when provided. Otherwise, initialize data to zero.
 
virtual double computeStableDtOnPatch (SAMRAI::hier::Patch< NDIM > &patch, bool initial_time, double dt_time) override
 Compute a stable time increment for patch using an explicit CFL condition and return the computed dt.
 
virtual void computeFluxesOnPatch (SAMRAI::hier::Patch< NDIM > &patch, double time, double dt) override
 Compute the time integral of the fluxes to be used in conservative difference for patch integration. More...
 
virtual void conservativeDifferenceOnPatch (SAMRAI::hier::Patch< NDIM > &patch, double time, double dt, bool at_synchronization) override
 Update solution variables by performing a conservative difference using the fluxes calculated by computeFluxesOnPatch().
 
virtual void preprocessAdvanceLevelState (const SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > &level, double current_time, double dt, bool first_step, bool last_step, bool regrid_advance) override
 Compute the values of any time-dependent source terms for use by the explicit predictor. More...
 
virtual void postprocessAdvanceLevelState (const SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > &level, double current_time, double dt, bool first_step, bool last_step, bool regrid_advance) override
 Add source terms to the updated solution. More...
 
virtual void tagGradientDetectorCells (SAMRAI::hier::Patch< NDIM > &patch, double regrid_time, bool initial_error, int tag_indexindx, bool uses_richardson_extrapolation_too) override
 Tag cells for refinement using a gradient detector.
 
virtual void setPhysicalBoundaryConditions (SAMRAI::hier::Patch< NDIM > &patch, double fill_time, const SAMRAI::hier::IntVector< NDIM > &ghost_width_to_fill) override
 Set the data in ghost cells corresponding to physical boundary conditions.
 
virtual void putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override
 Write state of AdvectorPredictorCorrectorHyperbolicPatchOps object to the given database for restart. More...
 
- Public Member Functions inherited from SAMRAI::algs::HyperbolicPatchStrategy< NDIM >
virtual void setupLoadBalancer (HyperbolicLevelIntegrator< NDIM > *integrator, mesh::GriddingAlgorithm< NDIM > *gridding_algorithm)
 
virtual void tagRichardsonExtrapolationCells (hier::Patch< NDIM > &patch, const int error_level_number, const tbox::Pointer< hier::VariableContext > coarsened_fine, const tbox::Pointer< hier::VariableContext > advanced_coarse, const double regrid_time, const double deltat, const int error_coarsen_ratio, const bool initial_error, const int tag_index, const bool uses_gradient_detector_too)
 
virtual hier::IntVector< NDIM > getRefineOpStencilWidth () const
 
virtual void preprocessRefine (hier::Patch< NDIM > &fine, const hier::Patch< NDIM > &coarse, const hier::Box< NDIM > &fine_box, const hier::IntVector< NDIM > &ratio)
 
virtual void postprocessRefine (hier::Patch< NDIM > &fine, const hier::Patch< NDIM > &coarse, const hier::Box< NDIM > &fine_box, const hier::IntVector< NDIM > &ratio)
 
virtual hier::IntVector< NDIM > getCoarsenOpStencilWidth () const
 
virtual void preprocessCoarsen (hier::Patch< NDIM > &coarse, const hier::Patch< NDIM > &fine, const hier::Box< NDIM > &coarse_box, const hier::IntVector< NDIM > &ratio)
 
virtual void postprocessCoarsen (hier::Patch< NDIM > &coarse, const hier::Patch< NDIM > &fine, const hier::Box< NDIM > &coarse_box, const hier::IntVector< NDIM > &ratio)
 
tbox::Pointer< hier::VariableContextgetDataContext () const
 
void setDataContext (tbox::Pointer< hier::VariableContext > context)
 
void clearDataContext ()
 
- Public Member Functions inherited from SAMRAI::xfer::RefinePatchStrategy< DIM >
virtual void setPhysicalBoundaryConditions (hier::Patch< DIM > &patch, const double fill_time, const hier::IntVector< DIM > &ghost_width_to_fill)=0
 
virtual void preprocessRefine (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::Box< DIM > &fine_box, const hier::IntVector< DIM > &ratio)=0
 
virtual void postprocessRefine (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::Box< DIM > &fine_box, const hier::IntVector< DIM > &ratio)=0
 
virtual void preprocessRefineBoxes (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::BoxList< DIM > &fine_boxes, const hier::IntVector< DIM > &ratio)
 
virtual void postprocessRefineBoxes (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::BoxList< DIM > &fine_boxes, const hier::IntVector< DIM > &ratio)
 
- Public Member Functions inherited from SAMRAI::xfer::CoarsenPatchStrategy< DIM >
virtual void preprocessCoarsen (hier::Patch< DIM > &coarse, const hier::Patch< DIM > &fine, const hier::Box< DIM > &coarse_box, const hier::IntVector< DIM > &ratio)=0
 
virtual void postprocessCoarsen (hier::Patch< DIM > &coarse, const hier::Patch< DIM > &fine, const hier::Box< DIM > &coarse_box, const hier::IntVector< DIM > &ratio)=0
 

Protected Member Functions

SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceData< NDIM, double > > getFluxIntegralData (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, SAMRAI::hier::Patch< NDIM > &patch, SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > context)
 Get a pointer to the requested flux integral patch data on the specified patch.
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceData< NDIM, double > > getQIntegralData (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, SAMRAI::hier::Patch< NDIM > &patch, SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > context)
 Get a pointer to the requested q integral patch data on the specified patch.
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceData< NDIM, double > > getUIntegralData (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, SAMRAI::hier::Patch< NDIM > &patch, SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > context)
 Get a pointer to the requested u integral patch data on the specified patch.
 

Protected Attributes

SAMRAI::algs::HyperbolicLevelIntegrator< NDIM > * d_integrator = nullptr
 
SAMRAI::tbox::Pointer< AdvectorExplicitPredictorPatchOpsd_explicit_predictor
 
std::set< SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > > d_u_var
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > >, bool > d_u_is_div_free
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > >, SAMRAI::tbox::Pointer< IBTK::CartGridFunction > > d_u_fcn
 
bool d_compute_init_velocity = true
 
bool d_compute_half_velocity = true
 
bool d_compute_final_velocity = true
 
std::set< SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > > d_F_var
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >, SAMRAI::tbox::Pointer< IBTK::CartGridFunction > > d_F_fcn
 
std::set< SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > > d_Q_var
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >, SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > > d_Q_u_map
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >, SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > > d_Q_F_map
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >, ConvectiveDifferencingTyped_Q_difference_form
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >, SAMRAI::tbox::Pointer< IBTK::CartGridFunction > > d_Q_init
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >, std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > > d_Q_bc_coef
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >, SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > > d_flux_integral_var
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >, SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > > d_q_integral_var
 
std::map< SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > >, SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > > d_u_integral_var
 
bool d_overwrite_tags = true
 

Detailed Description

Class AdvectorPredictorCorrectorHyperbolicPatchOps is a concrete SAMRAI::algs::HyperbolicPatchStrategy that makes use of class AdvectorExplicitPredictorPatchOps to solve the linear advection equation.

Class AdvectorPredictorCorrectorHyperbolicPatchOps provides numerical routines for solving the advection equation in conservative form,

\[ \frac{dQ}{dt} + \nabla \cdot (\vec{u}^{\mbox{\scriptsize ADV}} Q) = F - Q \nabla \cdot *\vec{u}^{\mbox{\scriptsize ADV}}, \]

where $ Q $ is a cell-centered scalar- or vector-valued quantity, $ \vec{u}^{\mbox{\scriptsize ADV}} $ is a specified face-centered advection velocity, and $ F $ is an optional source term. When $ \vec{u}^{\mbox{\scriptsize ADV}} $ is discretely divergence free, this is equivalent to solving

\[ \frac{dQ}{dt} + \nabla \cdot (\vec{u}^{\mbox{\scriptsize ADV}} Q) = F. \]

The class employs a predictor-corrector method to obtain a (formally) second-order accurate solution. The predicted fluxes are computed using the non-conservative advection equation, namely

\[ \frac{dQ}{dt} + (\vec{u}^{\mbox{\scriptsize ADV}} \cdot \nabla) Q = F. \]

Consequently, if the form of the source term $ F $ depends on whether $ Q $ is being conservatively or non-conservatively differenced, it is crucial that $ F $ correspond to the non-conservative form of the advection equation.

This class can also be used to solve the non-conservative form of the advection equation, i.e.,

\[ \frac{dQ}{dt} + (\vec{u}^{\mbox{\scriptsize ADV}} \cdot \nabla) Q = F. \]

Constructor & Destructor Documentation

◆ AdvectorPredictorCorrectorHyperbolicPatchOps()

IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::AdvectorPredictorCorrectorHyperbolicPatchOps ( const std::string object_name,
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database input_db,
SAMRAI::tbox::Pointer< AdvectorExplicitPredictorPatchOps explicit_predictor,
SAMRAI::tbox::Pointer< SAMRAI::geom::CartesianGridGeometry< NDIM > >  grid_geom,
bool  register_for_restart = true 
)

The constructor for AdvectorPredictorCorrectorHyperbolicPatchOps sets default parameters for the advection solver. The constructor also registers this object for restart with the restart manager using the object name.

After default values are set, this routine calls getFromRestart() if execution from a restart file is specified. Finally, getFromInput() is called to read values from the given input database (potentially overriding those found in the restart file).

◆ ~AdvectorPredictorCorrectorHyperbolicPatchOps()

IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::~AdvectorPredictorCorrectorHyperbolicPatchOps ( )
virtual

The destructor for AdvectorPredictorCorrectorHyperbolicPatchOps unregisters the patch strategy object with the restart manager when so registered.

Member Function Documentation

◆ computeFluxesOnPatch()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::computeFluxesOnPatch ( SAMRAI::hier::Patch< NDIM > &  patch,
double  time,
double  dt 
)
overridevirtual

Compute the time integral of the fluxes to be used in conservative difference for patch integration.

The conservative difference used to update the integrated quantities is implemented in conservativeDifferenceOnPatch().

Implements SAMRAI::algs::HyperbolicPatchStrategy< NDIM >.

◆ getName()

const std::string & IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::getName ( ) const

Return the name of the patch operations object.

◆ postprocessAdvanceLevelState()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::postprocessAdvanceLevelState ( const SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > &  level,
double  current_time,
double  dt,
bool  first_step,
bool  last_step,
bool  regrid_advance 
)
overridevirtual

Add source terms to the updated solution.

This routine is called after conservativeDifferenceOnPatch() is called and before computeStableDtOnPatch().

Note that when this routine is called, the scratch data is filled on all patches (i.e., ghost cells) and that data is the same as the new level data on all patch interiors. That is, both scratch and new data correspond to current_time + dt on patch interiors. The current data and ghost values correspond to the current_time.

Reimplemented from SAMRAI::algs::HyperbolicPatchStrategy< NDIM >.

Reimplemented in IBAMR::AdvDiffPredictorCorrectorHyperbolicPatchOps.

◆ preprocessAdvanceLevelState()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::preprocessAdvanceLevelState ( const SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > &  level,
double  current_time,
double  dt,
bool  first_step,
bool  last_step,
bool  regrid_advance 
)
overridevirtual

Compute the values of any time-dependent source terms for use by the explicit predictor.

This routine is called after patch boundary data is filled (i.e., ghosts) and before computeFluxesOnPatch().

Note that when this routine is called, the scratch data is filled on all patches (i.e., ghost cells) and that data is the same as the current level data on all patch interiors. That is, both scratch and current data correspond to current_time.

Reimplemented from SAMRAI::algs::HyperbolicPatchStrategy< NDIM >.

Reimplemented in IBAMR::AdvDiffPredictorCorrectorHyperbolicPatchOps.

◆ putToDatabase()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::putToDatabase ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db)
overridevirtual

Write state of AdvectorPredictorCorrectorHyperbolicPatchOps object to the given database for restart.

This routine is a concrete implementation of the function declared in the SAMRAI::tbox::Serializable abstract base class.

Implements SAMRAI::tbox::Serializable.

◆ registerAdvectionVelocity()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::registerAdvectionVelocity ( SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > >  u_var)

Register a face-centered advection velocity to be used to advect cell-centered quantities by the hierarchy integrator.

Note
By default, each registered advection velocity is assumed to be divergence free.

◆ registerModelVariables()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::registerModelVariables ( SAMRAI::algs::HyperbolicLevelIntegrator< NDIM > *  integrator)
overridevirtual

Register AdvectorPredictorCorrectorHyperbolicPatchOps model variables with the SAMRAI::algs::HyperbolicLevelIntegrator according to the variable registration function provided by the integrator.

In other words, variables are registered according to their role in the integration process (e.g. time-dependent, flux, etc.). This routine also registers variables for plotting with the VisIt writer.

Implements SAMRAI::algs::HyperbolicPatchStrategy< NDIM >.

◆ registerSourceTerm()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::registerSourceTerm ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  F_var)

Register a cell-centered source term.

◆ registerTransportedQuantity()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::registerTransportedQuantity ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var)

Register a cell-centered quantity to be advected and diffused by the hierarchy integrator.

◆ registerVisItDataWriter()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::registerVisItDataWriter ( SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > >  visit_writer)

Register a VisIt data writer so this class will write plot files that may be postprocessed with the VisIt visualization tool.

◆ setAdvectionVelocity()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::setAdvectionVelocity ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var,
SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > >  u_var 
)

Set the face-centered advection velocity to be used with a particular cell-centered quantity.

Note
The specified advection velocity must have been already registered with the hierarchy integrator.

◆ setAdvectionVelocityFunction()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::setAdvectionVelocityFunction ( SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > >  u_var,
SAMRAI::tbox::Pointer< IBTK::CartGridFunction u_fcn 
)

Set an IBTK::CartGridFunction object that specifies the value of a particular advection velocity.

◆ setAdvectionVelocityIsDivergenceFree()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::setAdvectionVelocityIsDivergenceFree ( SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > >  u_var,
bool  is_div_free 
)

Indicate whether a particular advection velocity is discretely divergence free.

◆ setConvectiveDifferencingType()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::setConvectiveDifferencingType ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var,
ConvectiveDifferencingType  difference_form 
)

Set the convective differencing form for a quantity that has been registered with the hierarchy integrator.

◆ setInitialConditions()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::setInitialConditions ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var,
SAMRAI::tbox::Pointer< IBTK::CartGridFunction Q_init 
)

Set a grid function to provide initial conditions for a quantity that has been registered with the hierarchy integrator.

◆ setPhysicalBcCoefs() [1/2]

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::setPhysicalBcCoefs ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var,
SAMRAI::solv::RobinBcCoefStrategy< NDIM > *  Q_bc_coef 
)

Set an object to provide boundary conditions for a scalar-valued quantity that has been registered with the hierarchy integrator.

◆ setPhysicalBcCoefs() [2/2]

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::setPhysicalBcCoefs ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var,
std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * >  Q_bc_coef 
)

Set objects to provide boundary conditions for a vector-valued quantity that has been registered with the hierarchy integrator.

◆ setSourceTerm()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::setSourceTerm ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var,
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  F_var 
)

Set the cell-centered source term to be used with a particular cell-centered quantity.

Note
The specified source term must have been already registered with the hierarchy integrator.

◆ setSourceTermFunction()

void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::setSourceTermFunction ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  F_var,
SAMRAI::tbox::Pointer< IBTK::CartGridFunction F_fcn 
)

Set an IBTK::CartGridFunction object that specifies the value of a particular source term.


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