LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Classes | Functions
projects::ipdg_stokes::assemble Namespace Reference

Classes

class  PiecewiseBoundaryNormalJumpAssembler
 Local assembler for a matrix mapping boundary basis function coefficients to jumps. More...
 
class  PiecewiseConstElementMatrixProvider
 Element matrix provider for the stokes system. More...
 
class  PiecewiseConstElementVectorProvider
 Element vector provider for the stokes system. More...
 

Functions

std::tuple< lf::assemble::COOMatrix< double >, Eigen::VectorXd > buildSystemMatrixNoFlow (const std::shared_ptr< const lf::mesh::Mesh > &mesh, const lf::assemble::DofHandler &dofh, const std::function< Eigen::Vector2d(const Eigen::Vector2d &)> &f, const std::function< Eigen::Vector2d(const lf::mesh::Entity &)> &dirichlet_data, double sigma, const lf::quad::QuadRule &quadrule, bool modified_penalty=false)
 Build the system matrix for the stokes system with no flow boundary conditions. More...
 
std::tuple< lf::assemble::COOMatrix< double >, Eigen::VectorXd, Eigen::VectorXd > buildSystemMatrixInOutFlow (const std::shared_ptr< const lf::mesh::Mesh > &mesh, const lf::assemble::DofHandler &dofh, const std::function< Eigen::Vector2d(const Eigen::Vector2d &)> &f, const std::function< Eigen::Vector2d(const lf::mesh::Entity &)> &dirichlet_data, double sigma, const lf::quad::QuadRule &quadrule, bool modified_penalty=false)
 Build the system matrix for the stokes system with in- and out flow boundary conditions. More...
 
Eigen::VectorXd createOffsetFunction (const std::shared_ptr< const lf::mesh::Mesh > &mesh, const lf::mesh::utils::MeshDataSet< bool > &boundary, const lf::assemble::DofHandler &dofh, const std::function< Eigen::Vector2d(const lf::mesh::Entity &)> &dirichlet_data, const Eigen::SparseMatrix< double > &A)
 Compute the offset function for given boundary conditions. More...
 

Function Documentation

◆ buildSystemMatrixInOutFlow()

std::tuple< lf::assemble::COOMatrix< double >, Eigen::VectorXd, Eigen::VectorXd > projects::ipdg_stokes::assemble::buildSystemMatrixInOutFlow ( const std::shared_ptr< const lf::mesh::Mesh > &  mesh,
const lf::assemble::DofHandler dofh,
const std::function< Eigen::Vector2d(const Eigen::Vector2d &)> &  f,
const std::function< Eigen::Vector2d(const lf::mesh::Entity &)> &  dirichlet_data,
double  sigma,
const lf::quad::QuadRule quadrule,
bool  modified_penalty = false 
)

Build the system matrix for the stokes system with in- and out flow boundary conditions.

Parameters
meshA shared pointer to some mesh instance
dofhThe dof handler for the given mesh
fA functor mapping 2D-coordinates to volumetric forces
dirichlet_dataA functor mapping entities on the mesh boundary to Dirichlet data for the velocity field
sigmaThe stabilization factor
quadruleThe quadrature rule used to integrate the volumetric forces over an element
modified_penaltyIf true, uses a modified penalty term where the jumps are not scaled by 1/|e_n|
Returns
A tuple with the first value containing the assembles COOMatrix of the stokes system, the second value being equal to the right hand side of the LSE and the third value contains the offset function

Definition at line 59 of file build_system_matrix.cc.

References lf::assemble::AssembleMatrixLocally(), lf::assemble::AssembleVectorLocally(), createOffsetFunction(), lf::assemble::DofHandler::Entity(), lf::assemble::FixFlaggedSolutionComponents(), lf::mesh::utils::flagEntitiesOnBoundary(), lf::base::kPoint, lf::assemble::COOMatrix< SCALAR >::makeSparse(), lf::assemble::DofHandler::NumDofs(), and lf::mesh::Entity::RefEl().

◆ buildSystemMatrixNoFlow()

std::tuple< lf::assemble::COOMatrix< double >, Eigen::VectorXd > projects::ipdg_stokes::assemble::buildSystemMatrixNoFlow ( const std::shared_ptr< const lf::mesh::Mesh > &  mesh,
const lf::assemble::DofHandler dofh,
const std::function< Eigen::Vector2d(const Eigen::Vector2d &)> &  f,
const std::function< Eigen::Vector2d(const lf::mesh::Entity &)> &  dirichlet_data,
double  sigma,
const lf::quad::QuadRule quadrule,
bool  modified_penalty = false 
)

Build the system matrix for the stokes system with no flow boundary conditions.

Parameters
meshA shared pointer to some mesh instance
dofhThe dof handler for the given mesh
fA functor mapping 2D-coordinates to volumetric forces
dirichlet_dataA functor mapping entities on the mesh boundary to Dirichlet data for the velocity field
sigmaThe stabilization factor
quadruleThe quadrature rule used to integrate the volumetric forces over an element
modified_penaltyIf true, uses a modified penalty term where the jumps are not scaled by 1/|e_n|
Returns
A tuple with the first value containing the assembles COOMatrix of the stokes system and the second value being equal to the right hand side of the LSE

Definition at line 20 of file build_system_matrix.cc.

References lf::assemble::AssembleMatrixLocally(), lf::assemble::AssembleVectorLocally(), lf::assemble::DofHandler::Entity(), lf::assemble::FixFlaggedSolutionComponents(), lf::mesh::utils::flagEntitiesOnBoundary(), lf::base::kPoint, lf::assemble::DofHandler::NumDofs(), and lf::mesh::Entity::RefEl().

◆ createOffsetFunction()

Eigen::VectorXd projects::ipdg_stokes::assemble::createOffsetFunction ( const std::shared_ptr< const lf::mesh::Mesh > &  mesh,
const lf::mesh::utils::MeshDataSet< bool > &  boundary,
const lf::assemble::DofHandler dofh,
const std::function< Eigen::Vector2d(const lf::mesh::Entity &)> &  dirichlet_data,
const Eigen::SparseMatrix< double > &  A 
)

Compute the offset function for given boundary conditions.

Parameters
meshA shared pointer the the mesh instance used
boundaryA MeshDataSet marking the boundary elements of the mesh
dofhThe DOF handler for the FE space used
dirichlet_dataA function mapping boundary entities to dirichlet values for the flow velocity
AThe full system matrix of the discretized problem
Returns
The offset function for the given boundary conditions

Definition at line 77 of file offset_function.cc.

References lf::assemble::AssembleMatrixLocally(), projects::ipdg_stokes::mesh::computeOutwardNormals(), lf::assemble::DofHandler::GlobalDofIndices(), lf::base::kPoint, lf::base::kSegment, lf::assemble::DofHandler::NumDofs(), and lf::mesh::Entity::RefEl().

Referenced by buildSystemMatrixInOutFlow().