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

Class CFINSForcing provides an interface for specifying a viscoelastic stress to be added to the Navier-Stokes equations. The class uses the advection diffusion integrator to update the viscoelastic stress. More...

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

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

Public Member Functions

 CFINSForcing (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, SAMRAI::tbox::Pointer< IBTK::CartGridFunction > u_fcn, SAMRAI::tbox::Pointer< SAMRAI::geom::CartesianGridGeometry< NDIM > > grid_geometry, SAMRAI::tbox::Pointer< IBAMR::AdvDiffSemiImplicitHierarchyIntegrator > adv_diff_integrator, SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > visit_data_writer)
 This constructor creates Variable and VariableContext objects for storing the viscoleastic stresses at the centers of the Cartesian grid. Sets up the advection diffusion solver to use the velocity function prescribed. The initial and boundary conditions should be specified on the quantity being solved for (e.g. Conformation tensor or square root or logarithm of the conformation tensor).
 
 CFINSForcing (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > app_initializer, const SAMRAI::tbox::Pointer< IBAMR::INSHierarchyIntegrator > fluid_solver, SAMRAI::tbox::Pointer< SAMRAI::geom::CartesianGridGeometry< NDIM > > grid_geometry, SAMRAI::tbox::Pointer< IBAMR::AdvDiffSemiImplicitHierarchyIntegrator > adv_diff_integrator, SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > visit_data_writer)
 This constructor creates Variable and VariableContext objects for storing the viscoleastic stresses at the centers of the Cartesian grid. Sets up the advection diffusion solver to use the fluid velocity from the fluid solver. Note that this function must be registered with the fluid solver as a forcing function for the viscoelastic stress to effect the fluid velocity. The initial and boundary conditions should be specified on the quantity being solved for (e.g. Conformation tensor or square root or logarithm of the conformation tensor).
 
 ~CFINSForcing ()
 Deallocates draw data and deletes boundary conditions.
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > getVariable ()
 This function returns a pointer to the cell variable that stores the viscoelastic stress.
 
int getVariableIdx ()
 This function returns the patch data index used to store the viscoelastic stress.
 
void registerCFStrategy (SAMRAI::tbox::Pointer< CFStrategy > rhs)
 This function registers a strategy with the class. More...
 
- Public Member Functions inherited from IBTK::CartGridFunction
 CartGridFunction (std::string object_name="")
 The default constructor sets the name of the strategy object.
 
virtual ~CartGridFunction ()=default
 Empty virtual destructor.
 

Methods to set patch data.

bool isTimeDependent () const override
 Indicates whether the concrete INSStaggeredStochasticForcing object is time-dependent.
 
void setDataOnPatchHierarchy (const int data_idx, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > var, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, const double data_time, const bool initial_time=false, const int coarsest_ln=IBTK::invalid_level_number, const int finest_ln=IBTK::invalid_level_number) override
 Compute the divergence of the stress tensor. Also sets up requested visualizations.
 
void setDataOnPatch (const int data_idx, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > var, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const double data_time, const bool initial_time=false, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > patch_level=SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > >(NULL)) override
 Evaluate the divergence on the patch interior.
 
void setDataOnPatchLevel (const int data_idx, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > var, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > level, const double data_time, const bool initial_time) override
 Evaluate the divergence on the specified patch level.
 
void checkPositiveDefinite (const int data_idx, const SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > var, const double data_time, const bool initial_time)
 Check whether the provided patch index stores a positive definite tensor.
 
void applyGradientDetector (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double error_data_time, int tag_index, bool initial_time, bool)
 Tag cells based on the specifications provided in the input database.
 
void projectTensor (const int data_idx, const SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > var, const double data_time, const bool initial_time, const bool extended_box)
 Projects the symmetric tensor stored in data_idx to the nearest non negative matrix in the L2 norm.
 
SAMRAI::tbox::Pointer< AdvDiffSemiImplicitHierarchyIntegratorgetAdvDiffHierarchyIntegrator ()
 Return the advection diffusion integrator used to evolve the conformation tensor.
 
static void apply_gradient_detector_callback (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double error_data_time, int tag_index, bool initial_time, bool richardson_extrapolation_too, void *ctx)
 
static void apply_project_tensor_callback (double current_time, double new_time, int cycle_num, void *ctx)
 

Additional Inherited Members

- Protected Attributes inherited from IBTK::CartGridFunction
std::string d_object_name
 

Detailed Description

Class CFINSForcing provides an interface for specifying a viscoelastic stress to be added to the Navier-Stokes equations. The class uses the advection diffusion integrator to update the viscoelastic stress.

Users must register a strategy operator CFStrategy that can compute (1) the relaxation of the stress and (2) the conversion from conformation tensor to stress tensor. One can choose from the pre-programmed models Oldroyd-B, Giesekus, or Rolie-Poly. Optionally, you can also register a your own strategy. The fluid model is specified through the database parameter "fluid_parameter". By specifying "USER_DEFINED", you can register your own strategy operator. This class currently solves for the conformation tensor or the square root or logarithm of the conformation tensor.

Input Database Parameters

There are also strategy specific items searched for in the database.

Member Function Documentation

◆ registerCFStrategy()

void IBAMR::CFINSForcing::registerCFStrategy ( SAMRAI::tbox::Pointer< CFStrategy rhs)

This function registers a strategy with the class.

Note
This function only should be called when fluid_model is set to USER_DEFINED in the input file.

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