Namespaces | |
| namespace | triangle |
| namespace | wedge |
| Features for wedge elements. | |
Functions | |
| template<typename ScalarType , linalg::OperatorLike OperatorType, typename FlagType > | |
| void | strong_algebraic_dirichlet_enforcement_poisson_like (OperatorType &A_neumann, OperatorType &A_neumann_diag, const linalg::VectorQ1Scalar< ScalarType > &g, linalg::VectorQ1Scalar< ScalarType > &tmp, linalg::VectorQ1Scalar< ScalarType > &b, const grid::Grid4DDataScalar< FlagType > &mask_data, const FlagType &dirichlet_boundary_mask) |
| Helper function to modify the right-hand side vector accordingly for strong Dirichlet boundary condition enforcement. | |
| template<typename ScalarType , linalg::OperatorLike OperatorType, typename FlagType > | |
| void | strong_algebraic_dirichlet_enforcement_vectorlaplace_like (OperatorType &A_neumann, OperatorType &A_neumann_diag, const linalg::VectorQ1Vec< ScalarType > &g, linalg::VectorQ1Vec< ScalarType > &tmp, linalg::VectorQ1Vec< ScalarType > &b, const grid::Grid4DDataScalar< FlagType > &mask_data, const FlagType &dirichlet_boundary_mask) |
| template<typename ScalarType , util::FlagLike FlagType> | |
| void | strong_algebraic_homogeneous_dirichlet_enforcement_poisson_like (linalg::VectorQ1Scalar< ScalarType > &b, const grid::Grid4DDataScalar< FlagType > &mask_data, const FlagType &dirichlet_boundary_mask) |
| Same as strong_algebraic_dirichlet_enforcement_poisson_like() for homogenous boundary conditions ( \( g = 0 \)). | |
| template<typename ScalarType , linalg::OperatorLike OperatorType, util::FlagLike FlagType> | |
| void | strong_algebraic_velocity_dirichlet_enforcement_stokes_like (OperatorType &K_neumann, OperatorType &K_neumann_diag, const linalg::VectorQ1IsoQ2Q1< ScalarType > &g, linalg::VectorQ1IsoQ2Q1< ScalarType > &tmp, linalg::VectorQ1IsoQ2Q1< ScalarType > &b, const grid::Grid4DDataScalar< FlagType > &mask_data, const FlagType &dirichlet_boundary_mask) |
| Same as strong_algebraic_dirichlet_enforcement_poisson_like() for Stokes-like systems (with strong enforcement of velocity boundary conditions). | |
| template<typename ScalarType , util::FlagLike FlagType> | |
| void | strong_algebraic_homogeneous_velocity_dirichlet_enforcement_stokes_like (linalg::VectorQ1IsoQ2Q1< ScalarType > &b, const grid::Grid4DDataScalar< FlagType > &mask_data, const FlagType &dirichlet_boundary_mask) |
| Same as strong_algebraic_homogeneous_dirichlet_enforcement_poisson_like() for Stokes-like systems (with strong enforcement of zero velocity boundary conditions). | |
| template<typename ScalarType , typename ScalarTypeGrid , util::FlagLike FlagType> | |
| void | strong_algebraic_freeslip_enforcement_in_place (linalg::VectorQ1IsoQ2Q1< ScalarType > &b, const grid::Grid3DDataVec< ScalarTypeGrid, 3 > &coords_shell, const grid::Grid4DDataScalar< FlagType > &mask_data, const FlagType &freeslip_boundary_mask) |
| Helper function to modify the right-hand side vector accordingly for strong free-slip boundary condition enforcement. | |
| void terra::fe::strong_algebraic_dirichlet_enforcement_poisson_like | ( | OperatorType & | A_neumann, |
| OperatorType & | A_neumann_diag, | ||
| const linalg::VectorQ1Scalar< ScalarType > & | g, | ||
| linalg::VectorQ1Scalar< ScalarType > & | tmp, | ||
| linalg::VectorQ1Scalar< ScalarType > & | b, | ||
| const grid::Grid4DDataScalar< FlagType > & | mask_data, | ||
| const FlagType & | dirichlet_boundary_mask | ||
| ) |
Helper function to modify the right-hand side vector accordingly for strong Dirichlet boundary condition enforcement.
| A_neumann | the "Neumann" operator (without any modification at the boundary) |
| A_neumann_diag | the diagonal of the "Neumann" operator |
| g | a coefficient vector with interpolated Dirichlet boundary conditions at the respective boundary nodes, zero elsewhere |
| tmp | a temporary vector |
| b | [in/out] RHS coefficient vector before boundary elimination (but including forcing etc.) - will be modified in this function to impose the Dirichlet BCs (after the function returns, this is what is called \(b_\mathrm{elim}\) in the documentation) |
| mask_data | the boundary mask data |
| dirichlet_boundary_mask | the flag that indicates where to apply the conditions |
| void terra::fe::strong_algebraic_dirichlet_enforcement_vectorlaplace_like | ( | OperatorType & | A_neumann, |
| OperatorType & | A_neumann_diag, | ||
| const linalg::VectorQ1Vec< ScalarType > & | g, | ||
| linalg::VectorQ1Vec< ScalarType > & | tmp, | ||
| linalg::VectorQ1Vec< ScalarType > & | b, | ||
| const grid::Grid4DDataScalar< FlagType > & | mask_data, | ||
| const FlagType & | dirichlet_boundary_mask | ||
| ) |
| void terra::fe::strong_algebraic_freeslip_enforcement_in_place | ( | linalg::VectorQ1IsoQ2Q1< ScalarType > & | b, |
| const grid::Grid3DDataVec< ScalarTypeGrid, 3 > & | coords_shell, | ||
| const grid::Grid4DDataScalar< FlagType > & | mask_data, | ||
| const FlagType & | freeslip_boundary_mask | ||
| ) |
Helper function to modify the right-hand side vector accordingly for strong free-slip boundary condition enforcement.
| b | [in/out] RHS coefficient vector before boundary elimination (but including forcing etc.) - will be modified in this function to impose the free-slip BCs (after the function returns, this is what is called \(b_\mathrm{elim}\) in the documentation) |
| coords_shell | the coordinates of the unit shell, obtained e.g. via terra::grid::shell::subdomain_unit_sphere_single_shell_coords |
| mask_data | the boundary mask data |
| freeslip_boundary_mask | the flag that indicates where to apply the conditions |
| void terra::fe::strong_algebraic_homogeneous_dirichlet_enforcement_poisson_like | ( | linalg::VectorQ1Scalar< ScalarType > & | b, |
| const grid::Grid4DDataScalar< FlagType > & | mask_data, | ||
| const FlagType & | dirichlet_boundary_mask | ||
| ) |
Same as strong_algebraic_dirichlet_enforcement_poisson_like() for homogenous boundary conditions ( \( g = 0 \)).
Does not require most of the steps since \( g = g_A = g_D = 0 \). Still requires solving \( A_\mathrm{elim} x = b_\mathrm{elim} \)after this.
| void terra::fe::strong_algebraic_homogeneous_velocity_dirichlet_enforcement_stokes_like | ( | linalg::VectorQ1IsoQ2Q1< ScalarType > & | b, |
| const grid::Grid4DDataScalar< FlagType > & | mask_data, | ||
| const FlagType & | dirichlet_boundary_mask | ||
| ) |
Same as strong_algebraic_homogeneous_dirichlet_enforcement_poisson_like() for Stokes-like systems (with strong enforcement of zero velocity boundary conditions).
| void terra::fe::strong_algebraic_velocity_dirichlet_enforcement_stokes_like | ( | OperatorType & | K_neumann, |
| OperatorType & | K_neumann_diag, | ||
| const linalg::VectorQ1IsoQ2Q1< ScalarType > & | g, | ||
| linalg::VectorQ1IsoQ2Q1< ScalarType > & | tmp, | ||
| linalg::VectorQ1IsoQ2Q1< ScalarType > & | b, | ||
| const grid::Grid4DDataScalar< FlagType > & | mask_data, | ||
| const FlagType & | dirichlet_boundary_mask | ||
| ) |
Same as strong_algebraic_dirichlet_enforcement_poisson_like() for Stokes-like systems (with strong enforcement of velocity boundary conditions).