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

Class PoissonUtilities provides utility functions for constructing Poisson solvers. More...

#include </home/runner/work/IBAMR/IBAMR/ibtk/include/ibtk/PoissonUtilities.h>

Static Public Member Functions

static void computeMatrixCoefficients (SAMRAI::pdat::CellData< NDIM, double > &matrix_coefficients, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const std::vector< SAMRAI::hier::Index< NDIM > > &stencil, const SAMRAI::solv::PoissonSpecifications &poisson_spec, SAMRAI::solv::RobinBcCoefStrategy< NDIM > *bc_coef, double data_time)
 
static void computeMatrixCoefficients (SAMRAI::pdat::CellData< NDIM, double > &matrix_coefficients, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const std::vector< SAMRAI::hier::Index< NDIM > > &stencil, const SAMRAI::solv::PoissonSpecifications &poisson_spec, const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &bc_coefs, double data_time)
 
static void computeMatrixCoefficients (SAMRAI::pdat::SideData< NDIM, double > &matrix_coefficients, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const std::vector< SAMRAI::hier::Index< NDIM > > &stencil, const SAMRAI::solv::PoissonSpecifications &poisson_spec, const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &bc_coefs, double data_time)
 
static void computeVCSCViscousOpMatrixCoefficients (SAMRAI::pdat::SideData< NDIM, double > &matrix_coefficients, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const std::vector< std::map< SAMRAI::hier::Index< NDIM >, int, IndexFortranOrder > > &stencil_map_vec, const SAMRAI::solv::PoissonSpecifications &poisson_spec, double alpha, double beta, const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &bc_coefs, double data_time, VCInterpType mu_interp_type=VC_HARMONIC_INTERP)
 
static void adjustRHSAtPhysicalBoundary (SAMRAI::pdat::CellData< NDIM, double > &rhs_data, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const SAMRAI::solv::PoissonSpecifications &poisson_spec, SAMRAI::solv::RobinBcCoefStrategy< NDIM > *bc_coef, double data_time, bool homogeneous_bc)
 
static void adjustRHSAtPhysicalBoundary (SAMRAI::pdat::CellData< NDIM, double > &rhs_data, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const SAMRAI::solv::PoissonSpecifications &poisson_spec, const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &bc_coefs, double data_time, bool homogeneous_bc)
 
static void adjustRHSAtPhysicalBoundary (SAMRAI::pdat::SideData< NDIM, double > &rhs_data, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const SAMRAI::solv::PoissonSpecifications &poisson_spec, const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &bc_coefs, double data_time, bool homogeneous_bc)
 
static void adjustVCSCViscousOpRHSAtPhysicalBoundary (SAMRAI::pdat::SideData< NDIM, double > &rhs_data, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const SAMRAI::solv::PoissonSpecifications &poisson_spec, double alpha, const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &bc_coefs, double data_time, bool homogeneous_bc, VCInterpType mu_interp_type=VC_HARMONIC_INTERP)
 
static void adjustRHSAtCoarseFineBoundary (SAMRAI::pdat::CellData< NDIM, double > &rhs_data, const SAMRAI::pdat::CellData< NDIM, double > &sol_data, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const SAMRAI::solv::PoissonSpecifications &poisson_spec, const SAMRAI::tbox::Array< SAMRAI::hier::BoundaryBox< NDIM > > &type1_cf_bdry)
 
static void adjustRHSAtCoarseFineBoundary (SAMRAI::pdat::SideData< NDIM, double > &rhs_data, const SAMRAI::pdat::SideData< NDIM, double > &sol_data, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const SAMRAI::solv::PoissonSpecifications &poisson_spec, const SAMRAI::tbox::Array< SAMRAI::hier::BoundaryBox< NDIM > > &type1_cf_bdry)
 
static void adjustVCSCViscousOpRHSAtCoarseFineBoundary (SAMRAI::pdat::SideData< NDIM, double > &rhs_data, const SAMRAI::pdat::SideData< NDIM, double > &sol_data, SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > > patch, const SAMRAI::solv::PoissonSpecifications &poisson_spec, double alpha, const SAMRAI::tbox::Array< SAMRAI::hier::BoundaryBox< NDIM > > &type1_cf_bdry, VCInterpType mu_interp_type=VC_HARMONIC_INTERP)
 

Detailed Description

Class PoissonUtilities provides utility functions for constructing Poisson solvers.

Member Function Documentation

◆ adjustRHSAtCoarseFineBoundary() [1/2]

void IBTK::PoissonUtilities::adjustRHSAtCoarseFineBoundary ( SAMRAI::pdat::CellData< NDIM, double > &  rhs_data,
const SAMRAI::pdat::CellData< NDIM, double > &  sol_data,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
const SAMRAI::tbox::Array< SAMRAI::hier::BoundaryBox< NDIM > > &  type1_cf_bdry 
)
static

Modify the right-hand side entries to account for coarse-fine interface boundary conditions corresponding to a cell-centered discretization of the Laplacian.

Note
This function simply uses ghost cell values in sol_data to provide Dirichlet boundary values at coarse-fine interfaces. A more complete implementation would employ the interpolation stencil used at coarse-fine interfaces to modify both the matrix coefficients and RHS values at coarse-fine interfaces.

◆ adjustRHSAtCoarseFineBoundary() [2/2]

void IBTK::PoissonUtilities::adjustRHSAtCoarseFineBoundary ( SAMRAI::pdat::SideData< NDIM, double > &  rhs_data,
const SAMRAI::pdat::SideData< NDIM, double > &  sol_data,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
const SAMRAI::tbox::Array< SAMRAI::hier::BoundaryBox< NDIM > > &  type1_cf_bdry 
)
static

Modify the right-hand side entries to account for coarse-fine interface boundary conditions corresponding to a side-centered discretization of the Laplacian.

Note
This function simply uses ghost cell values in sol_data to provide Dirichlet boundary values at coarse-fine interfaces. A more complete implementation would employ the interpolation stencil used at coarse-fine interfaces to modify both the matrix coefficients and RHS values at coarse-fine interfaces.

◆ adjustRHSAtPhysicalBoundary() [1/3]

void IBTK::PoissonUtilities::adjustRHSAtPhysicalBoundary ( SAMRAI::pdat::CellData< NDIM, double > &  rhs_data,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &  bc_coefs,
double  data_time,
bool  homogeneous_bc 
)
static

Modify the right-hand side entries to account for physical boundary conditions corresponding to a cell-centered discretization of the Laplacian.

◆ adjustRHSAtPhysicalBoundary() [2/3]

void IBTK::PoissonUtilities::adjustRHSAtPhysicalBoundary ( SAMRAI::pdat::CellData< NDIM, double > &  rhs_data,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
SAMRAI::solv::RobinBcCoefStrategy< NDIM > *  bc_coef,
double  data_time,
bool  homogeneous_bc 
)
static

Modify the right-hand side entries to account for physical boundary conditions corresponding to a cell-centered discretization of the Laplacian.

◆ adjustRHSAtPhysicalBoundary() [3/3]

void IBTK::PoissonUtilities::adjustRHSAtPhysicalBoundary ( SAMRAI::pdat::SideData< NDIM, double > &  rhs_data,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &  bc_coefs,
double  data_time,
bool  homogeneous_bc 
)
static

Modify the right-hand side entries to account for physical boundary conditions corresponding to a side-centered discretization of the Laplacian.

◆ adjustVCSCViscousOpRHSAtCoarseFineBoundary()

void IBTK::PoissonUtilities::adjustVCSCViscousOpRHSAtCoarseFineBoundary ( SAMRAI::pdat::SideData< NDIM, double > &  rhs_data,
const SAMRAI::pdat::SideData< NDIM, double > &  sol_data,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
double  alpha,
const SAMRAI::tbox::Array< SAMRAI::hier::BoundaryBox< NDIM > > &  type1_cf_bdry,
VCInterpType  mu_interp_type = VC_HARMONIC_INTERP 
)
static

Modify the right-hand side entries to account for coarse-fine interface boundary conditions corresponding to a side-centered discretization of the variable coefficient viscous operator.

Note
This function simply uses ghost cell values in sol_data to provide Dirichlet boundary values at coarse-fine interfaces. A more complete implementation would employ the interpolation stencil used at coarse-fine interfaces to modify both the matrix coefficients and RHS values at coarse-fine interfaces.
The scaling factors of $ D $ variable in the PoissonSpecification object is passed separately and is denoted $ \alpha $.

◆ adjustVCSCViscousOpRHSAtPhysicalBoundary()

void IBTK::PoissonUtilities::adjustVCSCViscousOpRHSAtPhysicalBoundary ( SAMRAI::pdat::SideData< NDIM, double > &  rhs_data,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
double  alpha,
const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &  bc_coefs,
double  data_time,
bool  homogeneous_bc,
VCInterpType  mu_interp_type = VC_HARMONIC_INTERP 
)
static

Modify the right-hand side entries to account for physical boundary conditions corresponding to a side-centered discretization of the variable-coefficient viscous operator.

Note
The scaling factors of $ D $ variable in the PoissonSpecification object is passed separately and is denoted $ \alpha $.

◆ computeMatrixCoefficients() [1/3]

void IBTK::PoissonUtilities::computeMatrixCoefficients ( SAMRAI::pdat::CellData< NDIM, double > &  matrix_coefficients,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const std::vector< SAMRAI::hier::Index< NDIM > > &  stencil,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &  bc_coefs,
double  data_time 
)
static

Compute the matrix coefficients corresponding to a cell-centered discretization of the Laplacian.

◆ computeMatrixCoefficients() [2/3]

void IBTK::PoissonUtilities::computeMatrixCoefficients ( SAMRAI::pdat::CellData< NDIM, double > &  matrix_coefficients,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const std::vector< SAMRAI::hier::Index< NDIM > > &  stencil,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
SAMRAI::solv::RobinBcCoefStrategy< NDIM > *  bc_coef,
double  data_time 
)
static

Compute the matrix coefficients corresponding to a cell-centered discretization of the Laplacian.

◆ computeMatrixCoefficients() [3/3]

void IBTK::PoissonUtilities::computeMatrixCoefficients ( SAMRAI::pdat::SideData< NDIM, double > &  matrix_coefficients,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const std::vector< SAMRAI::hier::Index< NDIM > > &  stencil,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &  bc_coefs,
double  data_time 
)
static

Compute the matrix coefficients corresponding to a side-centered discretization of the Laplacian.

◆ computeVCSCViscousOpMatrixCoefficients()

void IBTK::PoissonUtilities::computeVCSCViscousOpMatrixCoefficients ( SAMRAI::pdat::SideData< NDIM, double > &  matrix_coefficients,
SAMRAI::tbox::Pointer< SAMRAI::hier::Patch< NDIM > >  patch,
const std::vector< std::map< SAMRAI::hier::Index< NDIM >, int, IndexFortranOrder > > &  stencil_map_vec,
const SAMRAI::solv::PoissonSpecifications poisson_spec,
double  alpha,
double  beta,
const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &  bc_coefs,
double  data_time,
VCInterpType  mu_interp_type = VC_HARMONIC_INTERP 
)
static

Compute the matrix coefficients corresponding to a side-centered discretization of the divergence of the viscous stress tensor.

Note
The scaling factors of $ C $ and $ D $ variables in the PoissonSpecification object are passed separately and are denoted by $ \beta $ and $ \alpha $, respectively.

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