IBAMR  IBAMR version 0.19.
Public Member Functions | List of all members
IBAMR::MarangoniSurfaceTensionForceFunction Class Reference

Class MarangoniSurfaceTensionForceFunction provides Marangoni forcing due to temperature variations. This class computes the forcing term \( F = \frac{\mathrm{d} \sigma}{\mathrm{d} T}(\nabla T |\nabla C| - (\nabla T \cdot \nabla \phi) \nabla C)\) and add it to the surface tension forcing. More...

#include <ibamr/MarangoniSurfaceTensionForceFunction.h>

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

Public Member Functions

 MarangoniSurfaceTensionForceFunction (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, AdvDiffHierarchyIntegrator *adv_diff_solver, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > level_set_var, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > T_var, SAMRAI::solv::RobinBcCoefStrategy< NDIM > *&T_bc_coef)
 Constructor. More...
 
virtual ~MarangoniSurfaceTensionForceFunction ()=default
 Destructor. More...
 
virtual void setSmoother (const std::string &kernel_fcn)
 Set the smoother (kernel function) to mollify the Heaviside function. More...
 
virtual void setSurfaceTensionCoef (double sigma)
 Set the constant surface tension coefficient. More...
 
virtual void setNumberOfInterfaceCells (double m)
 Set the number of interface cells m over which the surface tension force will be applied. The surface tension will take effect in the band -m*h to m*h around the interface, where h = (dx*dy)^(1/2) in 2D and h = (dx*dy*dz)^(1/3) in 3D. More...
 
std::string getSmoother () const
 Get the smoother (kernel function) to mollify the Heaviside function. More...
 
double getSurfaceTensionCoef () const
 Get the constant surface tension coefficient. More...
 
double getNumberOfInterfaceCells () const
 Get the number of interface cells over which the surface tension force will be applied. More...
 

Methods to set the data.

using ComputeMarangoniCoefPtr = void(*)(int F_idx, SAMRAI::tbox::Pointer< IBTK::HierarchyMathOps > hier_math_ops, int cycle_num, double time, double current_time, double new_time, void *ctx)
 Callback function to compute the marangoni coefficient as a function of temperature and multiply it with the F_data as F_data = marangoni_coef*F_data. More...
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_T_var
 
int d_T_idx = IBTK::invalid_index
 
int d_F_cloned_idx = IBTK::invalid_index
 
double d_marangoni_coefficient = 0.0
 
SAMRAI::solv::RobinBcCoefStrategy< NDIM > * d_T_bc_coef = nullptr
 
ComputeMarangoniCoefPtr d_compute_marangoni_coef = nullptr
 
void * d_compute_marangoni_coef_ctx = nullptr
 
bool isTimeDependent () const override
 
void setDataOnPatchHierarchy (int data_idx, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > var, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, double data_time, bool initial_time=false, int coarsest_ln=-1, int finest_ln=-1) override
 Evaluate the function on the patch interiors on the specified levels of the patch hierarchy using the virtual function setDataOnPatch(). More...
 
void setDataOnPatch (int data_idx, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > var, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, double data_time, bool initial_time=false, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > level=SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > >(nullptr)) override
 
void registerMarangoniCoefficientFunction (ComputeMarangoniCoefPtr callback, void *ctx)
 Register callback function to compute the variable marangoni coefficient. More...
 
 MarangoniSurfaceTensionForceFunction ()=delete
 Default constructor. More...
 
 MarangoniSurfaceTensionForceFunction (const SurfaceTensionForceFunction &from)=delete
 Copy constructor. More...
 
MarangoniSurfaceTensionForceFunctionoperator= (const SurfaceTensionForceFunction &that)=delete
 Assignment operator. More...
 
void setDataOnPatchCell (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellData< NDIM, double > > F_data, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const double data_time, const bool initial_time, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > level)
 
void setDataOnPatchSide (SAMRAI::tbox::Pointer< SAMRAI::pdat::SideData< NDIM, double > > F_data, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const double data_time, const bool initial_time, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > level)
 

Methods to set the data.

using MaskSurfaceTensionForcePtr = void(*)(int F_idx, SAMRAI::tbox::Pointer< IBTK::HierarchyMathOps > hier_math_ops, int cycle_num, double time, double current_time, double new_time, void *ctx)
 Function to Mask surface tension force to act only on the liquid-gas interface. More...
 
using ComputeSurfaceTensionCoefficientPtr = void(*)(int F_idx, SAMRAI::tbox::Pointer< IBTK::HierarchyMathOps > hier_math_ops, int cycle_num, double time, double current_time, double new_time, void *ctx)
 Function to compute the variable surface tension coefficient. More...
 
MaskSurfaceTensionForcePtr d_mask_surface_tension_force = nullptr
 
void * d_mask_surface_tension_force_ctx = nullptr
 
ComputeSurfaceTensionCoefficientPtr d_compute_surface_tension_coef = nullptr
 
void * d_compute_surface_tension_coef_ctx = nullptr
 
void registerSurfaceTensionForceMasking (MaskSurfaceTensionForcePtr callback, void *ctx)
 Register function to limit the surface tension force. More...
 
void registerSurfaceTensionCoefficientFunction (ComputeSurfaceTensionCoefficientPtr callback, void *ctx)
 Register function to compute the variable surface tension coefficient. More...
 
void convertToHeaviside (int phi_idx, int coarsest_ln, int finest_ln, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > patch_hierarchy)
 
void mollifyData (int phi_idx, int coarsest_ln, int finest_ln, double data_time, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< IBTK::HierarchyGhostCellInterpolation > fill_op)
 
int getStencilSize (const std::string &kernel_fcn)
 
const AdvDiffHierarchyIntegrator *const d_adv_diff_solver
 
const SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > d_ls_var
 
TimeSteppingType d_ts_type
 
int d_C_idx = IBTK::invalid_index
 
int d_phi_idx = IBTK::invalid_index
 
std::string d_kernel_fcn
 
double d_sigma = std::numeric_limits<double>::signaling_NaN()
 
double d_num_interface_cells = std::numeric_limits<double>::signaling_NaN()
 
SAMRAI::tbox::Pointer< IBTK::HierarchyMathOpsd_hier_math_ops
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchySideDataOpsReal< NDIM, double > > d_hier_sc_data_ops
 
int getMinimumGhostWidth (const std::string &kernel_fcn)
 

Methods to set patch interior data.

virtual void setDataOnPatchLevel (int data_idx, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > var, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > patch_level, double data_time, bool initial_time=false)
 Evaluate the function on the patch interiors on the specified level of the patch hierarchy using the virtual function setDataOnPatch(). More...
 
std::string d_object_name
 

Detailed Description

This class uses the callback function to compute the variable Marangoni coefficient as a function of Temperature.

Note
Presently, this class assumes that the indicator function is a cell centered level-set variable that is maintained by the advection-diffusion integrator. In general, the indicator variable can either be a level set function, a volume fraction function, or a phase field function.

Member Typedef Documentation

◆ ComputeMarangoniCoefPtr

using IBAMR::MarangoniSurfaceTensionForceFunction::ComputeMarangoniCoefPtr = void (*)(int F_idx, SAMRAI::tbox::Pointer<IBTK::HierarchyMathOps> hier_math_ops, int cycle_num, double time, double current_time, double new_time, void* ctx)

◆ MaskSurfaceTensionForcePtr

using IBAMR::SurfaceTensionForceFunction::MaskSurfaceTensionForcePtr = void (*)(int F_idx, SAMRAI::tbox::Pointer<IBTK::HierarchyMathOps> hier_math_ops, int cycle_num, double time, double current_time, double new_time, void* ctx)
inherited

◆ ComputeSurfaceTensionCoefficientPtr

using IBAMR::SurfaceTensionForceFunction::ComputeSurfaceTensionCoefficientPtr = void (*)(int F_idx, SAMRAI::tbox::Pointer<IBTK::HierarchyMathOps> hier_math_ops, int cycle_num, double time, double current_time, double new_time, void* ctx)
inherited

Constructor & Destructor Documentation

◆ MarangoniSurfaceTensionForceFunction() [1/3]

IBAMR::MarangoniSurfaceTensionForceFunction::MarangoniSurfaceTensionForceFunction ( const std::string &  object_name,
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database input_db,
AdvDiffHierarchyIntegrator adv_diff_solver,
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > >  level_set_var,
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > >  T_var,
SAMRAI::solv::RobinBcCoefStrategy< NDIM > *&  T_bc_coef 
)

◆ ~MarangoniSurfaceTensionForceFunction()

virtual IBAMR::MarangoniSurfaceTensionForceFunction::~MarangoniSurfaceTensionForceFunction ( )
virtualdefault

◆ MarangoniSurfaceTensionForceFunction() [2/3]

IBAMR::MarangoniSurfaceTensionForceFunction::MarangoniSurfaceTensionForceFunction ( )
privatedelete
Note
This constructor is not implemented and should not be used.

◆ MarangoniSurfaceTensionForceFunction() [3/3]

IBAMR::MarangoniSurfaceTensionForceFunction::MarangoniSurfaceTensionForceFunction ( const SurfaceTensionForceFunction from)
privatedelete
Note
This constructor is not implemented and should not be used.
Parameters
fromThe value to copy to this object.

Member Function Documentation

◆ isTimeDependent()

bool IBAMR::MarangoniSurfaceTensionForceFunction::isTimeDependent ( ) const
overridevirtual
Note
This concrete IBTK::CartGridFunction is time-dependent.

Implements IBTK::CartGridFunction.

◆ setDataOnPatchHierarchy()

void IBAMR::MarangoniSurfaceTensionForceFunction::setDataOnPatchHierarchy ( int  data_idx,
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > >  var,
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  hierarchy,
double  data_time,
bool  initial_time = false,
int  coarsest_ln = -1,
int  finest_ln = -1 
)
overridevirtual
See also
setDataOnPatch

Reimplemented from IBTK::CartGridFunction.

◆ setDataOnPatch()

void IBAMR::MarangoniSurfaceTensionForceFunction::setDataOnPatch ( int  data_idx,
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > >  var,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
double  data_time,
bool  initial_time = false,
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > >  level = SAMRAI::tbox::PointerSAMRAI::hier::PatchLevel< NDIM > >(nullptr) 
)
overridevirtual

Set the data on the patch interior.

Implements IBTK::CartGridFunction.

◆ registerMarangoniCoefficientFunction()

void IBAMR::MarangoniSurfaceTensionForceFunction::registerMarangoniCoefficientFunction ( ComputeMarangoniCoefPtr  callback,
void *  ctx 
)

◆ operator=()

MarangoniSurfaceTensionForceFunction& IBAMR::MarangoniSurfaceTensionForceFunction::operator= ( const SurfaceTensionForceFunction that)
privatedelete
Note
This operator is not implemented and should not be used.
Parameters
thatThe value to assign to this object.
Returns
A reference to this object.

◆ setDataOnPatchCell()

void IBAMR::MarangoniSurfaceTensionForceFunction::setDataOnPatchCell ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellData< NDIM, double > >  F_data,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const double  data_time,
const bool  initial_time,
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > >  level 
)
private

Set the data on the patch interior.

◆ setDataOnPatchSide()

void IBAMR::MarangoniSurfaceTensionForceFunction::setDataOnPatchSide ( SAMRAI::tbox::Pointer< SAMRAI::pdat::SideData< NDIM, double > >  F_data,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const double  data_time,
const bool  initial_time,
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > >  level 
)
private

Set the data on the patch interior.

◆ setSmoother()

virtual void IBAMR::SurfaceTensionForceFunction::setSmoother ( const std::string &  kernel_fcn)
virtualinherited

◆ setSurfaceTensionCoef()

virtual void IBAMR::SurfaceTensionForceFunction::setSurfaceTensionCoef ( double  sigma)
virtualinherited

◆ setNumberOfInterfaceCells()

virtual void IBAMR::SurfaceTensionForceFunction::setNumberOfInterfaceCells ( double  m)
virtualinherited

◆ getSmoother()

std::string IBAMR::SurfaceTensionForceFunction::getSmoother ( ) const
inlineinherited

◆ getSurfaceTensionCoef()

double IBAMR::SurfaceTensionForceFunction::getSurfaceTensionCoef ( ) const
inlineinherited

◆ getNumberOfInterfaceCells()

double IBAMR::SurfaceTensionForceFunction::getNumberOfInterfaceCells ( ) const
inlineinherited

◆ registerSurfaceTensionForceMasking()

void IBAMR::SurfaceTensionForceFunction::registerSurfaceTensionForceMasking ( MaskSurfaceTensionForcePtr  callback,
void *  ctx 
)
inherited

◆ registerSurfaceTensionCoefficientFunction()

void IBAMR::SurfaceTensionForceFunction::registerSurfaceTensionCoefficientFunction ( ComputeSurfaceTensionCoefficientPtr  callback,
void *  ctx 
)
inherited

◆ getMinimumGhostWidth()

int IBAMR::SurfaceTensionForceFunction::getMinimumGhostWidth ( const std::string &  kernel_fcn)
protectedinherited

Get the ghost cell width of scratch data.

◆ convertToHeaviside()

void IBAMR::SurfaceTensionForceFunction::convertToHeaviside ( int  phi_idx,
int  coarsest_ln,
int  finest_ln,
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  patch_hierarchy 
)
privateinherited

Convert the level set variable to a smoothed heaviside function.

◆ mollifyData()

void IBAMR::SurfaceTensionForceFunction::mollifyData ( int  phi_idx,
int  coarsest_ln,
int  finest_ln,
double  data_time,
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  hierarchy,
SAMRAI::tbox::Pointer< IBTK::HierarchyGhostCellInterpolation fill_op 
)
privateinherited

Mollify data.

◆ getStencilSize()

int IBAMR::SurfaceTensionForceFunction::getStencilSize ( const std::string &  kernel_fcn)
privateinherited

Get the stencil size for the kernel.

◆ setDataOnPatchLevel()

virtual void IBTK::CartGridFunction::setDataOnPatchLevel ( int  data_idx,
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > >  var,
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > >  patch_level,
double  data_time,
bool  initial_time = false 
)
virtualinherited

Member Data Documentation

◆ d_T_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::MarangoniSurfaceTensionForceFunction::d_T_var
private

Temperature variable and its patch data index.

◆ d_T_idx

int IBAMR::MarangoniSurfaceTensionForceFunction::d_T_idx = IBTK::invalid_index
private

◆ d_F_cloned_idx

int IBAMR::MarangoniSurfaceTensionForceFunction::d_F_cloned_idx = IBTK::invalid_index
private

◆ d_marangoni_coefficient

double IBAMR::MarangoniSurfaceTensionForceFunction::d_marangoni_coefficient = 0.0
private

Marangoni coefficient.

◆ d_T_bc_coef

SAMRAI::solv::RobinBcCoefStrategy<NDIM>* IBAMR::MarangoniSurfaceTensionForceFunction::d_T_bc_coef = nullptr
private

Boundary condition object for temperature.

◆ d_compute_marangoni_coef

ComputeMarangoniCoefPtr IBAMR::MarangoniSurfaceTensionForceFunction::d_compute_marangoni_coef = nullptr
private

Call back function and the context to find marangoni coefficient as a function of temperature.

◆ d_compute_marangoni_coef_ctx

void* IBAMR::MarangoniSurfaceTensionForceFunction::d_compute_marangoni_coef_ctx = nullptr
private

◆ d_adv_diff_solver

const AdvDiffHierarchyIntegrator* const IBAMR::SurfaceTensionForceFunction::d_adv_diff_solver
protectedinherited

◆ d_ls_var

const SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > IBAMR::SurfaceTensionForceFunction::d_ls_var
protectedinherited

◆ d_ts_type

TimeSteppingType IBAMR::SurfaceTensionForceFunction::d_ts_type
protectedinherited

◆ d_C_idx

int IBAMR::SurfaceTensionForceFunction::d_C_idx = IBTK::invalid_index
protectedinherited

◆ d_phi_idx

int IBAMR::SurfaceTensionForceFunction::d_phi_idx = IBTK::invalid_index
protectedinherited

◆ d_kernel_fcn

std::string IBAMR::SurfaceTensionForceFunction::d_kernel_fcn
protectedinherited

◆ d_sigma

double IBAMR::SurfaceTensionForceFunction::d_sigma = std::numeric_limits<double>::signaling_NaN()
protectedinherited

◆ d_num_interface_cells

double IBAMR::SurfaceTensionForceFunction::d_num_interface_cells = std::numeric_limits<double>::signaling_NaN()
protectedinherited

◆ d_hier_math_ops

SAMRAI::tbox::Pointer<IBTK::HierarchyMathOps> IBAMR::SurfaceTensionForceFunction::d_hier_math_ops
protectedinherited

◆ d_hier_sc_data_ops

SAMRAI::tbox::Pointer<SAMRAI::math::HierarchySideDataOpsReal<NDIM, double> > IBAMR::SurfaceTensionForceFunction::d_hier_sc_data_ops
protectedinherited

◆ d_mask_surface_tension_force

MaskSurfaceTensionForcePtr IBAMR::SurfaceTensionForceFunction::d_mask_surface_tension_force = nullptr
privateinherited

◆ d_mask_surface_tension_force_ctx

void* IBAMR::SurfaceTensionForceFunction::d_mask_surface_tension_force_ctx = nullptr
privateinherited

◆ d_compute_surface_tension_coef

ComputeSurfaceTensionCoefficientPtr IBAMR::SurfaceTensionForceFunction::d_compute_surface_tension_coef = nullptr
privateinherited

◆ d_compute_surface_tension_coef_ctx

void* IBAMR::SurfaceTensionForceFunction::d_compute_surface_tension_coef_ctx = nullptr
privateinherited

◆ d_object_name

std::string IBTK::CartGridFunction::d_object_name
protectedinherited

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