Loading...
Searching...
No Matches
strong_algebraic_freeslip_enforcement.hpp
Go to the documentation of this file.
1
2#pragma once
3
6
7namespace terra::fe {
8
9/// @brief Helper function to modify the right-hand side vector accordingly for strong free-slip boundary condition
10/// enforcement.
11///
12/// \note The framework documentation features [a detailed description](#boundary-conditions)
13/// of the strong imposition of free-slip boundary conditions.
14///
15/// @param b [in/out] RHS coefficient vector before boundary elimination (but including forcing etc.) - will be modified
16/// in this function to impose the free-slip BCs (after the function returns, this is what is called \f$b_\mathrm{elim}\f$
17/// in the documentation)
18/// @param coords_shell the coordinates of the unit shell, obtained e.g. via
19/// @ref terra::grid::shell::subdomain_unit_sphere_single_shell_coords
20/// @param mask_data the boundary mask data
21/// @param freeslip_boundary_mask the flag that indicates where to apply the conditions
22template < typename ScalarType, typename ScalarTypeGrid, util::FlagLike FlagType >
27 const FlagType& freeslip_boundary_mask )
28{
29 // b <- R b trafo to n-t space
30 linalg::trafo::cartesian_to_normal_tangential_in_place( b.block_1(), coords_shell, mask_data, freeslip_boundary_mask );
31
32 // b <- 0 for normal components at FS boundary
33 kernels::common::assign_masked_else_keep_old<ScalarType, 3, FlagType>( b.block_1().grid_data(), 0, mask_data, freeslip_boundary_mask, 0 );
34
35 // b <- R^T b trafo back to carth space
36 linalg::trafo::normal_tangential_to_cartesian_in_place( b.block_1(), coords_shell, mask_data, freeslip_boundary_mask );
37
38}
39
40} // namespace terra::fe
Block vector consisting of a Q1 vector and a Q1 scalar vector on distributed shell grids.
Definition vector_q1isoq2_q1.hpp:20
const Block1Type & block_1() const
Get const reference to block 1 (vector field).
Definition vector_q1isoq2_q1.hpp:138
const grid::Grid4DDataVec< ScalarType, VecDim > & grid_data() const
Get const reference to grid data.
Definition vector_q1.hpp:280
Definition strong_algebraic_dirichlet_enforcement.hpp:8
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 condit...
Definition strong_algebraic_freeslip_enforcement.hpp:23
Kokkos::View< ScalarType ***[VecDim], Layout > Grid3DDataVec
Definition grid_types.hpp:40
Kokkos::View< ScalarType ****, Layout > Grid4DDataScalar
Definition grid_types.hpp:25
void cartesian_to_normal_tangential_in_place(VectorQ1Vec< ScalarType, 3 > &vec_cartesian, const grid::Grid3DDataVec< ScalarTypeGrid, 3 > &coords_shell, const grid::Grid4DDataScalar< FlagType > mask_data, const FlagType &flag)
Definition local_basis_trafo_normal_tangential.hpp:102
void normal_tangential_to_cartesian_in_place(VectorQ1Vec< ScalarType, 3 > &vec_normal_tangential, const grid::Grid3DDataVec< ScalarTypeGrid, 3 > &coords_shell, const grid::Grid4DDataScalar< FlagType > mask_data, const FlagType &flag)
Definition local_basis_trafo_normal_tangential.hpp:146