IBAMR  IBAMR version 0.19.
Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
IBAMR::AllenCahnHierarchyIntegrator Class Reference

Class AllenCahnHierarchyIntegrator is a concrete class that manages the time integration of the Allen-Cahn and energy equation. More...

#include <ibamr/AllenCahnHierarchyIntegrator.h>

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

Public Member Functions

 AllenCahnHierarchyIntegrator (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, bool register_for_restart=true)
 
 ~AllenCahnHierarchyIntegrator ()=default
 
void initializeHierarchyIntegrator (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) override
 
void preprocessIntegrateHierarchy (double current_time, double new_time, int num_cycles=1) override
 
void postprocessIntegrateHierarchy (double current_time, double new_time, bool skip_synchronize_new_state_data, int num_cycles=1) override
 
void registerLiquidFractionVariable (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > lf_var, const bool output_lf_var=true) override
 Register liquid fraction variable \( \varphi \). More...
 
void addTemporalAndLinearTermstoRHSOfEnergyEquation (int F_scratch_idx, const double dt) override
 
void computeDivergenceVelocitySourceTerm (int Div_U_F_idx, const double new_time) override
 
void setLiquidFractionPhysicalBcCoef (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > lf_var, SAMRAI::solv::RobinBcCoefStrategy< NDIM > *lf_bc_coef)
 
SAMRAI::solv::RobinBcCoefStrategy< NDIM > * getLiquidFractionPhysicalBcCoef ()
 
SAMRAI::tbox::Pointer< CellConvectiveOperatorgetAllenCahnEquationConvectiveOperator (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > lf_var, SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > H_var)
 
void putToDatabaseSpecialized (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override
 

Protected Member Functions

void integrateHierarchySpecialized (double current_time, double new_time, int cycle_num=0) override
 
virtual void regridHierarchyBeginSpecialized () override
 
virtual void regridHierarchyEndSpecialized () override
 

Private Member Functions

 AllenCahnHierarchyIntegrator ()=delete
 Default constructor. More...
 
 AllenCahnHierarchyIntegrator (const AllenCahnHierarchyIntegrator &from)=delete
 Copy constructor. More...
 
AllenCahnHierarchyIntegratoroperator= (const AllenCahnHierarchyIntegrator &that)=delete
 Assignment operator. More...
 
void computeDoubleWellPotential (int g_firstder_idx, int g_secondder_idx, const int liquid_fraction_idx)
 Compute the double-well potential on the cell-centers based on the liquid fraction value. More...
 
void computeInterpolationFunction (int p_firstder_idx, const int liquid_fraction_idx, const int T_idx)
 Compute an interpolation function on the cell-centers based on the liquid fraction value. More...
 
void computeLiquidFractionSourceTerm (int F_scratch_idx)
 
void getFromInput (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, bool is_from_restart)
 
void getFromRestart ()
 
SAMRAI::tbox::Pointer< IBTK::PoissonSolvergetAllenCahnEquationHelmholtzSolver (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > lf_var)
 
SAMRAI::tbox::Pointer< IBTK::LaplaceOperatorgetAllenCahnEquationHelmholtzRHSOperator (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > lf_var)
 

Private Attributes

SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_lf_F_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::SideVariable< NDIM, double > > d_lf_diffusion_coef_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::SideVariable< NDIM, double > > d_lf_diffusion_coef_rhs_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_lf_H_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_lf_N_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_lf_N_old_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_lf_rhs_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_lf_C_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_lf_temp_rhs_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_g_firstder_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_g_secondder_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_q_firstder_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_chemical_potential_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::SideVariable< NDIM, double > > d_grad_lf_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_Div_u_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_T_lf_N_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > d_lf_interp_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > d_H_interp_var
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > d_lf_flux_var
 
SAMRAI::tbox::Pointer< IBTK::CartGridFunctiond_lf_F_fcn
 
SAMRAI::tbox::Pointer< IBTK::CartGridFunctiond_lf_F_init
 
SAMRAI::solv::RobinBcCoefStrategy< NDIM > * d_lf_bc_coef = nullptr
 
ConvectiveDifferencingType d_lf_convective_difference_form = d_default_convective_difference_form
 
std::string d_lf_convective_op_type = d_default_convective_op_type
 
SAMRAI::tbox::Pointer< SAMRAI::tbox::Databased_lf_convective_op_input_db = d_default_convective_op_input_db
 
SAMRAI::tbox::Pointer< CellConvectiveOperatord_lf_convective_op = nullptr
 
bool d_lf_convective_op_needs_init = false
 
std::string d_lf_solver_type = IBTK::CCPoissonSolverManager::UNDEFINED
 
std::string d_lf_precond_type = IBTK::CCPoissonSolverManager::UNDEFINED
 
SAMRAI::tbox::Pointer< SAMRAI::tbox::Databased_lf_solver_db
 
SAMRAI::tbox::Pointer< SAMRAI::tbox::Databased_lf_precond_db
 
SAMRAI::tbox::Pointer< IBTK::LaplaceOperatord_lf_rhs_op = nullptr
 
SAMRAI::tbox::Pointer< IBTK::PoissonSolverd_lf_solver = nullptr
 
bool d_lf_solver_needs_init
 
bool d_lf_rhs_op_needs_init
 
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > d_lf_sol
 
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > d_lf_rhs
 
TimeSteppingType d_lf_diffusion_time_stepping_type = d_default_diffusion_time_stepping_type
 
TimeSteppingType d_lf_convective_time_stepping_type = d_default_convective_time_stepping_type
 
TimeSteppingType d_lf_init_convective_time_stepping_type = d_default_convective_time_stepping_type
 
bool d_solve_energy = false
 
int d_lf_F_current_idx = IBTK::invalid_index
 
int d_lf_F_scratch_idx = IBTK::invalid_index
 
int d_lf_F_new_idx = IBTK::invalid_index
 
int d_lf_diff_coef_current_idx = IBTK::invalid_index
 
int d_lf_diff_coef_scratch_idx = IBTK::invalid_index
 
int d_lf_diff_coef_new_idx = IBTK::invalid_index
 
int d_lf_N_old_current_idx = IBTK::invalid_index
 
int d_lf_N_old_new_idx = IBTK::invalid_index
 
int d_lf_N_old_scratch_idx = IBTK::invalid_index
 
int d_lf_temp_rhs_idx = IBTK::invalid_index
 
int d_lf_C_idx = IBTK::invalid_index
 
int d_g_firstder_idx = IBTK::invalid_index
 
int d_g_secondder_idx = IBTK::invalid_index
 
int d_q_firstder_idx = IBTK::invalid_index
 
int d_chemical_potential_idx = IBTK::invalid_index
 
int d_grad_lf_idx = IBTK::invalid_index
 
int d_H_sc_idx = IBTK::invalid_index
 
int d_T_lf_N_scratch_idx = IBTK::invalid_index
 
int d_lf_interp_idx = IBTK::invalid_index
 
int d_H_interp_idx = IBTK::invalid_index
 
int d_lf_flux_idx = IBTK::invalid_index
 
int d_lf_diffusion_coef_rhs_scratch_idx = IBTK::invalid_index
 
int d_lf_rhs_scratch_idx = IBTK::invalid_index
 
int d_lf_H_scratch_idx = IBTK::invalid_index
 
int d_Div_u_scratch_idx = IBTK::invalid_index
 
int d_lf_N_scratch_idx = IBTK::invalid_index
 
double d_mobility_lf
 
double d_mixing_energy_density_lf
 
double d_interface_thickness_lf
 
double d_numerical_diffusion = 1.0e-8
 
std::string d_interpolation_function_profile = "LINEAR_1"
 

Detailed Description

In this class, implementations are provided for phase change of a pure material using the phase-field method.

Constructor & Destructor Documentation

◆ AllenCahnHierarchyIntegrator() [1/3]

IBAMR::AllenCahnHierarchyIntegrator::AllenCahnHierarchyIntegrator ( const std::string &  object_name,
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database input_db,
bool  register_for_restart = true 
)

The constructor for class AllenCahnHierarchyIntegrator sets some default values, reads in configuration information from input and restart databases, and registers the integrator object with the restart manager when requested.

◆ ~AllenCahnHierarchyIntegrator()

IBAMR::AllenCahnHierarchyIntegrator::~AllenCahnHierarchyIntegrator ( )
default

The destructor for class AllenCahnHierarchyIntegrator unregisters the integrator object with the restart manager when the object is so registered.

◆ AllenCahnHierarchyIntegrator() [2/3]

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

◆ AllenCahnHierarchyIntegrator() [3/3]

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

Member Function Documentation

◆ initializeHierarchyIntegrator()

void IBAMR::AllenCahnHierarchyIntegrator::initializeHierarchyIntegrator ( SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  hierarchy,
SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > >  gridding_alg 
)
override

Initialize the variables, basic communications algorithms, solvers, and other data structures used by this time integrator object.

This method is called automatically by initializePatchHierarchy() prior to the construction of the patch hierarchy. It is also possible for users to make an explicit call to initializeHierarchyIntegrator() prior to calling initializePatchHierarchy().

◆ preprocessIntegrateHierarchy()

void IBAMR::AllenCahnHierarchyIntegrator::preprocessIntegrateHierarchy ( double  current_time,
double  new_time,
int  num_cycles = 1 
)
override

Prepare to advance the data from current_time to new_time.

◆ postprocessIntegrateHierarchy()

void IBAMR::AllenCahnHierarchyIntegrator::postprocessIntegrateHierarchy ( double  current_time,
double  new_time,
bool  skip_synchronize_new_state_data,
int  num_cycles = 1 
)
override

Clean up data following call(s) to integrateHierarchy().

◆ registerLiquidFractionVariable()

void IBAMR::AllenCahnHierarchyIntegrator::registerLiquidFractionVariable ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  lf_var,
const bool  output_lf_var = true 
)
override

◆ addTemporalAndLinearTermstoRHSOfEnergyEquation()

void IBAMR::AllenCahnHierarchyIntegrator::addTemporalAndLinearTermstoRHSOfEnergyEquation ( int  F_scratch_idx,
const double  dt 
)
override

Add the temporal terms to the RHS of the energy equation.

◆ computeDivergenceVelocitySourceTerm()

void IBAMR::AllenCahnHierarchyIntegrator::computeDivergenceVelocitySourceTerm ( int  Div_U_F_idx,
const double  new_time 
)
override

Compute the source term for the Div U equation.

◆ setLiquidFractionPhysicalBcCoef()

void IBAMR::AllenCahnHierarchyIntegrator::setLiquidFractionPhysicalBcCoef ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  lf_var,
SAMRAI::solv::RobinBcCoefStrategy< NDIM > *  lf_bc_coef 
)

Set an object to provide boundary conditions for \( \varphi \) variable, that has been registered with the hierarchy integrator.

◆ getLiquidFractionPhysicalBcCoef()

SAMRAI::solv::RobinBcCoefStrategy<NDIM>* IBAMR::AllenCahnHierarchyIntegrator::getLiquidFractionPhysicalBcCoef ( )

Return a \( \varphi \) boundary condition object.

◆ getAllenCahnEquationConvectiveOperator()

SAMRAI::tbox::Pointer<CellConvectiveOperator> IBAMR::AllenCahnHierarchyIntegrator::getAllenCahnEquationConvectiveOperator ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  lf_var,
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  H_var 
)

Get the convective operator being used by this solver class for Allen-Cahn equation.

If the convective operator has not already been constructed, then this function will initialize a default convective operator.

◆ putToDatabaseSpecialized()

void IBAMR::AllenCahnHierarchyIntegrator::putToDatabaseSpecialized ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db)
override

Write out specialized object state to the given database.

◆ integrateHierarchySpecialized()

void IBAMR::AllenCahnHierarchyIntegrator::integrateHierarchySpecialized ( double  current_time,
double  new_time,
int  cycle_num = 0 
)
overrideprotected

Synchronously advance each level in the hierarchy over the given time increment.

◆ regridHierarchyBeginSpecialized()

virtual void IBAMR::AllenCahnHierarchyIntegrator::regridHierarchyBeginSpecialized ( )
overrideprotectedvirtual

Reset cached hierarchy dependent data for solvers and operators before the regridding operation.

◆ regridHierarchyEndSpecialized()

virtual void IBAMR::AllenCahnHierarchyIntegrator::regridHierarchyEndSpecialized ( )
overrideprotectedvirtual

Reset cached hierarchy dependent data for solvers and operators before the regridding operation.

◆ operator=()

AllenCahnHierarchyIntegrator& IBAMR::AllenCahnHierarchyIntegrator::operator= ( const AllenCahnHierarchyIntegrator 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.

◆ computeDoubleWellPotential()

void IBAMR::AllenCahnHierarchyIntegrator::computeDoubleWellPotential ( int  g_firstder_idx,
int  g_secondder_idx,
const int  liquid_fraction_idx 
)
private

◆ computeInterpolationFunction()

void IBAMR::AllenCahnHierarchyIntegrator::computeInterpolationFunction ( int  p_firstder_idx,
const int  liquid_fraction_idx,
const int  T_idx 
)
private

◆ computeLiquidFractionSourceTerm()

void IBAMR::AllenCahnHierarchyIntegrator::computeLiquidFractionSourceTerm ( int  F_scratch_idx)
private

◆ getFromInput()

void IBAMR::AllenCahnHierarchyIntegrator::getFromInput ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database input_db,
bool  is_from_restart 
)
private

Read input values from a given database.

◆ getFromRestart()

void IBAMR::AllenCahnHierarchyIntegrator::getFromRestart ( )
private

Read object state from the restart file and initialize class data members.

◆ getAllenCahnEquationHelmholtzSolver()

SAMRAI::tbox::Pointer<IBTK::PoissonSolver> IBAMR::AllenCahnHierarchyIntegrator::getAllenCahnEquationHelmholtzSolver ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  lf_var)
private

Get the solver for the Allen-Cahn equation.

◆ getAllenCahnEquationHelmholtzRHSOperator()

SAMRAI::tbox::Pointer<IBTK::LaplaceOperator> IBAMR::AllenCahnHierarchyIntegrator::getAllenCahnEquationHelmholtzRHSOperator ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  lf_var)
private

Get the operator to use to evaluate the right-hand side of the Allen-Cahn equation.

Member Data Documentation

◆ d_lf_F_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_F_var
private

Solver variables.

◆ d_lf_diffusion_coef_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::SideVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_diffusion_coef_var
private

◆ d_lf_diffusion_coef_rhs_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::SideVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_diffusion_coef_rhs_var
private

◆ d_lf_H_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_H_var
private

◆ d_lf_N_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_N_var
private

◆ d_lf_N_old_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_N_old_var
private

◆ d_lf_rhs_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_rhs_var
private

◆ d_lf_C_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_C_var
private

◆ d_lf_temp_rhs_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_temp_rhs_var
private

◆ d_g_firstder_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_g_firstder_var
private

◆ d_g_secondder_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_g_secondder_var
private

◆ d_q_firstder_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_q_firstder_var
private

◆ d_chemical_potential_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_chemical_potential_var
private

◆ d_grad_lf_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::SideVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_grad_lf_var
private

◆ d_Div_u_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_Div_u_var
private

◆ d_T_lf_N_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_T_lf_N_var
private

◆ d_lf_interp_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::FaceVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_interp_var
private

◆ d_H_interp_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::FaceVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_H_interp_var
private

◆ d_lf_flux_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::FaceVariable<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_flux_var
private

◆ d_lf_F_fcn

SAMRAI::tbox::Pointer<IBTK::CartGridFunction> IBAMR::AllenCahnHierarchyIntegrator::d_lf_F_fcn
private

Cartgrid functions to be used to set the Allen-Cahn equation source term.

◆ d_lf_F_init

SAMRAI::tbox::Pointer<IBTK::CartGridFunction> IBAMR::AllenCahnHierarchyIntegrator::d_lf_F_init
private

◆ d_lf_bc_coef

SAMRAI::solv::RobinBcCoefStrategy<NDIM>* IBAMR::AllenCahnHierarchyIntegrator::d_lf_bc_coef = nullptr
private

Boundary condition object for \( \varphi \).

◆ d_lf_convective_difference_form

ConvectiveDifferencingType IBAMR::AllenCahnHierarchyIntegrator::d_lf_convective_difference_form = d_default_convective_difference_form
private

Convective operator and difference type for the Allen-Cahn equation.

◆ d_lf_convective_op_type

std::string IBAMR::AllenCahnHierarchyIntegrator::d_lf_convective_op_type = d_default_convective_op_type
private

◆ d_lf_convective_op_input_db

SAMRAI::tbox::Pointer<SAMRAI::tbox::Database> IBAMR::AllenCahnHierarchyIntegrator::d_lf_convective_op_input_db = d_default_convective_op_input_db
private

◆ d_lf_convective_op

SAMRAI::tbox::Pointer<CellConvectiveOperator> IBAMR::AllenCahnHierarchyIntegrator::d_lf_convective_op = nullptr
private

◆ d_lf_convective_op_needs_init

bool IBAMR::AllenCahnHierarchyIntegrator::d_lf_convective_op_needs_init = false
private

◆ d_lf_solver_type

std::string IBAMR::AllenCahnHierarchyIntegrator::d_lf_solver_type = IBTK::CCPoissonSolverManager::UNDEFINED
private

Hierarchy operators, solvers and related data for the Allen-Cahn equation.

◆ d_lf_precond_type

std::string IBAMR::AllenCahnHierarchyIntegrator::d_lf_precond_type = IBTK::CCPoissonSolverManager::UNDEFINED
private

◆ d_lf_solver_db

SAMRAI::tbox::Pointer<SAMRAI::tbox::Database> IBAMR::AllenCahnHierarchyIntegrator::d_lf_solver_db
private

◆ d_lf_precond_db

SAMRAI::tbox::Pointer<SAMRAI::tbox::Database> IBAMR::AllenCahnHierarchyIntegrator::d_lf_precond_db
private

◆ d_lf_rhs_op

SAMRAI::tbox::Pointer<IBTK::LaplaceOperator> IBAMR::AllenCahnHierarchyIntegrator::d_lf_rhs_op = nullptr
private

◆ d_lf_solver

SAMRAI::tbox::Pointer<IBTK::PoissonSolver> IBAMR::AllenCahnHierarchyIntegrator::d_lf_solver = nullptr
private

◆ d_lf_solver_needs_init

bool IBAMR::AllenCahnHierarchyIntegrator::d_lf_solver_needs_init
private

◆ d_lf_rhs_op_needs_init

bool IBAMR::AllenCahnHierarchyIntegrator::d_lf_rhs_op_needs_init
private

◆ d_lf_sol

SAMRAI::tbox::Pointer<SAMRAI::solv::SAMRAIVectorReal<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_sol
private

◆ d_lf_rhs

SAMRAI::tbox::Pointer<SAMRAI::solv::SAMRAIVectorReal<NDIM, double> > IBAMR::AllenCahnHierarchyIntegrator::d_lf_rhs
private

◆ d_lf_diffusion_time_stepping_type

TimeSteppingType IBAMR::AllenCahnHierarchyIntegrator::d_lf_diffusion_time_stepping_type = d_default_diffusion_time_stepping_type
private

◆ d_lf_convective_time_stepping_type

TimeSteppingType IBAMR::AllenCahnHierarchyIntegrator::d_lf_convective_time_stepping_type = d_default_convective_time_stepping_type
private

◆ d_lf_init_convective_time_stepping_type

TimeSteppingType IBAMR::AllenCahnHierarchyIntegrator::d_lf_init_convective_time_stepping_type = d_default_convective_time_stepping_type
private

◆ d_solve_energy

bool IBAMR::AllenCahnHierarchyIntegrator::d_solve_energy = false
private

Boolean to identify whether the energy equation is to be solved.

◆ d_lf_F_current_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_F_current_idx = IBTK::invalid_index
private

Patch data descriptor indices for all "state" variables managed by the integrator.

State variables have three contexts: current, scratch, and new.

◆ d_lf_F_scratch_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_F_scratch_idx = IBTK::invalid_index
private

◆ d_lf_F_new_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_F_new_idx = IBTK::invalid_index
private

◆ d_lf_diff_coef_current_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_diff_coef_current_idx = IBTK::invalid_index
private

◆ d_lf_diff_coef_scratch_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_diff_coef_scratch_idx = IBTK::invalid_index
private

◆ d_lf_diff_coef_new_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_diff_coef_new_idx = IBTK::invalid_index
private

◆ d_lf_N_old_current_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_N_old_current_idx = IBTK::invalid_index
private

◆ d_lf_N_old_new_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_N_old_new_idx = IBTK::invalid_index
private

◆ d_lf_N_old_scratch_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_N_old_scratch_idx = IBTK::invalid_index
private

◆ d_lf_temp_rhs_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_temp_rhs_idx = IBTK::invalid_index
private

Patch data descriptor indices for all "scratch" variables managed by the integrator.

Scratch variables have only one context: scratch.

◆ d_lf_C_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_C_idx = IBTK::invalid_index
private

◆ d_g_firstder_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_g_firstder_idx = IBTK::invalid_index
private

◆ d_g_secondder_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_g_secondder_idx = IBTK::invalid_index
private

◆ d_q_firstder_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_q_firstder_idx = IBTK::invalid_index
private

◆ d_chemical_potential_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_chemical_potential_idx = IBTK::invalid_index
private

◆ d_grad_lf_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_grad_lf_idx = IBTK::invalid_index
private

◆ d_H_sc_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_H_sc_idx = IBTK::invalid_index
private

◆ d_T_lf_N_scratch_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_T_lf_N_scratch_idx = IBTK::invalid_index
private

◆ d_lf_interp_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_interp_idx = IBTK::invalid_index
private

◆ d_H_interp_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_H_interp_idx = IBTK::invalid_index
private

◆ d_lf_flux_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_flux_idx = IBTK::invalid_index
private

◆ d_lf_diffusion_coef_rhs_scratch_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_diffusion_coef_rhs_scratch_idx = IBTK::invalid_index
private

◆ d_lf_rhs_scratch_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_rhs_scratch_idx = IBTK::invalid_index
private

◆ d_lf_H_scratch_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_H_scratch_idx = IBTK::invalid_index
private

◆ d_Div_u_scratch_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_Div_u_scratch_idx = IBTK::invalid_index
private

◆ d_lf_N_scratch_idx

int IBAMR::AllenCahnHierarchyIntegrator::d_lf_N_scratch_idx = IBTK::invalid_index
private

◆ d_mobility_lf

double IBAMR::AllenCahnHierarchyIntegrator::d_mobility_lf
private

Allen-Cahn equation parameters.

◆ d_mixing_energy_density_lf

double IBAMR::AllenCahnHierarchyIntegrator::d_mixing_energy_density_lf
private

◆ d_interface_thickness_lf

double IBAMR::AllenCahnHierarchyIntegrator::d_interface_thickness_lf
private

◆ d_numerical_diffusion

double IBAMR::AllenCahnHierarchyIntegrator::d_numerical_diffusion = 1.0e-8
private

To smoothly extend the liquid fraction to H=0 region.

◆ d_interpolation_function_profile

std::string IBAMR::AllenCahnHierarchyIntegrator::d_interpolation_function_profile = "LINEAR_1"
private

String to identify the profile that we want to use for interpolation function q'.


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