|
IBAMR
IBAMR version 0.19.
|
Class AdvectorPredictorCorrectorHyperbolicPatchOps is a concrete SAMRAI::algs::HyperbolicPatchStrategy that makes use of class AdvectorExplicitPredictorPatchOps to solve the linear advection equation. More...
#include <ibamr/AdvectorPredictorCorrectorHyperbolicPatchOps.h>

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::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. More... | |
| 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. More... | |
| 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(). More... | |
| 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. More... | |
| 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. More... | |
| virtual void | putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override |
| Write state of AdvectorPredictorCorrectorHyperbolicPatchOps object to the given database for restart. More... | |
| virtual void | registerModelVariables (HyperbolicLevelIntegrator< DIM > *integrator)=0 |
| virtual void | setupLoadBalancer (HyperbolicLevelIntegrator< DIM > *integrator, mesh::GriddingAlgorithm< DIM > *gridding_algorithm) |
| virtual void | initializeDataOnPatch (hier::Patch< DIM > &patch, const double data_time, const bool initial_time)=0 |
| virtual double | computeStableDtOnPatch (hier::Patch< DIM > &patch, const bool initial_time, const double dt_time)=0 |
| virtual void | computeFluxesOnPatch (hier::Patch< DIM > &patch, const double time, const double dt)=0 |
| virtual void | conservativeDifferenceOnPatch (hier::Patch< DIM > &patch, const double time, const double dt, bool at_syncronization)=0 |
| virtual void | preprocessAdvanceLevelState (const tbox::Pointer< hier::PatchLevel< DIM > > &level, double current_time, double dt, bool first_step, bool last_step, bool regrid_advance) |
| virtual void | postprocessAdvanceLevelState (const tbox::Pointer< hier::PatchLevel< DIM > > &level, double current_time, double dt, bool first_step, bool last_step, bool regrid_advance) |
| virtual void | tagGradientDetectorCells (hier::Patch< DIM > &patch, const double regrid_time, const bool initial_error, const int tag_index, const bool uses_richardson_extrapolation_too) |
| virtual void | tagRichardsonExtrapolationCells (hier::Patch< DIM > &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 void | setPhysicalBoundaryConditions (hier::Patch< DIM > &patch, const double fill_time, const hier::IntVector< DIM > &ghost_width_to_fill)=0 |
| virtual hier::IntVector< DIM > | getRefineOpStencilWidth () const |
| virtual hier::IntVector< DIM > | getRefineOpStencilWidth () const =0 |
| virtual void | preprocessRefine (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::Box< DIM > &fine_box, const hier::IntVector< DIM > &ratio) |
| virtual void | postprocessRefine (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::Box< DIM > &fine_box, const hier::IntVector< DIM > &ratio) |
| virtual hier::IntVector< DIM > | getCoarsenOpStencilWidth () const |
| virtual hier::IntVector< DIM > | getCoarsenOpStencilWidth () const =0 |
| virtual void | preprocessCoarsen (hier::Patch< DIM > &coarse, const hier::Patch< DIM > &fine, const hier::Box< DIM > &coarse_box, const hier::IntVector< DIM > &ratio) |
| virtual void | postprocessCoarsen (hier::Patch< DIM > &coarse, const hier::Patch< DIM > &fine, const hier::Box< DIM > &coarse_box, const hier::IntVector< DIM > &ratio) |
| tbox::Pointer< hier::VariableContext > | getDataContext () const |
| void | setDataContext (tbox::Pointer< hier::VariableContext > context) |
| void | clearDataContext () |
| 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) |
Private Member Functions | |
| AdvectorPredictorCorrectorHyperbolicPatchOps ()=delete | |
| Default constructor. More... | |
| AdvectorPredictorCorrectorHyperbolicPatchOps (const AdvectorPredictorCorrectorHyperbolicPatchOps &from)=delete | |
| Copy constructor. More... | |
| AdvectorPredictorCorrectorHyperbolicPatchOps & | operator= (const AdvectorPredictorCorrectorHyperbolicPatchOps &that)=delete |
| Assignment operator. More... | |
| void | setInflowBoundaryConditions (SAMRAI::hier::Patch< NDIM > &patch, double fill_time) |
| void | getFromInput (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db, bool is_from_restart) |
| void | getFromRestart () |
Private Attributes | |
| std::string | d_object_name |
| bool | d_registered_for_restart |
| SAMRAI::tbox::Pointer< SAMRAI::geom::CartesianGridGeometry< NDIM > > | d_grid_geometry |
| SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > | d_visit_writer |
| IBTK::CartExtrapPhysBdryOp | d_extrap_bc_helper |
| SAMRAI::hier::IntVector< NDIM > | d_ghosts |
| SAMRAI::hier::IntVector< NDIM > | d_flux_ghosts |
| std::string | d_extrap_type = "CONSTANT" |
| SAMRAI::tbox::Array< std::string > | d_refinement_criteria |
| SAMRAI::tbox::Array< double > | d_dev_tol |
| SAMRAI::tbox::Array< double > | d_dev |
| SAMRAI::tbox::Array< double > | d_dev_time_max |
| SAMRAI::tbox::Array< double > | d_dev_time_min |
| SAMRAI::tbox::Array< double > | d_grad_tol |
| SAMRAI::tbox::Array< double > | d_grad_time_max |
| SAMRAI::tbox::Array< double > | d_grad_time_min |
| tbox::Pointer< hier::VariableContext > | d_data_context |
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. \]
| 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).
|
virtual |
The destructor for AdvectorPredictorCorrectorHyperbolicPatchOps unregisters the patch strategy object with the restart manager when so registered.
|
privatedelete |
|
privatedelete |
| from | The value to copy to this object. |
| const std::string& IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::getName | ( | ) | const |
Return the name of the patch operations object.
| 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.
| 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.
| 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.
| 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.
| void IBAMR::AdvectorPredictorCorrectorHyperbolicPatchOps::registerSourceTerm | ( | SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > | F_var | ) |
Register a cell-centered source term.
| 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.
| 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.
| 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.
| 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.
| 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.
| 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.
| 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.
| 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.
|
overridevirtual |
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.
|
overridevirtual |
|
overridevirtual |
|
overridevirtual |
The conservative difference used to update the integrated quantities is implemented in conservativeDifferenceOnPatch().
|
overridevirtual |
Reimplemented in IBAMR::AdvDiffPredictorCorrectorHyperbolicPatchOps.
|
overridevirtual |
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 in IBAMR::AdvDiffPredictorCorrectorHyperbolicPatchOps.
|
overridevirtual |
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 in IBAMR::AdvDiffPredictorCorrectorHyperbolicPatchOps.
|
overridevirtual |
|
overridevirtual |
|
overridevirtual |
This routine is a concrete implementation of the function declared in the SAMRAI::tbox::Serializable abstract base class.
Implements SAMRAI::tbox::Serializable.
|
protected |
|
protected |
|
protected |
|
privatedelete |
| that | The value to assign to this object. |
|
private |
|
private |
|
private |
|
pure virtualinherited |
Register specific variables needed in the numerical routines with the hyperbolic level integrator using the registerVariable() function in that class. The integrator manipulates storage for the data and this registration defines the way in which data for each quantity will be manipulated on the patches. Typically, the derived data quantities for plotting are registered with a visualization data writer in this routine as well, since the hyperbolic level integrator provides the variable context for plotting (i.e., which data is available when a plot file is generated). The integrator pointer cannot be null in most cases.
The gridding algorithm pointer is provided so that patch data objects may be registered with the load balancer object (owned by the gridding algorithm) for non-uniform load balancing, if needed.
|
inlinevirtualinherited |
Set up parameters in the load balancer object (owned by the gridding algorithm) if needed. This function is called immediately after the registerModelVariables() function is called. The hyperbolic level integrator pointer is provided so that the integrator can be used to manage data for the load balancer if needed (e.g., when using non-uniform load balancing).
Note that this function is not pure virtual. It is given a dummy implementation here so that users may ignore it when inheriting from this class.
|
pure virtualinherited |
Set the initial data on a patch interior only. Note that no ghost cells need to be set in this routine regardless of whether the patch data corresponding to the context requires ghost cells. The data_time is the simulation time when the routine is called. The boolean initial_time is true if the routine is called at the initial time when the hierarchy is initially constructed, otherwise it is false.
|
pure virtualinherited |
Compute the stable time increment for a patch on the level with the given number. The boolean flag initial_time is true if the routine is called at the initial simulation time; otherwise it is false. The double argument dt_time is the simulation time.
|
pure virtualinherited |
Compute TIME INTEGRALS of fluxes to be used in conservative difference for patch integration. That is, it is assumed that this numerical routine will compute the fluxes corresponding to the cell faces multiplied by the time increment. Typically, the numerical flux is the normal flux at the cell face. The flux integrals will be used in the conservative difference that updates the conserved quantities.
Note that the numerical routines in this method generally require ghost cells. Ghost cells data is filled before this routine is called.
|
pure virtualinherited |
Update patch data with a conservative difference (approximating the divergence theorem) using the flux integrals computed in computeFluxesOnPatch() routine. The boolean flag is true when this routine is called during a flux synchronization step. Otherwise, it is false. Note that the computeFluxesOnPatch() routine computes TIME INTEGRALs of the numerical fluxes (e.g., they have been multiplied by dt). So the conservative difference routine should be consistent with this.
|
inlinevirtualinherited |
This is an optional routine for user to process any application-specific patch strategy data BEFORE patches are advanced on the given level. This routine is called after patch boundary data is filled (i.e., ghosts) and before computeFluxesOnPatch(). The arguments are: level – level that will be advanced, current_time – current integration time, dt – current time increment, first_step – boolean flag that is true if advance is first in time step sequence on level (i.e., previous advance step was on another level, false otherwise, last_step – boolean flag that is true if advance is last in time step sequence on level (i.e., synchronization with coarser level will occur immediately after this advance), regrid_advance – boolean flag that is true if the advance is during a regridding phase (i.e., the advance is not used to integrate data on the hierarchy) in which case the results of the advance will be discarded.
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.
Note that this function is not pure virtual. It is given a dummy implementation here so that users may ignore it when inheriting from this class.
|
inlinevirtualinherited |
This is an optional routine for user to process any application-specific patch strategy data AFTER patches are advanced on the given level. This routine is called after conservativeDifferenceOnPatch() is called and before computeStableDtOnPatch(). The arguments are: level – level that will be advanced, current_time – current integration time, dt – current time increment, first_step – boolean flag that is true if advance is first in time step sequence on level (i.e., previous advance step was on another level, false otherwise, last_step – boolean flag that is true if advance is last in time step sequence on level (i.e., synchronization with coarser level will occur immediately after this advance), regrid_advance – boolean flag that is true if the advance is during a regridding phase (i.e., the advance is not used to integrate data on the hierarchy) in which case the results of the advance will be discarded.
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.
Note that this function is not pure virtual. It is given a dummy implementation here so that users may ignore it when inheriting from this class.
|
virtualinherited |
Tag cells on the given patch that require refinement based on application-specific numerical quantities. The tag index argument indicates the index of the tag data on the patch data array. The boolean argument initial_error is true if tagging is being done at the initial simulation time; otherwise, it is false. The other boolean flag uses_richardson_extrapolation_too is true when Richardson extrapolation is used in addition to the gradient detector. This flag helps users manage multiple regridding criteria.
Note that this function is not pure virtual. It is given a dummy implementation here so that users may ignore it when inheriting from this class.
|
virtualinherited |
Tag cells based from differences computed in the Richardson extrapolation. The Richardson extrapolation algorithm creates a coarsened version of some hierarchy patch level and advances data in time on both the coarsened patch level and the hierarchy level. This routine takes the data resulting from the advance on both the coarse and fine levels, compares them, and tags cells according to the difference.
* (2) * n+1 ^------->x finish (1) advanced_coarse * | ^ (2) coarsened_fine * time n - | * ^ |(1) * | | * <--------o start * fine coarse *
The patch supplied to this routine is on the coarsened level. However, the error_level_number corresponds to the actual hierarchy level from which it was coarsened. Data resides on this patch in two contexts - advanced_coarse'' andcoarsened_fine''. Advanced coarse is data advanced on the coarsened version of the level, while coarsened fine is the data advanced on the fine level and then coarsened to the coarse level. The regrid time and the time increment are given for the actual hierarchy level. The error coarsen ratio argument is the ratio between the index spaces on the hierarchy level and the coarsened hierarchy level. The boolean flag `‘initial_error’' is true when the error estimation is performed at the initial simulation time; i.e., when the hierarchy levels are being constructed for the first time. The tag index argument is the index of the tag data on the patch data array. The other boolean flag uses_gradient_detector_too is true when a gradient detector scheme is used in addition to Richardson extrapolation. This flag helps users manage multiple regridding criteria.
Note that this function is not pure virtual. It is given a dummy implementation here so that users may ignore it when inheriting from this class.
|
pure virtualinherited |
Set user-defined boundary conditions at the physical domain boundary.
Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.
|
inlinevirtualinherited |
Return maximum stencil width needed for user-defined data interpolation operations. Default is to return zero, assuming no user-defined operations provided.
Note that this function is not pure virtual. It is given a dummy implementation here so that users may ignore it when inheriting from this class.
|
pure virtualinherited |
Pure virtual function to return maximum stencil width needed over user-defined data interpolation operations. This is needed to determine the correct interpolation data dependencies.
Implemented in SAMRAI::xfer::MultiblockRefinePatchStrategy< DIM >, SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >, SAMRAI::algs::HyperbolicPatchStrategy< DIM >, SAMRAI::algs::MethodOfLinesPatchStrategy< DIM >, SAMRAI::mesh::MultiblockGriddingTagger< DIM >, and SAMRAI::solv::CartesianRobinBcHelper< DIM >.
|
inlinevirtualinherited |
Pre- and post-processing routines for implementing user-defined spatial interpolation routines applied to variables. The interpolation routines are used in the hyperbolic AMR algorithm for filling patch ghost cells before advancing data on a level and after regridding a level to fill portions of the new level from some coarser level. These routines are called automatically from within patch boundary filling schedules; thus, some concrete function matching these signatures must be provided in the user's patch routines. However, the routines only need to perform some operations when "USER_DEFINED_REFINE" is given as the interpolation method for some variable when the patch routines register variables with the hyperbolic level integration algorithm, typically. If the user does not provide operations that refine such variables in either of these routines, then they will not be refined.
The order in which these operations are used in each patch boundary filling schedule is:
Note that these functions are not pure virtual. They are given dummy implementations here so that users may ignore them when inheriting from this class.
Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.
|
inlinevirtualinherited |
Pure virtual function to perform user-defined postprocess data refine operations. This member function is called after standard refine operations (expressed using concrete subclasses of the RefineOperator<DIM> base class). The postprocess function must refine data from the scratch components of the coarse patch into the scratch components of the fine patch on the specified fine box region. Recall that the scratch components are specified in calls to the registerRefine() function in the RefineAlgorithm<DIM> class.
| fine | Fine patch containing destination data. |
| coarse | Coarse patch containing source data. |
| fine_box | hier::Box region on fine patch into which data is refined. |
| ratio | Integer vector containing ratio relating index space between coarse and fine patches. |
Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.
|
inlinevirtualinherited |
Return maximum stencil width needed for user-defined data coarsen operations. Default is to return zero, assuming no user-defined operations provided.
Note that this function is not pure virtual. It is given a dummy implementation here so that users may ignore it when inheriting from this class.
|
pure virtualinherited |
Return maximum stencil width needed over all user-defined data coarsening operations. This is needed to determine the correct coarsening data dependencies.
Implemented in SAMRAI::algs::HyperbolicPatchStrategy< DIM >, and SAMRAI::algs::MethodOfLinesPatchStrategy< DIM >.
|
inlinevirtualinherited |
Pre- and post-processing routines for implementing user-defined spatial coarsening routines applied to variables. The coarsening routines are used in the hyperbolic AMR algorithm synchronizing coarse and fine levels when they have been integrated to the same point. These routines are called automatically from within the data synchronization coarsen schedules; thus, some concrete function matching these signatures must be provided in the user's patch routines. However, the routines only need to perform some operations when "USER_DEFINED_COARSEN" is given as the coarsening method for some variable when the patch routines register variables with the hyperbolic level integration algorithm, typically. If the user does not provide operations that coarsen such variables in either of these routines, then they will not be coarsened.
The order in which these operations are used in each coarsening schedule is:
Note that these functions are not pure virtual. They are given dummy implementations here so that users may ignore them when inheriting from this class.
Implements SAMRAI::xfer::CoarsenPatchStrategy< DIM >.
|
inlinevirtualinherited |
Perform user-defined coarsening operations. This member function is called after standard coarsening operations (expressed using concrete subclasses of the CoarsenOperator<DIM> base class).
The postprocess function should move data from the source components on the fine patch into the source components on the coarse patch in the specified coarse box region. Recall that the source components are specified in calls to the registerCoarsen() function in the CoarsenAlgorithm<DIM> class.
| coarse | Coarse patch containing destination data. |
| fine | Fine patch containing source data. |
| coarse_box | hier::Box region on coarse patch into which data is copied. |
| ratio | Integer vector containing ratio |
Implements SAMRAI::xfer::CoarsenPatchStrategy< DIM >.
|
inlineinherited |
Return pointer to patch data context.
|
inlineinherited |
The hyperbolic integrator controls the context for the data to be used in the numerical routines implemented in the concrete patch strategy. The setDataContext() allows the integrator to set the context for data on a patch on which to operate.
|
inlineinherited |
The clearDataContext() routine resets the data context to be null.
|
virtualinherited |
Virtual function to perform user-defined refine operations. This member function is called before standard refining operations (expressed using concrete subclasses of the RefineOperator<DIM> base class). The preprocess function must refine data from the scratch components of the coarse patch into the scratch components of the fine patch on the specified fine box regions.
Typically, only the pure virtual members of this class are implemented in user-defined subclasses of this base class. This version of the preprocess function operates on an entire box list. By default, this version simply loops over the box list and calls the single-box version, which is a pure virtual function.
| fine | Fine patch containing destination data. |
| coarse | Coarse patch containing source data. |
| fine_boxes | tbox::List of box regions on fine patch into which data is refined. |
| ratio | Integer vector containing ratio relating index space between coarse and fine patches. |
Reimplemented in SAMRAI::solv::CartesianRobinBcHelper< DIM >.
|
virtualinherited |
Virtual function to perform user-defined refine operations. This member function is called after standard refining operations (expressed using concrete subclasses of the RefineOperator<DIM> base class). The postprocess function must refine data from the scratch components of the coarse patch into the scratch components of the fine patch on the specified fine box regions.
Typically, only the pure virtual members of this class are implemented in user-defined subclasses of this base class. This version of the postprocess function operates on an entire box list. By default, this version simply loops over the box list and calls the single-box version, which is a pure virtual function.
| fine | Fine patch containing destination data. |
| coarse | Coarse patch containing source data. |
| fine_boxes | tbox::List of box regions on fine patch into which data is refined. |
| ratio | Integer vector containing ratio relating index space between coarse and fine patches. |
Reimplemented in SAMRAI::solv::CartesianRobinBcHelper< DIM >.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
privateinherited |
1.8.17