LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
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... | |
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.
mesh | A shared pointer to some mesh instance |
dofh | The dof handler for the given mesh |
f | A functor mapping 2D-coordinates to volumetric forces |
dirichlet_data | A functor mapping entities on the mesh boundary to Dirichlet data for the velocity field |
sigma | The stabilization factor |
quadrule | The quadrature rule used to integrate the volumetric forces over an element |
modified_penalty | If true, uses a modified penalty term where the jumps are not scaled by 1/|e_n| |
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().
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.
mesh | A shared pointer to some mesh instance |
dofh | The dof handler for the given mesh |
f | A functor mapping 2D-coordinates to volumetric forces |
dirichlet_data | A functor mapping entities on the mesh boundary to Dirichlet data for the velocity field |
sigma | The stabilization factor |
quadrule | The quadrature rule used to integrate the volumetric forces over an element |
modified_penalty | If true, uses a modified penalty term where the jumps are not scaled by 1/|e_n| |
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().
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.
mesh | A shared pointer the the mesh instance used |
boundary | A MeshDataSet marking the boundary elements of the mesh |
dofh | The DOF handler for the FE space used |
dirichlet_data | A function mapping boundary entities to dirichlet values for the flow velocity |
A | The full system matrix of the discretized problem |
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().