IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
Classes | Typedefs | Enumerations | Functions
IBTK Namespace Reference

Classes

class  AppInitializer
 Class AppInitializer provides functionality to simplify the initialization code in an application code. More...
 
class  BGaussSeidelPreconditioner
 Class BGaussSeidelPreconditioner is a block Gauss-Seidel preconditioner which implements the abstract LinearSolver interface. More...
 
class  BJacobiPreconditioner
 Class BJacobiPreconditioner is a block Jacobi preconditioner which implements the abstract LinearSolver interface. More...
 
class  CartCellDoubleBoundsPreservingConservativeLinearRefine
 Class CartCellDoubleBoundsPreservingConservativeLinearRefine is a concrete SAMRAI::xfer::RefineOperator object which prolongs cell-centered double precision patch data via conservative linear interpolation with an additional bounds preservation repair step. More...
 
class  CartCellDoubleCubicCoarsen
 Class CartCellDoubleCubicCoarsen is a concrete SAMRAI::xfer::CoarsenOperator for restricting cell-centered double precision patch data via cubic interpolation. More...
 
class  CartCellDoubleLinearCFInterpolation
 Class CartCellDoubleLinearCFInterpolation is a concrete SAMRAI::xfer::RefinePatchStrategy which sets coarse-fine interface ghost cell values for cell-centered double precision patch data via linear interpolation in the normal and tangential directions at coarse-fine interfaces. More...
 
class  CartCellDoubleQuadraticCFInterpolation
 Class CartCellDoubleQuadraticCFInterpolation is a concrete SAMRAI::xfer::RefinePatchStrategy which sets coarse-fine interface ghost cell values for cell-centered double precision patch data via quadratic interpolation in the normal and tangential directions at coarse-fine interfaces. More...
 
class  CartCellDoubleQuadraticRefine
 Class CartCellDoubleQuadraticRefine is a concrete SAMRAI::xfer::RefineOperator object which prolongs cell-centered double precision patch data via quadratic interpolation. More...
 
class  CartCellRobinPhysBdryOp
 Class CartCellRobinPhysBdryOp is a concrete SAMRAI::xfer::RefinePatchStrategy for setting Robin boundary conditions at physical boundaries for cell-centered scalar- and vector-valued quantities. More...
 
class  CartExtrapPhysBdryOp
 Class CartExtrapPhysBdryOp is a concrete SAMRAI::xfer::RefinePatchStrategy for setting ghost cell values at physical boundaries via constant, linear, or quadratic extrapolation from interior values. More...
 
class  CartGridFunction
 Class CartGridFunction provides an abstract interface for objects for evaluating functions to set values in SAMRAI::hier::PatchData objects. More...
 
class  CartGridFunctionSet
 Class CartGridFunctionSet is a concrete CartGridFunction that is used to allow multiple CartGridFunction objects to act as a single function. More...
 
class  CartSideDoubleCubicCoarsen
 Class CartSideDoubleCubicCoarsen is a concrete SAMRAI::xfer::CoarsenOperator for restricting side-centered double precision patch data via cubic interpolation. More...
 
class  CartSideDoubleDivPreservingRefine
 Class CartSideDoubleDivPreservingRefine is a concrete SAMRAI::xfer::RefinePatchStrategy which prolongs side-centered double precision patch data via conservative linear interpolation with divergence- and curl-preserving corrections. More...
 
class  CartSideDoubleQuadraticCFInterpolation
 Class CartSideDoubleQuadraticCFInterpolation is a concrete SAMRAI::xfer::RefinePatchStrategy which sets coarse-fine interface ghost cell values for side-centered double precision patch data via quadratic interpolation in the normal and tangential directions at coarse-fine interfaces. More...
 
class  CartSideDoubleRT0Coarsen
 Class CartSideDoubleRT0Coarsen is a concrete SAMRAI::xfer::CoarsenOperator for restricting side-centered double precision patch data via the adjoint of RT0 interpolation. More...
 
class  CartSideDoubleRT0Refine
 Class CartSideDoubleRT0Refine is a concrete SAMRAI::xfer::RefineOperator object that prolongs side-centered double precision patch data via RT0-based interpolation. More...
 
class  CartSideDoubleSpecializedLinearRefine
 Class CartSideDoubleSpecializedLinearRefine is a concrete SAMRAI::xfer::RefineOperator object that prolongs side-centered double precision patch data via linear interpolation in the normal direction and MC-limited piecewise-linear interpolation in the tangential direction. More...
 
class  CartSideRobinPhysBdryOp
 Class CartSideRobinPhysBdryOp is a concrete SAMRAI::xfer::RefinePatchStrategy for setting Robin boundary conditions at physical boundaries for side-centered vector-valued quantities. More...
 
class  CCLaplaceOperator
 Class CCLaplaceOperator is a concrete LaplaceOperator which implements a globally second-order accurate cell-centered finite difference discretization of a scalar elliptic operator of the form $ L = C I + \nabla \cdot D \nabla$. More...
 
class  CCPoissonBoxRelaxationFACOperator
 Class CCPoissonBoxRelaxationFACOperator is a concrete PoissonFACPreconditionerStrategy for solving elliptic equations of the form $ \mbox{$L u$} = \mbox{$(C I + \nabla \cdot D \nabla) u$} = f $ using a globally second-order accurate cell-centered finite-volume discretization, with support for Robin and periodic boundary conditions. More...
 
class  CCPoissonHypreLevelSolver
 Class CCPoissonHypreLevelSolver is a concrete LinearSolver for solving elliptic equations of the form $ \mbox{$L u$} = \mbox{$(C I + \nabla \cdot D \nabla) u$} = f $ on a single SAMRAI::hier::PatchLevel using hypre. More...
 
class  CCPoissonLevelRelaxationFACOperator
 Class CCPoissonLevelRelaxationFACOperator is a concrete PoissonFACPreconditionerStrategy for solving elliptic equations of the form $ \mbox{$L u$} = \mbox{$(C I + \nabla \cdot D \nabla) u$} = f $ using a globally second-order accurate cell-centered finite-volume discretization, with support for Robin and periodic boundary conditions. More...
 
class  CCPoissonPETScLevelSolver
 Class CCPoissonPETScLevelSolver is a concrete PETScLevelSolver for solving elliptic equations of the form $ \mbox{$L u$} = \mbox{$(C I + \nabla \cdot D \nabla) u$} = f $ on a single SAMRAI::hier::PatchLevel. More...
 
class  CCPoissonPointRelaxationFACOperator
 Class CCPoissonPointRelaxationFACOperator is a concrete PoissonFACPreconditionerStrategy for solving elliptic equations of the form $ \mbox{$L u$} = \mbox{$(C I + \nabla \cdot D \nabla) u$} = f $ using a globally second-order accurate cell-centered finite-volume discretization, with support for Robin and periodic boundary conditions. More...
 
class  CCPoissonSolverManager
 Class CCPoissonSolverManager is a singleton manager class to provide access to generic cell-centered PoissonSolver implementations. More...
 
class  CellNoCornersFillPattern
 Class CellCellNoCornersFillPattern is a concrete implementation of the abstract base class SAMRAI::xfer::VariableFillPattern. It is used to calculate overlaps according to a pattern which limits overlaps to the cell-centered ghost region surrounding a patch, excluding all corners. In 3D, it is also possible to configure this fill pattern object also to exclude all edges. More...
 
class  CoarseFineBoundaryRefinePatchStrategy
 Class CoarseFineBoundaryRefinePatchStrategy is a subclass of the abstract base class SAMRAI::xfer::RefinePatchStrategy that extends the functionality of SAMRAI::xfer::RefinePatchStrategy to facilitate the implementation of coarse-fine interface discretizations. More...
 
class  CoarsenPatchStrategySet
 Class CoarsenPatchStrategySet is a utility class that allows multiple SAMRAI::xfer::CoarsenPatchStrategy objects to be employed by a single SAMRAI::xfer::CoarsenSchedule. More...
 
class  CopyToRootSchedule
 Class CopyToRootSchedule is used to communicate distributed patch data to a unified patch data object on a root MPI process. More...
 
class  CopyToRootTransaction
 Class CopyToRootTransaction is a concrete implementation of the abstract base class SAMRAI::tbox::Transaction. It is used to communicate distributed patch data to a unified patch data object on a root MPI process. More...
 
class  DebuggingUtilities
 Class DebuggingUtilities provides debugging functionality. More...
 
class  EdgeDataSynchronization
 Class EdgeDataSynchronization encapsulates the operations required to "synchronize" edge-centered values defined at patch boundaries. More...
 
class  EdgeSynchCopyFillPattern
 Class EdgeCellSynchCopyFillPattern is a concrete implementation of the abstract base class SAMRAI::xfer::VariableFillPattern. It is used to calculate overlaps according to a pattern which limits overlaps to the edge-centered ghost region surrounding a patch appropriate for "synchronizing" edge-centered values in an axis-by-axis manner at patch boundaries. More...
 
class  ExtendedRobinBcCoefStrategy
 Class ExtendedRobinBcCoefStrategy is a subclass of the abstract base class SAMRAI::solv::RobinBcCoefStrategy that extends the functionality of SAMRAI::solv::RobinBcCoefStrategy to allow for the specification of patch data descriptor indices that are required for filling, and the specification of whether homogeneous or inhomogeneous boundary data should be set. More...
 
class  FaceDataSynchronization
 Class FaceDataSynchronization encapsulates the operations required to "synchronize" face-centered values defined at patch boundaries. More...
 
class  FaceSynchCopyFillPattern
 Class FaceCellSynchCopyFillPattern is a concrete implementation of the abstract base class SAMRAI::xfer::VariableFillPattern. It is used to calculate overlaps according to a pattern which limits overlaps to the face-centered ghost region surrounding a patch appropriate for "synchronizing" face-centered values at patch boundaries. More...
 
class  FACPreconditioner
 Class FACPreconditioner is a concrete LinearSolver for implementing FAC (multilevel multigrid) preconditioners. More...
 
class  FACPreconditionerStrategy
 Class FACPreconditionerStrategy provides an interface for specifying the problem-specific operations needed to implement a specific FAC preconditioner. More...
 
class  FixedSizedStream
 Class FixedSizedStream provides a fixed-size message buffer used by various communication routines. More...
 
class  GeneralOperator
 Class GeneralOperator provides an abstract interface for the specification of general operators to compute $ y=F[x] $ and $ z=F[x]+y $. More...
 
class  GeneralSolver
 Class GeneralSolver provides an abstract interface for the implementation of linear or nonlinear solvers for systems of equations defined on an AMR patch hierarchy. More...
 
class  HierarchyGhostCellInterpolation
 Class HierarchyGhostCellInterpolation encapsulates the operations required to set ghost cell values at physical and coarse-fine boundaries across a range of levels of a locally refined patch hierarchy. More...
 
class  HierarchyIntegrator
 Class HierarchyIntegrator provides an abstract interface for a time integrator for a system of equations defined on an AMR grid hierarchy, along with basic data management for variables defined on that hierarchy. More...
 
class  HierarchyMathOps
 Class HierarchyMathOps provides functionality to perform "composite-grid" mathematical operations on a range of levels in a SAMRAI::hier::PatchHierarchy object. More...
 
struct  IBTK_MPI
 Provides C++ wrapper around MPI routines. More...
 
class  IBTKInit
 Initialization for IBAMR programs. More...
 
struct  IndexOrder
 
class  IndexUtilities
 Class IndexUtilities is a utility class that defines simple functions such as conversion routines between physical coordinates and Cartesian index space. More...
 
class  JacobianOperator
 Class JacobianOperator provides an abstract interface for the specification of general operators to compute Jacobian-vector products, i.e., $ F'[x]v $. More...
 
class  KrylovLinearSolver
 Class KrylovLinearSolver provides an abstract interface for the implementation of Krylov subspace solvers for linear problems of the form $Ax=b$. More...
 
class  KrylovLinearSolverManager
 Class KrylovLinearSolverManager is a singleton manager class to provide access to generic KrylovLinearSolver implementations. More...
 
class  KrylovLinearSolverPoissonSolverInterface
 Class KrylovLinearSolverPoissonSolverInterface provides an interface for KrylovLinearSolvers that are to be used as Poisson solvers. More...
 
class  LaplaceOperator
 Class LaplaceOperator is an abstract base class for a Laplace-type operators. More...
 
class  LData
 Class LData provides storage for a single scalar- or vector-valued Lagrangian quantity. More...
 
class  LDataManager
 Class LDataManager coordinates the irregular distribution of LNode and LData on the patch hierarchy. More...
 
class  LIndexSetData
 Class LIndexSetData is a specialization of the templated class LSetData that is intended to be used with Lagrangian data objects that provide Lagrangian and PETSc indexing information. More...
 
class  LEInteractor
 Class LEInteractor is a utility class that defines several functions to interpolate data from Eulerian grid patches onto Lagrangian meshes and to spread values (not densities) from Lagrangian meshes to Eulerian grid patches. More...
 
class  LibMeshSystemIBVectors
 Class LibMeshSystemIBVectors is a convenience class that manages access to libMesh vectors for the same system defined on multiple parts. It extends the base class LibMeshSystemVectors to provide access to vectors ghosted with both the Lagrangian partitioning (i.e., libMesh's computed partitioning, as in the base class) as well as the IB partitioning (i.e., the partitioning based on the distribution of SAMRAI data). More...
 
class  LIndexSetDataFactory
 Class LIndexSetPatchDataFactory provides a SAMRAI::hier::PatchDataFactory class corresponding to patch data of type LIndexSetData. More...
 
class  LIndexSetVariable
 Class LIndexSetVariable provides a SAMRAI::hier::Variable class corresponding to patch data of type LIndexSetData. More...
 
class  LinearOperator
 Class LinearOperator provides an abstract interface for the specification of linear operators to compute $ y=Ax $ and $ z=Ax+y $ and, optionally, $ y=A^{T} x $ and $ z=A^{T}x+y $. More...
 
class  LinearSolver
 Class LinearSolver provides an abstract interface for the implementation of solvers for linear problems of the form $Ax=b$. More...
 
class  LInitStrategy
 Class LInitStrategy provides a mechanism for specifying the initial configuration of the curvilinear mesh. More...
 
class  LMesh
 Class LMesh is a collection of LNode objects. More...
 
class  LNode
 Class LNode is the basic element of an LMesh. More...
 
class  LNodeIndex
 Class LNodeIndex provides Lagrangian and PETSc indexing information for a single node of a Lagrangian mesh. More...
 
class  LNodeIndexPosnComp
 Comparison functor to order on the physical location of the Lagrangian node. More...
 
struct  LNodeIndexLagrangianIndexComp
 Comparison functor to order on the Lagrangian index of the Lagrangian node. More...
 
struct  LNodeIndexGlobalPETScIndexComp
 Comparison functor to order on the global PETSc index of the Lagrangian node. More...
 
struct  LNodeIndexLocalPETScIndexComp
 Comparison functor to order on the local PETSc index of the Lagrangian node. More...
 
class  LNodeIndexPosnEqual
 Comparison functor to check for equality between LNodeIndex objects based on their positions. More...
 
struct  LNodeIndexLagrangianIndexEqual
 Comparison functor to check for equality between LNodeIndex objects based on their Lagrangian indices. More...
 
struct  LNodeIndexGlobalPETScIndexEqual
 Comparison functor to check for equality between between LNodeIndex objects based on their global PETSc indices. More...
 
struct  LNodeIndexLocalPETScIndexEqual
 Comparison functor to check for equality between LNodeIndex objects based on their local PETSc indices. More...
 
class  LSet
 Class LSet provides inter-processor communications and database access functionality to a collection of Lagrangian objects. More...
 
class  LSetData
 Class LSetData is a specialization of the templated class SAMRAI::pdat::IndexData that provides access to Lagrangian objects that are embedded in the a Cartesian grid patch. More...
 
class  LSetDataFactory
 Class LSetPatchDataFactory provides a SAMRAI::hier::PatchDataFactory class corresponding to patch data of type LSetData. More...
 
class  LSetDataIterator
 Class LSetDataIterator is an iterator class which may be used to iterate through LSet objects associated with a specified box in cell-centered index space. More...
 
class  LSetVariable
 Class LSetVariable provides a SAMRAI::hier::Variable class corresponding to patch data of type LSetData. More...
 
class  LSiloDataWriter
 Class LSiloDataWriter provides functionality to output Lagrangian data for visualization via the VisIt visualization tool in the Silo data format. More...
 
class  LTransaction
 Class LTransaction represents a communication transaction between two processors or a local data copy for communicating or copying Lagrangian objects. More...
 
class  MarkerPatch
 
class  MarkerPatchHierarchy
 
class  MergingLoadBalancer
 Class MergingLoadBalancer merges the boxes generated by a load balancer in a final step to decrease the total number of boxes. In essence, it postprocesses the list of boxes generated by its parent class to try and coalesce the set of boxes on each process. More...
 
class  muParserCartGridFunction
 Class muParserCartGridFunction is an implementation of the strategy class CartGridFunction that allows for the run-time specification of (possibly spatially- and temporally-varying) functions which are used to set double precision values on standard SAMRAI SAMRAI::hier::PatchData objects. More...
 
class  muParserRobinBcCoefs
 Class muParserRobinBcCoefs is an implementation of the strategy class SAMRAI::solv::RobinBcCoefStrategy that allows for the run-time specification of (possibly spatially- and temporally-varying) Robin boundary conditions. More...
 
class  NewtonKrylovSolver
 Class NewtonKrylovSolver provides an abstract interface for the implementation of inexact Newton-Krylov solvers for nonlinear problems of the form $ F[x]=b $. More...
 
class  NewtonKrylovSolverManager
 Class NewtonKrylovSolverManager is a singleton manager class to provide access to generic NewtonKrylovSolver implementations. More...
 
class  NodeDataSynchronization
 Class NodeDataSynchronization encapsulates the operations required to "synchronize" node-centered values defined at patch boundaries. More...
 
class  NodeSynchCopyFillPattern
 Class NodeCellSynchCopyFillPattern is a concrete implementation of the abstract base class SAMRAI::xfer::VariableFillPattern. It is used to calculate overlaps according to a pattern which limits overlaps to the node-centered ghost region surrounding a patch appropriate for "synchronizing" node-centered values in an axis-by-axis manner at patch boundaries. More...
 
class  NormOps
 Class NormOps provides functionality for computing discrete vector norms. More...
 
class  ParallelEdgeMap
 Class ParallelEdgeMap is a utility class for managing edge maps (i.e., maps from vertices to links between vertices) in parallel. More...
 
class  ParallelMap
 Class ParallelMap is a utility class for associating integer keys with arbitrary data items in parallel. More...
 
class  ParallelSet
 Class ParallelSet is a utility class for storing collections of integer keys in parallel. More...
 
class  PartitioningBox
 Class PartitioningBox implements an NDIM-dimensional bounding box defined by two points. Unlike a standard bounding box, a PartitioningBox is an NDIM-dimensional tensor product of half-open intervals: i.e., it is a half-open box. This property allows one to partition a domain into a set of boxes. More...
 
class  PartitioningBoxes
 Class PartitioningBoxes stores a set of bounding boxes and can check if a point is in the set of partitioning boxes or not in a more optimized way than just looping over a std::vector<PartitioningBox>. More...
 
class  PatchMathOps
 Class PatchMathOps provides functionality to perform mathematical operations on individual SAMRAI::hier::Patch objects. More...
 
class  PETScKrylovLinearSolver
 Class PETScKrylovLinearSolver provides a KrylovLinearSolver interface for a PETSc Krylov subspace iterative linear solver (KSP). More...
 
class  PETScKrylovPoissonSolver
 Class PETScKrylovPoissonSolver is an extension of class PETScKrylovLinearSolver that provides an implementation of the PoissonSolver interface. More...
 
class  PETScLevelSolver
 Class PETScLevelSolver is an abstract LinearSolver for solving systems of linear equations on a single SAMRAI::hier::PatchLevel using PETSc. More...
 
class  PETScMatLOWrapper
 Class PETScMatLOWrapper provides a LinearOperator interface for a PETSc Mat object. More...
 
class  PETScMatUtilities
 Class PETScMatUtilities provides utility functions for PETSc Mat objects. More...
 
class  PETScMFFDJacobianOperator
 Class PETScMFFDJacobianOperator provides a method for computing Jacobian-vector products, i.e., $ F'[x]v $, via a matrix-free finite-difference approach. More...
 
class  PETScNewtonKrylovSolver
 Class PETScNewtonKrylovSolver provides a NewtonKrylovSolver interface for a PETSc inexact Newton-Krylov iterative nonlinear solver (SNES). More...
 
class  PETScPCLSWrapper
 Class PETScPCLSWrapper provides a LinearSolver interface for a PETSc PC object. More...
 
class  PETScSAMRAIVectorReal
 Class PETScSAMRAIVectorReal is a class for wrapping SAMRAI::solv::SAMRAIVectorReal objects so that they may be used with the PETSc solver package. More...
 
class  PETScSNESFunctionGOWrapper
 Class PETScSNESFunctionGOWrapper provides a GeneralOperator interface for a PETSc SNES nonlinear function. More...
 
class  PETScSNESJacobianJOWrapper
 Class PETScSNESJacobianJOWrapper provides a JacobianOperator interface for a PETSc SNES Jacobian. More...
 
class  PETScVecUtilities
 Class PETScVecUtilities provides utility functions for PETSc Vec objects. More...
 
class  PhysicalBoundaryUtilities
 Class PhysicalBoundaryUtilities is a utility class to organize functions related to setting values at physical boundaries. More...
 
class  PoissonFACPreconditioner
 Class PoissonFACPreconditioner is a FACPreconditioner that has been specialized for Poisson problems. More...
 
class  PoissonFACPreconditionerStrategy
 Class PoissonFACPreconditionerStrategy is an abstract FACPreconditionerStrategy implementing many of the operations required by smoothers for the Poisson equation and related problems. More...
 
class  PoissonSolver
 Class PoissonSolver is an abstract base class for Poisson solvers. More...
 
class  PoissonUtilities
 Class PoissonUtilities provides utility functions for constructing Poisson solvers. More...
 
class  RefinePatchStrategySet
 Class RefinePatchStrategySet is a utility class that allows multiple SAMRAI::xfer::RefinePatchStrategy objects to be employed by a single SAMRAI::xfer::RefineSchedule. More...
 
class  RobinPhysBdryPatchStrategy
 Class RobinPhysBdryPatchStrategy is an abstract strategy class that extends the SAMRAI::xfer::RefinePatchStrategy interface to provide routines specific for setting Robin boundary conditions at physical boundaries. This class also provides default implementations of some methods defined in SAMRAI::xfer::RefinePatchStrategy that are generally not needed for filling ghost cell values at physical boundaries. More...
 
class  SAMRAIDataCache
 Class SAMRAIDataCache is a utility class for caching cloned SAMRAI patch data. Patch data are allocated as needed and should not be deallocated by the caller. More...
 
class  SAMRAIGhostDataAccumulator
 Class that can accumulate data summed into ghost regions on a patch hierarchy into their correct locations. More...
 
class  SAMRAIScopedVectorCopy
 
class  SAMRAIScopedVectorDuplicate
 
class  SCLaplaceOperator
 Class SCLaplaceOperator is a concrete LaplaceOperator which implements a globally second-order accurate side-centered finite difference discretization of a scalar elliptic operator of the form $ L = C I + \nabla \cdot D \nabla$. More...
 
class  SCPoissonHypreLevelSolver
 Class SCPoissonHypreLevelSolver is a concrete LinearSolver for solving elliptic equations of the form $ \mbox{$L u$} = \mbox{$(C I + \nabla \cdot D \nabla) u$} = f $ on a single SAMRAI::hier::PatchLevel using hypre. More...
 
class  SCPoissonPETScLevelSolver
 Class SCPoissonPETScLevelSolver is a concrete PETScLevelSolver for solving elliptic equations of the form $ \mbox{$L u$} = \mbox{$(C I + \nabla \cdot D \nabla) u$} = f $ on a single SAMRAI::hier::PatchLevel using PETSc. More...
 
class  SCPoissonPointRelaxationFACOperator
 Class SCPoissonPointRelaxationFACOperator is a concrete PoissonFACPreconditionerStrategy for solving elliptic equations of the form $ \mbox{$L u$} = \mbox{$(C I + \nabla \cdot D \nabla) u$} = f $ using a globally second-order accurate side-centered finite-difference discretization, with support for Robin and periodic boundary conditions. More...
 
class  SCPoissonSolverManager
 Class SCPoissonSolverManager is a singleton manager class to provide access to generic side-centered PoissonSolver implementations. More...
 
class  SecondaryHierarchy
 
class  SideDataSynchronization
 Class SideDataSynchronization encapsulates the operations required to "synchronize" side-centered values defined at patch boundaries. More...
 
class  SideNoCornersFillPattern
 Class SideCellNoCornersFillPattern is a concrete implementation of the abstract base class SAMRAI::xfer::VariableFillPattern. It is used to calculate overlaps according to a pattern which limits overlaps to the cell-centered ghost region surrounding a patch, excluding all corners. In 3D, it is also possible to configure this fill pattern object also to exclude all edges. More...
 
class  SideSynchCopyFillPattern
 Class SideCellSynchCopyFillPattern is a concrete implementation of the abstract base class SAMRAI::xfer::VariableFillPattern. It is used to calculate overlaps according to a pattern which limits overlaps to the side-centered ghost region surrounding a patch appropriate for "synchronizing" side-centered values at patch boundaries. More...
 
class  SnapshotCache
 Class SnapshotCache provides a method of storing snapshots of patch data on a snapshot of the patch hierarchy. More...
 
class  StaggeredPhysicalBoundaryHelper
 Class StaggeredPhysicalBoundaryHelper provides helper functions to handle physical boundary conditions for a staggered grid discretizations. More...
 
class  StandardTagAndInitStrategySet
 Class StandardTagAndInitStrategySet is a utility class that allows multiple SAMRAI::mesh::StandardTagAndInitStrategy objects to be employed by a single SAMRAI::mesh::StandardTagAndInitialize object. More...
 
class  Streamable
 Class Streamable is an abstract interface for objects that can be packed into SAMRAI::tbox::AbstractStream data streams. More...
 
class  StreamableFactory
 Class StreamableFactory is an abstract interface for classes that can unpack particular concrete Streamable objects from SAMRAI::tbox::AbstractStream data streams. More...
 
class  StreamableManager
 Class StreamableManager is a singleton manager class that organizes the actual packing and unpacking of concrete Streamable objects for SAMRAI::tbox::AbstractStream based communication. More...
 
class  VCSCViscousOperator
 Class VCSCViscousOperator is a subclass of SCLaplaceOperator which implements a globally second-order accurate side-centered finite difference discretization of a vector elliptic operator of the form $ L = \beta C I + \alpha \nabla \cdot \mu ( (\nabla u) + (\nabla u)^T ) $. More...
 
class  VCSCViscousOpPointRelaxationFACOperator
 Class VCSCViscousOpPointRelaxationFACOperator is a specialization of SCPoissonPointRelaxationFACOperator for solving vector elliptic equation of the form $ \mbox{$L u$} = C u + \nabla \cdot \mu (\nabla u + (\nabla u)^T) = f $ using a globally second-order accurate side-centered finite-difference discretization, with support for Robin and periodic boundary conditions. More...
 
class  VCSCViscousPETScLevelSolver
 Class VCSCViscousPETScLevelSolver is a subclass of SCPoissonPETScLevelSolver class which solves vector-valued elliptic equation of the form $ \mbox{$L u$} = C u + \nabla \cdot \mu (\nabla u + (\nabla u)^T) = f $ on a single SAMRAI::hier::PatchLevel using PETSc. More...
 

Typedefs

template<typename T >
using EigenAlignedVector = std::vector< T, Eigen::aligned_allocator< T > >
 
using Matrix2d = Eigen::Matrix< double, 2, 2 >
 
using Vector2d = Eigen::Matrix< double, 2, 1 >
 
using ColumnVector2d = Eigen::Matrix< double, 2, 1 >
 
using RowVector2d = Eigen::Matrix< double, 1, 2 >
 
using Matrix3d = Eigen::Matrix< double, 3, 3 >
 
using Vector3d = Eigen::Matrix< double, 3, 1 >
 
using ColumnVector3d = Eigen::Matrix< double, 3, 1 >
 
using RowVector3d = Eigen::Matrix< double, 1, 3 >
 
using MatrixNd = Eigen::Matrix< double, NDIM, NDIM >
 
using VectorNd = Eigen::Matrix< double, NDIM, 1 >
 
using ColumnVectorNd = Eigen::Matrix< double, NDIM, 1 >
 
using RowVectorNd = Eigen::Matrix< double, 1, NDIM >
 
using MatrixXd = Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic >
 
using VectorXd = Eigen::Matrix< double, Eigen::Dynamic, 1 >
 
using ColumnVectorXd = Eigen::Matrix< double, Eigen::Dynamic, 1 >
 
using RowVectorXd = Eigen::Matrix< double, 1, Eigen::Dynamic >
 
using Matrix = MatrixNd
 
using Point = VectorNd
 
using Vector = VectorNd
 
using RigidDOFVector = Eigen::Matrix< double, s_max_free_dofs, 1 >
 
using FreeRigidDOFVector = Eigen::Matrix< int, s_max_free_dofs, 1 >
 
using RDV = RigidDOFVector
 
using FRDV = FreeRigidDOFVector
 
using IndexFortranOrder = struct IndexOrder< SAMRAI::hier::Index< NDIM > >
 
using CellIndexFortranOrder = struct IndexOrder< SAMRAI::pdat::CellIndex< NDIM > >
 
using LNodeIndexSet = LSet< LNodeIndex >
 
using LNodeIndexSetData = LIndexSetData< LNodeIndex >
 
using LNodeIndexSetDataFactory = LIndexSetDataFactory< LNodeIndex >
 
using LNodeIndexSetDataIterator = LSetDataIterator< LNodeIndex >
 
using LNodeIndexSetVariable = LIndexSetVariable< LNodeIndex >
 
using LNodeIndexTransaction = LTransaction< LNodeIndex >
 
using LNodeSet = LSet< LNode >
 
using LNodeSetData = LIndexSetData< LNode >
 
using LNodeSetDataFactory = LIndexSetDataFactory< LNode >
 
using LNodeSetDataIterator = LSetDataIterator< LNode >
 
using LNodeSetVariable = LIndexSetVariable< LNode >
 
using LNodeTransaction = LTransaction< LNode >
 

Enumerations

enum  MGCycleType {
  F_CYCLE , FMG_CYCLE , V_CYCLE , W_CYCLE ,
  UNKNOWN_MG_CYCLE_TYPE = -1
}
 Enumerated type for different multigrid cycle types.
 
enum  RegridMode { STANDARD , AGGRESSIVE , UNKNOWN_REGRID_MODE = -1 }
 Enumerated type for different regridding modes.
 
enum  VariableContextType { CURRENT_DATA , NEW_DATA , SCRATCH_DATA , UNKNOWN_VARIABLE_CONTEXT_TYPE = -1 }
 Enumerated type for different standard data contexts.
 
enum  VCInterpType { VC_AVERAGE_INTERP = 1 , VC_HARMONIC_INTERP = 2 , UNKNOWN_VC_INTERP_TYPE = -1 }
 Enumerated type for different interpolation types for the material properties of the viscous solver.
 
enum  NodeOutsidePatchCheckType { NODE_OUTSIDE_PERMIT = 1 , NODE_OUTSIDE_WARN = 2 , NODE_OUTSIDE_ERROR = 3 , UNKNOWN_NODE_OUTSIDE_PATCH_CHECK_TYPE = -1 }
 

Functions

std::vector< SAMRAI::hier::Box< NDIM > > merge_boxes_by_longest_edge (const std::vector< SAMRAI::hier::Box< NDIM > > &boxes)
 
template<typename T >
string_to_enum (const std::string &)
 Routine for converting strings to enums.
 
template<typename T >
std::string enum_to_string (T)
 Routine for converting enums to strings.
 
template<>
MGCycleType string_to_enum< MGCycleType > (const std::string &val)
 
template<>
std::string enum_to_string< MGCycleType > (MGCycleType val)
 
template<>
RegridMode string_to_enum< RegridMode > (const std::string &val)
 
template<>
std::string enum_to_string< RegridMode > (RegridMode val)
 
template<>
VariableContextType string_to_enum< VariableContextType > (const std::string &val)
 
template<>
std::string enum_to_string< VariableContextType > (VariableContextType val)
 
template<>
VCInterpType string_to_enum< VCInterpType > (const std::string &val)
 
template<>
std::string enum_to_string< VCInterpType > (VCInterpType val)
 
template<>
NodeOutsidePatchCheckType string_to_enum< NodeOutsidePatchCheckType > (const std::string &val)
 
template<>
std::string enum_to_string< NodeOutsidePatchCheckType > (NodeOutsidePatchCheckType val)
 
MPI_Datatype mpi_type_id (const double)
 
MPI_Datatype mpi_type_id (const std::pair< double, int >)
 
MPI_Datatype mpi_type_id (const int)
 
MPI_Datatype mpi_type_id (const std::pair< int, int >)
 
MPI_Datatype mpi_type_id (const float)
 
MPI_Datatype mpi_type_id (const std::pair< float, int >)
 
MPI_Datatype mpi_type_id (const char)
 
MPI_Datatype mpi_type_id (const unsigned int)
 
template<typename T >
MPI_Datatype mpi_type_id (const T &)
 
bool rel_equal_eps (double a, double b, double eps=std::sqrt(std::numeric_limits< double >::epsilon()))
 
bool abs_equal_eps (double a, double b, double eps=std::sqrt(std::numeric_limits< double >::epsilon()))
 Check whether the absolute difference between a and b are within the threshold eps. More...
 
std::string get_data_time_str (const double data_time, const double current_time, const double new_time)
 
double get_min_patch_dx (const SAMRAI::hier::PatchLevel< NDIM > &patch_level)
 
template<class T , unsigned N>
std::array< T, N > array_constant (const T &v)
 
template<class T , unsigned N>
std::array< T, N > array_one ()
 
template<class T , unsigned N>
std::array< T, N > array_zero ()
 
bool level_can_be_refined (int level_number, int max_levels)
 
std::pair< int, int > voigt_to_tensor_idx (const int k)
 
int tensor_idx_to_voigt (const std::pair< int, int > &idx)
 
double smooth_heaviside (const double &phi, const double &alpha)
 
double smooth_delta (const double &phi, const double &alpha)
 
double discontinuous_heaviside (const double &phi)
 
void deallocate_vector_data (SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number)
 
void free_vector_components (SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number)
 
template<class T >
const T & checked_dereference (const SAMRAI::tbox::Pointer< T > &p)
 
template<class T >
T & checked_dereference (SAMRAI::tbox::Pointer< T > &p)
 
Vec setup_petsc_vector (const unsigned int num_local_nodes, const unsigned int depth, const std::vector< int > &nonlocal_petsc_indices)
 
void update_snapshot (SnapshotCache &cache, int u_idx, double time, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > current_hierarchy, double tol=1.0e-8)
 
void fill_snapshot_on_hierarchy (SnapshotCache &cache, int u_idx, double time, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > current_hierarchy, const std::string &snapshot_refine_type, double tol=1.0e-8)
 
void fill_snapshot_at_time (SnapshotCache &snapshot_cache, int u_idx, double t, int scr_idx, SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, const std::string &refine_type, SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyDataOpsReal< NDIM, double > > hier_data_ops=nullptr, double period=std::numeric_limits< double >::quiet_NaN())
 
void reportPETScKSPConvergedReason (const std::string &object_name, const KSPConvergedReason &reason, std::ostream &os)
 Report the KSPConvergedReason.
 
void reportPETScSNESConvergedReason (const std::string &object_name, const SNESConvergedReason &reason, std::ostream &os)
 Report the SNESConvergedReason.
 
std::array< HYPRE_Int, NDIM > hypre_array (const SAMRAI::hier::Index< NDIM > &index)
 Helper function to convert SAMRAI indices to Hypre integers. More...
 
void copyFromHypre (SAMRAI::pdat::CellData< NDIM, double > &dst_data, const std::vector< HYPRE_StructVector > &vectors, const SAMRAI::hier::Box< NDIM > &box)
 Copy data from a vector of Hypre vectors to SAMRAI cell centered data with depth equal to number of Hypre vectors. More...
 
void copyFromHypre (SAMRAI::pdat::SideData< NDIM, double > &dst_data, HYPRE_SStructVector vector, const SAMRAI::hier::Box< NDIM > &box)
 Copy data from a Hypre vector to SAMRAI side centered data. More...
 
void copyToHypre (const std::vector< HYPRE_StructVector > &vectors, SAMRAI::pdat::CellData< NDIM, double > &src_data, const SAMRAI::hier::Box< NDIM > &box)
 Copy data from SAMRAI cell centered data to Hypre vectors. More...
 
void copyToHypre (HYPRE_SStructVector &vector, SAMRAI::pdat::SideData< NDIM, double > &src_data, const SAMRAI::hier::Box< NDIM > &box)
 Copy data from SAMRAI side centered data to a Hypre vector. More...
 
int add_42 (int a)
 
template<int n_vars, int n_basis, bool weights_are_unity = false>
void sum_weighted_elem_solution_n_vars_n_basis (const int qp_offset, const std::vector< std::vector< double > > &phi_F, const std::vector< double > &weights, const boost::multi_array< double, 2 > &F_node, std::vector< double > &F_w_qp)
 Compute the product of the finite element representation of the force and (optionally) a set of weights (e.g., JxW values) at all quadrature points on an element. See integrate_elem_rhs for details on the first two template arguments. More...
 
template<int n_vars, bool weights_are_unity = false>
void sum_weighted_elem_solution_n_vars (const int n_basis, const int qp_offset, const std::vector< std::vector< double > > &phi_F, const std::vector< double > &weights, const boost::multi_array< double, 2 > &F_node, std::vector< double > &F_w_qp)
 
template<bool weights_are_unity = false>
void sum_weighted_elem_solution (const int n_vars, const int n_basis, const int qp_offset, const std::vector< std::vector< double > > &phi_F, const std::vector< double > &weights, const boost::multi_array< double, 2 > &F_node, std::vector< double > &F_w_qp)
 
template<int n_vars, int n_basis>
void integrate_elem_rhs_n_vars_n_basis (const int qp_offset, const std::vector< std::vector< double > > &phi_F, const std::vector< double > &JxW_F, const std::vector< double > &F_qp, std::vector< double > &F_rhs_concatenated)
 Assemble the element contribution to a load vector. More...
 
template<int n_vars>
void integrate_elem_rhs_n_vars (const int n_basis, const int qp_offset, const std::vector< std::vector< double > > &phi_F, const std::vector< double > &JxW_F, const std::vector< double > &F_qp, std::vector< double > &F_rhs_concatenated)
 
void integrate_elem_rhs (const int n_vars, const int n_basis, const int qp_offset, const std::vector< std::vector< double > > &phi_F, const std::vector< double > &JxW_F, const std::vector< double > &F_qp, std::vector< double > &F_rhs_concatenated)
 
template<int dim, int spacedim>
Eigen::Matrix< double, spacedim, dim > getCovariant (const Eigen::Matrix< double, spacedim, dim > &contravariant)
 
template<int dim>
Eigen::Matrix< double, dim, dim > getCovariant (const Eigen::Matrix< double, dim, dim > &contravariant)
 
void setup_system_vectors (libMesh::EquationSystems *equation_systems, const std::vector< std::string > &system_names, const std::vector< std::string > &vector_names, const bool from_restart)
 
void setup_system_vector (libMesh::System &system, const std::string &vector_name, const bool from_restart)
 
void apply_transposed_constraint_matrix (const libMesh::DofMap &dof_map, libMesh::PetscVector< double > &rhs)
 
quadrature_key_type getQuadratureKey (const libMesh::QuadratureType quad_type, libMesh::Order order, const bool use_adaptive_quadrature, const double point_density, const bool allow_rules_with_negative_weights, const libMesh::Elem *const elem, const boost::multi_array< double, 2 > &X_node, const double dx_min)
 
void write_elem_partitioning (const std::string &file_name, const libMesh::System &position_system)
 
void write_node_partitioning (const std::string &file_name, const libMesh::System &position_system)
 
std::vector< libMeshWrappers::BoundingBox > get_local_element_bounding_boxes (const libMesh::MeshBase &mesh, const libMesh::System &X_system, const libMesh::QuadratureType quad_type, const libMesh::Order quad_order, const bool use_adaptive_quadrature, const double point_density, bool allow_rules_with_negative_weights, const double patch_dx_min)
 
std::vector< libMeshWrappers::BoundingBox > get_local_element_bounding_boxes (const libMesh::MeshBase &mesh, const libMesh::System &X_system)
 
std::vector< libMeshWrappers::BoundingBox > get_global_element_bounding_boxes (const libMesh::MeshBase &mesh, const std::vector< libMeshWrappers::BoundingBox > &local_bboxes)
 
std::vector< libMeshWrappers::BoundingBox > get_global_element_bounding_boxes (const libMesh::MeshBase &mesh, const libMesh::System &X_system)
 

Detailed Description

Defines "using" declarations for all IBTK namespaces. This header file may be included in application codes, but it MUST NOT be included in any other header (.h) or inline (-inl.h) file in the library.

Typedef Documentation

◆ EigenAlignedVector

template<typename T >
using IBTK::EigenAlignedVector = typedef std::vector<T, Eigen::aligned_allocator<T> >

Eigen types have special alignment requirements and require a specific memory allocator. This is a convenience type alias for a std::vector with the correct allocator.

Function Documentation

◆ abs_equal_eps()

bool IBTK::abs_equal_eps ( double  a,
double  b,
double  eps = std::sqrt(std::numeric_limits<double>::epsilon()) 
)
inline

Check whether the absolute difference between a and b are within the threshold eps.

Note
This function should be used with caution to check numbers that have large magnitudes. In these cases, consider using the rel_equal_eps function.

◆ checked_dereference() [1/2]

template<class T >
const T& IBTK::checked_dereference ( const SAMRAI::tbox::Pointer< T > &  p)

Utility function which asserts that the SAMRAI pointer is not null.

This is useful for writing generic code in which we might want to assert that a pointer is not null in the initialization list where we cannot use the normal assertion macro.

◆ checked_dereference() [2/2]

template<class T >
T& IBTK::checked_dereference ( SAMRAI::tbox::Pointer< T > &  p)

Same idea, but for a non-const pointer.

◆ copyFromHypre() [1/2]

void IBTK::copyFromHypre ( SAMRAI::pdat::CellData< NDIM, double > &  dst_data,
const std::vector< HYPRE_StructVector > &  vectors,
const SAMRAI::hier::Box< NDIM > &  box 
)

Copy data from a vector of Hypre vectors to SAMRAI cell centered data with depth equal to number of Hypre vectors.

Parameters
[out]dst_dataReference to destination for data to be copied.
[in]vectorsVector of Hypre data to copy.
[in]boxBox over which to copy.

◆ copyFromHypre() [2/2]

void IBTK::copyFromHypre ( SAMRAI::pdat::SideData< NDIM, double > &  dst_data,
HYPRE_SStructVector  vector,
const SAMRAI::hier::Box< NDIM > &  box 
)

Copy data from a Hypre vector to SAMRAI side centered data.

Note
This function is specialized for cases when the Hypre vector has one part with number of variables equal to the spatial dimension.
Parameters
[out]dst_dataReference to destination for data to be copied.
[in]vectorVector of Hypre data to copy
[in]boxBox over which to copy.

◆ copyToHypre() [1/2]

void IBTK::copyToHypre ( const std::vector< HYPRE_StructVector > &  vectors,
SAMRAI::pdat::CellData< NDIM, double > &  src_data,
const SAMRAI::hier::Box< NDIM > &  box 
)

Copy data from SAMRAI cell centered data to Hypre vectors.

Parameters
[out]vectorsReference to vector of Hypre vectors to be copied to.
[in]src_dataReference to cell centered data to be copied.
[in]boxBox over which to copy.

◆ copyToHypre() [2/2]

void IBTK::copyToHypre ( HYPRE_SStructVector &  vector,
SAMRAI::pdat::SideData< NDIM, double > &  src_data,
const SAMRAI::hier::Box< NDIM > &  box 
)

Copy data from SAMRAI side centered data to a Hypre vector.

Note
This function is specialized for cases when the Hypre vector has one part with number of variables equal to the spatial dimension.
Parameters
[out]vectorsReference to Hypre vector to be copied to.
[in]src_dataReference to side centered data to be copied.
[in]boxBox over which to copy.

◆ deallocate_vector_data()

void IBTK::deallocate_vector_data ( SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  x,
int  coarsest_ln = invalid_level_number,
int  finest_ln = invalid_level_number 
)
inline

Deallocate a SAMRAIVectorReal.

◆ discontinuous_heaviside()

double IBTK::discontinuous_heaviside ( const double &  phi)
inline

Discontinuous heaviside function.

◆ fill_snapshot_at_time()

void IBTK::fill_snapshot_at_time ( SnapshotCache snapshot_cache,
int  u_idx,
double  t,
int  scr_idx,
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  hierarchy,
const std::string refine_type,
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyDataOpsReal< NDIM, double > >  hier_data_ops = nullptr,
double  period = std::numeric_limits<double>::quiet_NaN() 
)

This function interpolates in time between two snapshots. We find the two closest time points stored in the snapshot cache and linear interpolate in time between them. Each snapshot will be transferred onto the provided hierarchy. A scratch patch index with the same patch layout as u_idx must be provided, and data for both must be already allocated. The refine_type must be a valid refinement operation for the data layout.

This assumes that the time is between two snapshots stored in the cache. An error will occur if two time points can not be found. If there is a single snapshot, this function returns that snapshot.

The hier_data_ops object, if provided, must match the same variable type used in for the snapshots.

If the optional period parameter is provided, this function will treat the first snapshot time point t_1 as also the last snapshot time point t_1 + period.

◆ fill_snapshot_on_hierarchy()

void IBTK::fill_snapshot_on_hierarchy ( SnapshotCache cache,
int  u_idx,
double  time,
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  current_hierarchy,
const std::string snapshot_refine_type,
double  tol = 1.0e-8 
)

Transfer the snapshot at the specified time value within a given tolerance to the patch index u_idx on the supplied patch hierarchy. The patch index u_idx must contain sufficient ghost cell width to perform the operations used by the refinement operator.

Note this can require a significant amount of communication if the supplied patch hierarchy has a different configuration of patches than the snapshot patch hierarchy.

This function does not synchronize the data on the current hierarchy (i.e. no coarsening is performed). It transfers the snapshot using the refinement operator.

◆ free_vector_components()

void IBTK::free_vector_components ( SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  x,
int  coarsest_ln = invalid_level_number,
int  finest_ln = invalid_level_number 
)
inline

Free the components of a SAMRAIVectorReal.

◆ get_min_patch_dx()

double IBTK::get_min_patch_dx ( const SAMRAI::hier::PatchLevel< NDIM > &  patch_level)

Get the smallest cell width on the specified level. This operation is collective.

◆ hypre_array()

std::array< HYPRE_Int, NDIM > IBTK::hypre_array ( const SAMRAI::hier::Index< NDIM > &  index)

Helper function to convert SAMRAI indices to Hypre integers.

Note
Hypre can use 64 bit indices, but SAMRAI IntVectors are always 32.

◆ integrate_elem_rhs_n_vars_n_basis()

template<int n_vars, int n_basis>
void IBTK::integrate_elem_rhs_n_vars_n_basis ( const int  qp_offset,
const std::vector< std::vector< double > > &  phi_F,
const std::vector< double > &  JxW_F,
const std::vector< double > &  F_qp,
std::vector< double > &  F_rhs_concatenated 
)

Assemble the element contribution to a load vector.

The three functions below allow compile-time selection of the number of basis functions and number of variables. Fixing the loop bounds via templates and switch statements makes most IBFE codes a few percent faster.

This function (the last of the three) takes the number of variables and number of basis functions as template argument which allows the compiler to unroll the two inner loops when appropriate. Setting either value to -1 results in run-time selection of loop bounds (which is considerably slower).

Parameters
[in]qp_offsetOffset used to calculate an index into F_qp. The relevant part of the array is assumed to start at n_vars * qp_offset.
[in]phi_FValues of test functions evaluated at quadrature points, indexed by test function number and then quadrature point number.
[in]JxW_FProducts of Jacobian and quadrature weight at each quadrature point.
[in]F_qpGlobally indexed array containing function values at quadrature points. The array is indexed by quadrature point and then by variable: i.e., if n_vars is greater than one then components of the vector-valued function being projected at a certain quadrature point are contiguous.
[out]F_rhs_concatenatedVector containing element integrals of products of test functions and values of the interpolated function. Unlike F_qp, the values in this vector are first indexed by variable and then by quadrature point.

◆ merge_boxes_by_longest_edge()

std::vector< SAMRAI::hier::Box< NDIM > > IBTK::merge_boxes_by_longest_edge ( const std::vector< SAMRAI::hier::Box< NDIM > > &  boxes)

Given a set of boxes (which describe a region in index space), return another set of boxes whose union covers the same index space. The returned set of boxes is formed by merging boxes in boxes along their longest edges.

◆ rel_equal_eps()

bool IBTK::rel_equal_eps ( double  a,
double  b,
double  eps = std::sqrt(std::numeric_limits<double>::epsilon()) 
)
inline

Check whether the relative difference between a and b are within the threshold eps.

Note
This function should be used with caution to check numbers close to zero. In this case, consider using the abs_equal_eps function.

◆ setup_petsc_vector()

Vec IBTK::setup_petsc_vector ( const unsigned int  num_local_nodes,
const unsigned int  depth,
const std::vector< int > &  nonlocal_petsc_indices 
)
inline

Internal function for setting up a PETSc function of given size and depth. See the documentation of LData for more information.

◆ smooth_delta()

double IBTK::smooth_delta ( const double &  phi,
const double &  alpha 
)
inline

Smooth delta function.

◆ smooth_heaviside()

double IBTK::smooth_heaviside ( const double &  phi,
const double &  alpha 
)
inline

Smooth heaviside function.

◆ sum_weighted_elem_solution_n_vars_n_basis()

template<int n_vars, int n_basis, bool weights_are_unity = false>
void IBTK::sum_weighted_elem_solution_n_vars_n_basis ( const int  qp_offset,
const std::vector< std::vector< double > > &  phi_F,
const std::vector< double > &  weights,
const boost::multi_array< double, 2 > &  F_node,
std::vector< double > &  F_w_qp 
)

Compute the product of the finite element representation of the force and (optionally) a set of weights (e.g., JxW values) at all quadrature points on an element. See integrate_elem_rhs for details on the first two template arguments.

Template Parameters
weights_are_unityAssume that all values in weights are 1. This value should be true when we are only interpolating values at quadrature points and otherwise false.
Parameters
[in]qp_offsetOffset used to calculate an index into F_w_qp. The relevant part of the array is assumed to start at n_vars * qp_offset.
[in]phi_FValues of test functions evaluated at quadrature points, indexed by test function number and then quadrature point number.
[in]weightsVector containing weights at quadrature points: in the case of force spreading this is the standard JxW array. If weights_are_unity is true then its values are assumed to be 1.
[in]F_nodeValues of the finite element solution (i.e., the multipliers on each trial function) on the current element.
[out]F_w_qpGlobally indexed array containing the product of the finite element representation of the force and weights at the quadrature points on the current element. The array is indexed by quadrature point and then by variable: i.e., if n_vars is greater than one then variables of the vector-valued force being evaluated at a certain quadrature point are contiguous.

◆ tensor_idx_to_voigt()

int IBTK::tensor_idx_to_voigt ( const std::pair< int, int > &  idx)
inline

Convert a symmetric tensor index to the corresponding Voigt notation index.

◆ update_snapshot()

void IBTK::update_snapshot ( SnapshotCache cache,
int  u_idx,
double  time,
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  current_hierarchy,
double  tol = 1.0e-8 
)

Update a snapshot at the given time. This function copies the data stored in u_idx into the snapshot index in the cache object. The snapshotted hierarchy is also updated to the current hierarchy.

This outputs an error if a snapshot with that time point within the specified tolerance can not be found.

Calls to this function may use a different patch index than the one used in setSnapshot(), but the underlying data layout must be consistent. This data layout is specified by the variable used to construct the cache.

◆ voigt_to_tensor_idx()

std::pair<int, int> IBTK::voigt_to_tensor_idx ( const int  k)
inline

Convert a Voigt notation index to the corresponding symmetric tensor index. This function only returns the upper triangular index.