IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
Class AdvDiffPredictorCorrectorHyperbolicPatchOps is a specialization of class AdvectorPredictorCorrectorHyperbolicPatchOps for use with a linearly-implicit time integrator for the advection-diffusion equation. More...
#include </home/runner/work/IBAMR/IBAMR/include/ibamr/AdvDiffPredictorCorrectorHyperbolicPatchOps.h>
Public Member Functions | |
AdvDiffPredictorCorrectorHyperbolicPatchOps (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) | |
~AdvDiffPredictorCorrectorHyperbolicPatchOps ()=default | |
void | conservativeDifferenceOnPatch (SAMRAI::hier::Patch< NDIM > &patch, double time, double dt, bool at_synchronization) override |
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 |
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 |
Public Member Functions inherited from 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) | |
virtual | ~AdvectorPredictorCorrectorHyperbolicPatchOps () |
const std::string & | getName () 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 | 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::VariableContext > | getDataContext () 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 |
Class AdvDiffPredictorCorrectorHyperbolicPatchOps is a specialization of class AdvectorPredictorCorrectorHyperbolicPatchOps for use with a linearly-implicit time integrator for the advection-diffusion equation.
IBAMR::AdvDiffPredictorCorrectorHyperbolicPatchOps::AdvDiffPredictorCorrectorHyperbolicPatchOps | ( | 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 AdvDiffPredictorCorrectorHyperbolicPatchOps sets default parameters for the patch strategy. The constructor also registers this object for restart with the restart manager using the object name when so requested.
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).
|
default |
The destructor for AdvDiffPredictorCorrectorHyperbolicPatchOps unregisters the patch strategy object with the restart manager when so registered.
|
overridevirtual |
Update solution variables by performing a conservative difference using the fluxes calculated in computeFluxesOnPatch().
Reimplemented from IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps.
|
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 IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps.
|
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 IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps.