FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
FEAT::VoxelAssembly::Kernel Namespace Reference

Namespace for the actual assembly kernels. More...

Functions

template<typename SpaceHelp_ , typename LocVecType_ , int dim_, int num_verts_>
CUDA_HOST_DEVICE void burgers_defect_assembly_kernel (LocVecType_ &loc_vec, const LocVecType_ &local_prim_dofs, const LocVecType_ &local_conv_dofs, const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &local_coeffs, const typename SpaceHelp_::DomainPointType *cub_pt, const typename SpaceHelp_::DataType *cub_wg, const int num_cubs, const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &burgers_params, const bool need_convection)
 Burgers Vector assembly kernel. More...
 
template<typename SpaceHelp_ , typename LocMatType_ , typename LocVecType_ , int dim_, int num_verts_>
CUDA_HOST_DEVICE void burgers_mat_assembly_kernel (LocMatType_ &loc_mat, const LocVecType_ &local_conv_dofs, const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &local_coeffs, const typename SpaceHelp_::DomainPointType *cub_pt, const typename SpaceHelp_::DataType *cub_wg, const int num_cubs, const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &burgers_params, const bool need_streamline, const bool need_convection, const typename SpaceHelp_::DataType tol_eps)
 Burgers Matrix assembly kernel. More...
 
template<typename SpaceHelp_ , typename LocVecType_ , int dim_, int num_verts_, typename ViscFunc_ >
CUDA_HOST_DEVICE void burgers_velo_material_defect_assembly_kernel (LocVecType_ &loc_vec, const LocVecType_ &local_prim_dofs, const LocVecType_ &local_conv_dofs, const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &local_coeffs, const typename SpaceHelp_::DomainPointType *cub_pt, const typename SpaceHelp_::DataType *cub_wg, const int num_cubs, const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &burgers_params, const bool need_convection, ViscFunc_ visco_func)
 Burgers Vector assembly kernel. More...
 
template<typename SpaceHelp_ , typename LocMatType_ , typename LocVecType_ , int dim_, int num_verts_, typename ViscFunc_ , typename ViscDerFunc_ , bool need_stream_diff_>
CUDA_HOST_DEVICE void burgers_velo_material_mat_assembly_kernel (LocMatType_ &loc_mat, const LocVecType_ &local_conv_dofs, const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &local_coeffs, const typename SpaceHelp_::DomainPointType *cub_pt, const typename SpaceHelp_::DataType *cub_wg, const int num_cubs, const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &burgers_params, const VoxelAssembly::AssemblyMaterialData< typename SpaceHelp_::DataType > &material_params, const bool need_convection, const typename SpaceHelp_::DataType tol_eps, ViscFunc_ visco_func, ViscDerFunc_ visco_d_func)
 Burgers Velo Material Matrix assembly kernel. More...
 
template<typename Space_ , typename DT_ , typename IT_ , FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
void defo_assembler_matrix1_csr_host (DT_ *matrix_data, const IT_ *matrix_row_ptr, const IT_ *matrix_col_idx, Index matrix_num_rows, Index matrix_num_cols, const Tiny::Vector< DT_, Space_::world_dim > *cub_pt, const DT_ *cub_wg, int num_cubs, DT_ alpha, const IT_ *cell_to_dof, Index cell_num, const Tiny::Vector< DT_, Space_::world_dim > *nodes, Index node_size, const int *coloring_map, Index coloring_size, const IT_ *cell_to_dof_sorter, DT_ nu)
 
template<typename SpaceHelp_ , typename LocMatType_ , int dim_, int num_verts_>
CUDA_HOST_DEVICE void defo_assembly_kernel (LocMatType_ &loc_mat, const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &local_coeffs, const typename SpaceHelp_::DomainPointType *cub_pt, const typename SpaceHelp_::DataType *cub_wg, const int num_cubs, const typename SpaceHelp_::DataType nu)
 
template<typename Space_ , typename DT_ , typename IT_ , FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
void full_burgers_assembler_matrix1_bcsr_host (DT_ *matrix_data, const DT_ *conv_data, const IT_ *matrix_row_ptr, const IT_ *matrix_col_idx, Index matrix_num_rows, Index matrix_num_cols, const Tiny::Vector< DT_, Space_::world_dim > *cub_pt, const DT_ *cub_wg, int num_cubs, DT_ alpha, const IT_ *cell_to_dof, Index cell_num, const Tiny::Vector< DT_, Space_::world_dim > *nodes, Index node_size, const int *coloring_map, Index coloring_size, const IT_ *cell_to_dof_sorter, const VoxelAssembly::AssemblyBurgersData< DT_ > &burgers_params)
 
template<typename Space_ , typename DT_ , typename IT_ >
void full_burgers_assembler_vector_bd_host (DT_ *vector_data, const DT_ *conv_data, const DT_ *primal_data, Index vec_size, const Tiny::Vector< DT_, Space_::world_dim > *cub_pt, const DT_ *cub_wg, int num_cubs, DT_ alpha, const IT_ *cell_to_dof, Index cell_num, const Tiny::Vector< DT_, Space_::world_dim > *nodes, Index node_size, const int *coloring_map, Index coloring_size, const VoxelAssembly::AssemblyBurgersData< DT_ > &burgers_params)
 
template<typename Space_ , typename DT_ , typename IT_ , typename ViscFunc , typename ViscDFunc , FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
void full_burgers_vm_assembler_matrix1_bcsr_host (DT_ *matrix_data, const DT_ *conv_data, const IT_ *matrix_row_ptr, const IT_ *matrix_col_idx, Index matrix_num_rows, Index matrix_num_cols, const Tiny::Vector< DT_, Space_::world_dim > *cub_pt, const DT_ *cub_wg, int num_cubs, DT_ alpha, const IT_ *cell_to_dof, Index cell_num, const Tiny::Vector< DT_, Space_::world_dim > *nodes, Index node_size, const int *coloring_map, Index coloring_size, const IT_ *cell_to_dof_sorter, const VoxelAssembly::AssemblyBurgersData< DT_ > &burgers_params, const VoxelAssembly::AssemblyMaterialData< DT_ > &material_params, ViscFunc visco_func, ViscDFunc visco_d_func)
 
template<typename Space_ , typename DT_ , typename IT_ , typename ViscFunc >
void full_burgers_vm_assembler_vector_bd_host (DT_ *vector_data, const DT_ *conv_data, const DT_ *primal_data, Index vec_size, const Tiny::Vector< DT_, Space_::world_dim > *cub_pt, const DT_ *cub_wg, int num_cubs, DT_ alpha, const IT_ *cell_to_dof, Index cell_num, const Tiny::Vector< DT_, Space_::world_dim > *nodes, Index node_size, const int *coloring_map, Index coloring_size, const VoxelAssembly::AssemblyBurgersData< DT_ > &burgers_params, const ViscFunc &visco_func)
 
template<typename ThreadGroup_ , typename SpaceHelp_ >
CUDA_DEVICE __forceinline__ void grouped_burgers_mat_alt_assembly_kernel (const ThreadGroup_ &tg, BurgersSharedDataKernelWrapper< SpaceHelp_ > *shared_wrapper, int loc_assemble_size, int assemble_offset, typename SpaceHelp_::DataType *loc_mat, const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &burgers_params, const bool need_streamline, const bool need_convection, const typename SpaceHelp_::DataType tol_eps)
 appears to not work better... compiler optimizes this sufficiently? More...
 
template<typename ThreadGroup_ , typename SpaceHelp_ >
CUDA_DEVICE __forceinline__ void grouped_burgers_mat_alt_prepare_assembly_kernel (const ThreadGroup_ &tg, BurgersSharedDataKernelWrapper< SpaceHelp_ > *shared_wrapper, int loc_assemble_size, int assemble_offset, typename SpaceHelp_::DataType *loc_mat, const typename SpaceHelp_::DataType *local_conv_dofs, const Tiny::Matrix< typename SpaceHelp_::DataType, SpaceHelp_::dim, SpaceHelp_::num_verts > &local_coeffs, const typename SpaceHelp_::DomainPointType *cub_pt, const typename SpaceHelp_::DataType *cub_wg, const int cub_ind, const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &burgers_params, const bool need_streamline, const bool need_convection, const typename SpaceHelp_::DataType tol_eps)
 
template<typename ThreadGroup_ , typename SpaceHelp_ >
CUDA_DEVICE __forceinline__ void grouped_burgers_mat_assembly_kernel (const ThreadGroup_ &tg, BurgersSharedDataKernelWrapper< SpaceHelp_ > *shared_wrapper, int loc_assemble_size, int assemble_offset, typename SpaceHelp_::DataType *loc_mat, const typename SpaceHelp_::DataType *local_conv_dofs, const Tiny::Matrix< typename SpaceHelp_::DataType, SpaceHelp_::dim, SpaceHelp_::num_verts > &local_coeffs, const typename SpaceHelp_::DomainPointType *cub_pt, const typename SpaceHelp_::DataType *cub_wg, const int num_cubs, const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &burgers_params, const bool need_streamline, const bool need_convection, const typename SpaceHelp_::DataType tol_eps)
 Burgers Matrix assembly kernel. More...
 
template<typename Space_ , typename DT_ , typename IT_ , FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
void poisson_assembler_matrix1_csr_host (DT_ *matrix_data, const IT_ *matrix_row_ptr, const IT_ *matrix_col_idx, Index matrix_num_rows, Index matrix_num_cols, const Tiny::Vector< DT_, Space_::world_dim > *cub_pt, const DT_ *cub_wg, int num_cubs, DT_ alpha, const IT_ *cell_to_dof, Index cell_num, const Tiny::Vector< DT_, Space_::world_dim > *nodes, Index node_size, const int *coloring_map, Index coloring_size, const IT_ *cell_to_dof_sorter)
 
template<typename SpaceHelp_ , typename LocMatType_ , int dim_, int num_verts_>
CUDA_HOST_DEVICE void poisson_assembly_kernel (LocMatType_ &loc_mat, const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &local_coeffs, const typename SpaceHelp_::DomainPointType *cub_pt, const typename SpaceHelp_::DataType *cub_wg, const int num_cubs)
 
template<typename DT_ , int dim_>
void set_sd_v_norm_host (const Tiny::Vector< DT_, dim_ > *convect, DT_ *result, Index vec_size)
 Reduces the max local vector norm of a convection vector. More...
 

Detailed Description

Namespace for the actual assembly kernels.

Function Documentation

◆ burgers_defect_assembly_kernel()

template<typename SpaceHelp_ , typename LocVecType_ , int dim_, int num_verts_>
CUDA_HOST_DEVICE void FEAT::VoxelAssembly::Kernel::burgers_defect_assembly_kernel ( LocVecType_ &  loc_vec,
const LocVecType_ &  local_prim_dofs,
const LocVecType_ &  local_conv_dofs,
const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &  local_coeffs,
const typename SpaceHelp_::DomainPointType *  cub_pt,
const typename SpaceHelp_::DataType *  cub_wg,
const int  num_cubs,
const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &  burgers_params,
const bool  need_convection 
)

Burgers Vector assembly kernel.

Outsources the the most inner execution of the assembly of the full burgers term for a cell previously gathered. Local always refers to the local dofs of the current cell.

Template Parameters
SpaceHelp_The psacehelper to be used. Holds most information about the underlying space.
LocMatType_The Matrixtype of the local matrix dofs.
LocVecType_The Vectortype of the local convection dofs.
dim_Dimension of our space. For ease of access.
num_verts_The number of vertices of one cell.
Parameters
[out]loc_vecA reference to the local defect to be assembled.
[in]local_prim_dofsA reference to the local primal values the matrix is implicitly applied to.
[in]local_conv_dofsA reference to the local convection dofs. If convection is not needed values are ignored.
[in]local_coeffsThe local trafo coefficients.
[in]cub_ptA pointer to the cubature point coordinate vector.
[in]cub_wgA pointer to the cubature weights. Same order as cub_pt array.
[in]burgers_paramsA struct holding the burgers parameter configuration.
[in]need_convectionDo we need convection?

Definition at line 901 of file burgers_assembler.hpp.

References FEAT::Math::abs(), FEAT::dom_point, FEAT::Tiny::dot(), and FEAT::value.

◆ burgers_mat_assembly_kernel()

template<typename SpaceHelp_ , typename LocMatType_ , typename LocVecType_ , int dim_, int num_verts_>
CUDA_HOST_DEVICE void FEAT::VoxelAssembly::Kernel::burgers_mat_assembly_kernel ( LocMatType_ &  loc_mat,
const LocVecType_ &  local_conv_dofs,
const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &  local_coeffs,
const typename SpaceHelp_::DomainPointType *  cub_pt,
const typename SpaceHelp_::DataType *  cub_wg,
const int  num_cubs,
const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &  burgers_params,
const bool  need_streamline,
const bool  need_convection,
const typename SpaceHelp_::DataType  tol_eps 
)

Burgers Matrix assembly kernel.

Outsources the the most inner execution of the assembly of the full burgers term for a cell previously gathered. Local always refers to the local dofs of the current cell.

Template Parameters
SpaceHelp_The psacehelper to be used. Holds most information about the underlying space.
LocMatType_The Matrixtype of the local matrix dofs.
LocVecType_The Vectortype of the local convection dofs.
dim_Dimension of our space. For ease of access.
num_verts_The number of vertices of one cell.
Parameters
[out]loc_matA reference to the local matrix to be assembled.
[in]local_conv_dofsA reference to the local convection dofs. If convection is not needed values are ignored.
[in]local_coeffsThe local trafo coefficients.
[in]cub_ptA pointer to the cubature point coordinate vector.
[in]cub_wgA pointer to the cubature weights. Same order as cub_pt array.
[in]burgers_paramsA struct holding the burgers parameter configuration.
[in]need_streamlineDo we need streamline?
[in]need_convectionDo we need convection?
[in]tol_epsTolerance for local stream norm to be regarded as zero.

Definition at line 105 of file burgers_assembler.hpp.

References FEAT::Math::abs(), FEAT::dom_point, FEAT::Tiny::dot(), and FEAT::value.

◆ burgers_velo_material_defect_assembly_kernel()

template<typename SpaceHelp_ , typename LocVecType_ , int dim_, int num_verts_, typename ViscFunc_ >
CUDA_HOST_DEVICE void FEAT::VoxelAssembly::Kernel::burgers_velo_material_defect_assembly_kernel ( LocVecType_ &  loc_vec,
const LocVecType_ &  local_prim_dofs,
const LocVecType_ &  local_conv_dofs,
const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &  local_coeffs,
const typename SpaceHelp_::DomainPointType *  cub_pt,
const typename SpaceHelp_::DataType *  cub_wg,
const int  num_cubs,
const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &  burgers_params,
const bool  need_convection,
ViscFunc_  visco_func 
)

Burgers Vector assembly kernel.

Outsources the the most inner execution of the assembly of the full burgers term for a cell previously gathered. Local always refers to the local dofs of the current cell.

Template Parameters
SpaceHelp_The psacehelper to be used. Holds most information about the underlying space.
LocMatType_The Matrixtype of the local matrix dofs.
LocVecType_The Vectortype of the local convection dofs.
dim_Dimension of our space. For ease of access.
num_verts_The number of vertices of one cell.
ViscFunc_Function mapping from R to R.
Parameters
[out]loc_vecA reference to the local defect to be assembled.
[in]local_prim_dofsA reference to the local primal values the matrix is implicitly applied to.
[in]local_conv_dofsA reference to the local convection dofs. If convection is not needed values are ignored.
[in]local_coeffsThe local trafo coefficients.
[in]cub_ptA pointer to the cubature point coordinate vector.
[in]cub_wgA pointer to the cubature weights. Same order as cub_pt array.
[in]burgers_paramsA struct holding the burgers parameter configuration.
[in]need_convectionDo we need convection?
[in]visco_funcFunction mapping invariants of velocity jacobian to viscosity.

Definition at line 704 of file burgers_velo_material_assembler.hpp.

References FEAT::Math::abs(), FEAT::dom_point, FEAT::Tiny::dot(), FEAT::Math::sqrt(), and FEAT::value.

◆ burgers_velo_material_mat_assembly_kernel()

template<typename SpaceHelp_ , typename LocMatType_ , typename LocVecType_ , int dim_, int num_verts_, typename ViscFunc_ , typename ViscDerFunc_ , bool need_stream_diff_>
CUDA_HOST_DEVICE void FEAT::VoxelAssembly::Kernel::burgers_velo_material_mat_assembly_kernel ( LocMatType_ &  loc_mat,
const LocVecType_ &  local_conv_dofs,
const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &  local_coeffs,
const typename SpaceHelp_::DomainPointType *  cub_pt,
const typename SpaceHelp_::DataType *  cub_wg,
const int  num_cubs,
const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &  burgers_params,
const VoxelAssembly::AssemblyMaterialData< typename SpaceHelp_::DataType > &  material_params,
const bool  need_convection,
const typename SpaceHelp_::DataType  tol_eps,
ViscFunc_  visco_func,
ViscDerFunc_  visco_d_func 
)

Burgers Velo Material Matrix assembly kernel.

Outsources the the most inner execution of the assembly of the full burgers term for a cell previously gathered. Local always refers to the local dofs of the current cell.

Template Parameters
SpaceHelp_The psacehelper to be used. Holds most information about the underlying space.
LocMatType_The Matrixtype of the local matrix dofs.
LocVecType_The Vectortype of the local convection dofs.
dim_Dimension of our space. For ease of access.
num_verts_The number of vertices of one cell.
ViscFunc_Function mapping from R to R.
ViscDerFunc_Function mapping from R to R.
need_stream_diff_Streamdiff required?
Parameters
[out]loc_matA reference to the local matrix to be assembled.
[in]local_conv_dofsA reference to the local convection dofs. If convection is not needed values are ignored.
[in]local_coeffsThe local trafo coefficients.
[in]cub_ptA pointer to the cubature point coordinate vector.
[in]cub_wgA pointer to the cubature weights. Same order as cub_pt array.
[in]burgers_paramsA struct holding the burgers parameter configuration.
[in]material_paramsA struct holding the additional material parameters.
[in]need_convectionDo we need convection?
[in]tol_epsTolerance for local stream norm to be regarded as zero.
[in]visco_funcFunction mapping invariants of velocity jacobian to viscosity.
[in]visco_d_funcFunction representing the derivative of visc_func.

Definition at line 131 of file burgers_velo_material_assembler.hpp.

References FEAT::Math::abs(), FEAT::dom_point, FEAT::Tiny::dot(), FEAT::Tiny::Vector< T_, n_, s_ >::set_mat_vec_mult(), FEAT::Math::sqrt(), and FEAT::value.

◆ defo_assembler_matrix1_csr_host()

template<typename Space_ , typename DT_ , typename IT_ , FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
void FEAT::VoxelAssembly::Kernel::defo_assembler_matrix1_csr_host ( DT_ *  matrix_data,
const IT_ *  matrix_row_ptr,
const IT_ *  matrix_col_idx,
Index  matrix_num_rows,
Index  matrix_num_cols,
const Tiny::Vector< DT_, Space_::world_dim > *  cub_pt,
const DT_ *  cub_wg,
int  num_cubs,
DT_  alpha,
const IT_ *  cell_to_dof,
Index  cell_num,
const Tiny::Vector< DT_, Space_::world_dim > *  nodes,
Index  node_size,
const int *  coloring_map,
Index  coloring_size,
const IT_ *  cell_to_dof_sorter,
DT_  nu 
)

Definition at line 17 of file defo_assembler.cpp.

◆ defo_assembly_kernel()

template<typename SpaceHelp_ , typename LocMatType_ , int dim_, int num_verts_>
CUDA_HOST_DEVICE void FEAT::VoxelAssembly::Kernel::defo_assembly_kernel ( LocMatType_ &  loc_mat,
const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &  local_coeffs,
const typename SpaceHelp_::DomainPointType *  cub_pt,
const typename SpaceHelp_::DataType *  cub_wg,
const int  num_cubs,
const typename SpaceHelp_::DataType  nu 
)

Definition at line 42 of file defo_assembler.hpp.

◆ full_burgers_assembler_matrix1_bcsr_host()

template<typename Space_ , typename DT_ , typename IT_ , FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
void FEAT::VoxelAssembly::Kernel::full_burgers_assembler_matrix1_bcsr_host ( DT_ *  matrix_data,
const DT_ *  conv_data,
const IT_ *  matrix_row_ptr,
const IT_ *  matrix_col_idx,
Index  matrix_num_rows,
Index  matrix_num_cols,
const Tiny::Vector< DT_, Space_::world_dim > *  cub_pt,
const DT_ *  cub_wg,
int  num_cubs,
DT_  alpha,
const IT_ *  cell_to_dof,
Index  cell_num,
const Tiny::Vector< DT_, Space_::world_dim > *  nodes,
Index  node_size,
const int *  coloring_map,
Index  coloring_size,
const IT_ *  cell_to_dof_sorter,
const VoxelAssembly::AssemblyBurgersData< DT_ > &  burgers_params 
)

Definition at line 21 of file burgers_assembler.cpp.

◆ full_burgers_assembler_vector_bd_host()

template<typename Space_ , typename DT_ , typename IT_ >
void FEAT::VoxelAssembly::Kernel::full_burgers_assembler_vector_bd_host ( DT_ *  vector_data,
const DT_ *  conv_data,
const DT_ *  primal_data,
Index  vec_size,
const Tiny::Vector< DT_, Space_::world_dim > *  cub_pt,
const DT_ *  cub_wg,
int  num_cubs,
DT_  alpha,
const IT_ *  cell_to_dof,
Index  cell_num,
const Tiny::Vector< DT_, Space_::world_dim > *  nodes,
Index  node_size,
const int *  coloring_map,
Index  coloring_size,
const VoxelAssembly::AssemblyBurgersData< DT_ > &  burgers_params 
)

Definition at line 93 of file burgers_assembler.cpp.

◆ full_burgers_vm_assembler_matrix1_bcsr_host()

template<typename Space_ , typename DT_ , typename IT_ , typename ViscFunc , typename ViscDFunc , FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
void FEAT::VoxelAssembly::Kernel::full_burgers_vm_assembler_matrix1_bcsr_host ( DT_ *  matrix_data,
const DT_ *  conv_data,
const IT_ *  matrix_row_ptr,
const IT_ *  matrix_col_idx,
Index  matrix_num_rows,
Index  matrix_num_cols,
const Tiny::Vector< DT_, Space_::world_dim > *  cub_pt,
const DT_ *  cub_wg,
int  num_cubs,
DT_  alpha,
const IT_ *  cell_to_dof,
Index  cell_num,
const Tiny::Vector< DT_, Space_::world_dim > *  nodes,
Index  node_size,
const int *  coloring_map,
Index  coloring_size,
const IT_ *  cell_to_dof_sorter,
const VoxelAssembly::AssemblyBurgersData< DT_ > &  burgers_params,
const VoxelAssembly::AssemblyMaterialData< DT_ > &  material_params,
ViscFunc  visco_func,
ViscDFunc  visco_d_func 
)

Definition at line 17 of file burgers_velo_material_assembler.cpp.

◆ full_burgers_vm_assembler_vector_bd_host()

template<typename Space_ , typename DT_ , typename IT_ , typename ViscFunc >
void FEAT::VoxelAssembly::Kernel::full_burgers_vm_assembler_vector_bd_host ( DT_ *  vector_data,
const DT_ *  conv_data,
const DT_ *  primal_data,
Index  vec_size,
const Tiny::Vector< DT_, Space_::world_dim > *  cub_pt,
const DT_ *  cub_wg,
int  num_cubs,
DT_  alpha,
const IT_ *  cell_to_dof,
Index  cell_num,
const Tiny::Vector< DT_, Space_::world_dim > *  nodes,
Index  node_size,
const int *  coloring_map,
Index  coloring_size,
const VoxelAssembly::AssemblyBurgersData< DT_ > &  burgers_params,
const ViscFunc &  visco_func 
)

Definition at line 95 of file burgers_velo_material_assembler.cpp.

◆ grouped_burgers_mat_alt_assembly_kernel()

template<typename ThreadGroup_ , typename SpaceHelp_ >
CUDA_DEVICE __forceinline__ void FEAT::VoxelAssembly::Kernel::grouped_burgers_mat_alt_assembly_kernel ( const ThreadGroup_ &  tg,
BurgersSharedDataKernelWrapper< SpaceHelp_ > *  shared_wrapper,
int  loc_assemble_size,
int  assemble_offset,
typename SpaceHelp_::DataType *  loc_mat,
const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &  burgers_params,
const bool  need_streamline,
const bool  need_convection,
const typename SpaceHelp_::DataType  tol_eps 
)

appears to not work better... compiler optimizes this sufficiently?

Definition at line 500 of file burgers_assembler.hpp.

References FEAT::Math::abs(), FEAT::dom_point, FEAT::Tiny::dot(), and FEAT::value.

◆ grouped_burgers_mat_alt_prepare_assembly_kernel()

template<typename ThreadGroup_ , typename SpaceHelp_ >
CUDA_DEVICE __forceinline__ void FEAT::VoxelAssembly::Kernel::grouped_burgers_mat_alt_prepare_assembly_kernel ( const ThreadGroup_ &  tg,
BurgersSharedDataKernelWrapper< SpaceHelp_ > *  shared_wrapper,
int  loc_assemble_size,
int  assemble_offset,
typename SpaceHelp_::DataType *  loc_mat,
const typename SpaceHelp_::DataType *  local_conv_dofs,
const Tiny::Matrix< typename SpaceHelp_::DataType, SpaceHelp_::dim, SpaceHelp_::num_verts > &  local_coeffs,
const typename SpaceHelp_::DomainPointType *  cub_pt,
const typename SpaceHelp_::DataType *  cub_wg,
const int  cub_ind,
const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &  burgers_params,
const bool  need_streamline,
const bool  need_convection,
const typename SpaceHelp_::DataType  tol_eps 
)

Definition at line 341 of file burgers_assembler.hpp.

◆ grouped_burgers_mat_assembly_kernel()

template<typename ThreadGroup_ , typename SpaceHelp_ >
CUDA_DEVICE __forceinline__ void FEAT::VoxelAssembly::Kernel::grouped_burgers_mat_assembly_kernel ( const ThreadGroup_ &  tg,
BurgersSharedDataKernelWrapper< SpaceHelp_ > *  shared_wrapper,
int  loc_assemble_size,
int  assemble_offset,
typename SpaceHelp_::DataType *  loc_mat,
const typename SpaceHelp_::DataType *  local_conv_dofs,
const Tiny::Matrix< typename SpaceHelp_::DataType, SpaceHelp_::dim, SpaceHelp_::num_verts > &  local_coeffs,
const typename SpaceHelp_::DomainPointType *  cub_pt,
const typename SpaceHelp_::DataType *  cub_wg,
const int  num_cubs,
const VoxelAssembly::AssemblyBurgersData< typename SpaceHelp_::DataType > &  burgers_params,
const bool  need_streamline,
const bool  need_convection,
const typename SpaceHelp_::DataType  tol_eps 
)

Burgers Matrix assembly kernel.

Outsources the the most inner execution of the assembly of the full burgers term for a cell previously gathered. Local always refers to the local dofs of the current cell.

Template Parameters
SpaceHelp_The psacehelper to be used. Holds most information about the underlying space.
LocMatType_The Matrixtype of the local matrix dofs.
LocVecType_The Vectortype of the local convection dofs.
dim_Dimension of our space. For ease of access.
num_verts_The number of vertices of one cell.
Parameters
[out]loc_matA reference to the local matrix to be assembled.
[in]local_conv_dofsA reference to the local convection dofs. If convection is not needed values are ignored.
[in]local_coeffsThe local trafo coefficients.
[in]cub_ptA pointer to the cubature point coordinate vector.
[in]cub_wgA pointer to the cubature weights. Same order as cub_pt array.
[in]burgers_paramsA struct holding the burgers parameter configuration.
[in]need_streamlineDo we need streamline?
[in]need_convectionDo we need convection?
[in]tol_epsTolerance for local stream norm to be regarded as zero.

Definition at line 652 of file burgers_assembler.hpp.

References FEAT::Math::abs(), FEAT::dom_point, FEAT::Tiny::dot(), and FEAT::value.

◆ poisson_assembler_matrix1_csr_host()

template<typename Space_ , typename DT_ , typename IT_ , FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
void FEAT::VoxelAssembly::Kernel::poisson_assembler_matrix1_csr_host ( DT_ *  matrix_data,
const IT_ *  matrix_row_ptr,
const IT_ *  matrix_col_idx,
Index  matrix_num_rows,
Index  matrix_num_cols,
const Tiny::Vector< DT_, Space_::world_dim > *  cub_pt,
const DT_ *  cub_wg,
int  num_cubs,
DT_  alpha,
const IT_ *  cell_to_dof,
Index  cell_num,
const Tiny::Vector< DT_, Space_::world_dim > *  nodes,
Index  node_size,
const int *  coloring_map,
Index  coloring_size,
const IT_ *  cell_to_dof_sorter 
)

Definition at line 17 of file poisson_assembler.cpp.

◆ poisson_assembly_kernel()

template<typename SpaceHelp_ , typename LocMatType_ , int dim_, int num_verts_>
CUDA_HOST_DEVICE void FEAT::VoxelAssembly::Kernel::poisson_assembly_kernel ( LocMatType_ &  loc_mat,
const Tiny::Matrix< typename SpaceHelp_::DataType, dim_, num_verts_ > &  local_coeffs,
const typename SpaceHelp_::DomainPointType *  cub_pt,
const typename SpaceHelp_::DataType *  cub_wg,
const int  num_cubs 
)

Definition at line 47 of file poisson_assembler.hpp.

◆ set_sd_v_norm_host()

template<typename DT_ , int dim_>
void FEAT::VoxelAssembly::Kernel::set_sd_v_norm_host ( const Tiny::Vector< DT_, dim_ > *  convect,
DT_ *  result,
Index  vec_size 
)

Reduces the max local vector norm of a convection vector.

This uses multiple reductions with red black sumations on a single block and an atomic add at the end to reduce the vector value.

Attention
Only one dimensional kernel and blockDim.x * num_kernels >= vec_size required.
Warning
Will fail if if blockdim is not a multiple of 2 (which you should do in any case)
Parameters
[in]convectDevice memory to the convection vector in native view.
[out]resultPtr to device global variable.
[in]vec_sizeThe size of the convection vector in native view.

Definition at line 174 of file burgers_assembler.cpp.

References FEAT::Math::max().