1#include <kernel/voxel_assembly/burgers_velo_material_assembler.hpp>
2#include <kernel/lafem/matrix_gather_scatter_helper.hpp>
3#include <kernel/lafem/vector_gather_scatter_helper.hpp>
4#include <kernel/util/math.hpp>
5#include <kernel/lafem/vector_mirror.hpp>
6#include <kernel/global/vector.hpp>
10 namespace VoxelAssembly
15 template<
typename Space_,
typename DT_,
typename IT_,
typename ViscFunc,
typename ViscDFunc,
16 FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
17 void full_burgers_vm_assembler_matrix1_bcsr_host(DT_* matrix_data,
const DT_* conv_data,
18 const IT_* matrix_row_ptr,
const IT_* matrix_col_idx,
Index matrix_num_rows,
Index matrix_num_cols,
20 const DT_* cub_wg,
int num_cubs, DT_ alpha,
21 const IT_* cell_to_dof, [[maybe_unused]]
Index cell_num,
23 const int* coloring_map,
Index coloring_size,
const IT_* cell_to_dof_sorter,
26 ViscFunc visco_func, ViscDFunc visco_d_func)
29 typedef Space_ SpaceType;
31 typedef IT_ IndexType;
33 const DataType& beta{burgers_params.beta};
34 const DataType& sd_delta{burgers_params.sd_delta};
35 const DataType& sd_v_norm{burgers_params.sd_v_norm};
37 const DataType tol_eps =
Math::sqrt(Math::eps<DataType>());
40 const bool need_streamline = (
Math::abs(sd_delta) > DataType(0)) && (sd_v_norm > tol_eps);
41 const bool need_convection =
Math::abs(beta) > DataType(0);
45 constexpr int dim = SpaceHelp::dim;
48 constexpr int num_loc_dofs = SpaceType::DofMappingType::dof_count;
50 constexpr int num_loc_verts = SpaceType::MeshType::template IndexSet<dim, 0>::Type::num_indices;
58 FEAT_PRAGMA_OMP(parallel
for)
59 for(
Index idx = 0; idx < coloring_size; ++idx)
64 LocalMatrixType loc_mat;
65 LocalVectorType local_conv_dofs(DataType(0));
67 IndexType cell = IndexType(coloring_map[idx]);
70 const IndexType* local_dofs = cell_to_dof + cell*num_loc_dofs;
71 const IndexType* local_dof_sorter = cell_to_dof_sorter + cell*num_loc_dofs;
73 SpaceHelp::set_coefficients(local_coeffs, local_dofs_w, nodes, cell);
75 LAFEM::template VectorGatherScatterHelper<SpaceType, DataType, IndexType>::gather_vector_dense(local_conv_dofs,
76 (
const VecValueType*)conv_data, IndexType(matrix_num_rows), local_dofs,DataType(1));
79 VoxelAssembly::Kernel::burgers_velo_material_mat_assembly_kernel<SpaceHelp, LocalMatrixType, LocalVectorType, dim, num_loc_verts, ViscFunc, ViscDFunc, true>(
80 loc_mat, local_conv_dofs, local_coeffs, cub_pt, cub_wg, num_cubs, burgers_params, material_params, need_convection, tol_eps, visco_func, visco_d_func);
84 VoxelAssembly::Kernel::burgers_velo_material_mat_assembly_kernel<SpaceHelp, LocalMatrixType, LocalVectorType, dim, num_loc_verts, ViscFunc, ViscDFunc, false>(
85 loc_mat, local_conv_dofs, local_coeffs, cub_pt, cub_wg, num_cubs, burgers_params, material_params, need_convection, tol_eps, visco_func, visco_d_func);
89 LAFEM::template MatrixGatherScatterHelper<SpaceType, DataType, IndexType, pol_>::scatter_matrix_csr(loc_mat, (MatValueType*)matrix_data, local_dofs, local_dofs, IndexType(matrix_num_rows), IndexType(matrix_num_cols), matrix_row_ptr, matrix_col_idx, alpha, local_dof_sorter);
94 template<
typename Space_,
typename DT_,
typename IT_,
typename ViscFunc>
95 void full_burgers_vm_assembler_vector_bd_host(DT_* vector_data,
96 const DT_* conv_data,
const DT_* primal_data,
Index vec_size,
98 const DT_* cub_wg,
int num_cubs, DT_ alpha,
99 const IT_* cell_to_dof, [[maybe_unused]]
Index cell_num,
101 const int* coloring_map,
Index coloring_size,
103 const ViscFunc& visco_func)
106 typedef Space_ SpaceType;
107 typedef DT_ DataType;
108 typedef IT_ IndexType;
110 const DataType& beta{burgers_params.beta};
112 const bool need_convection =
Math::abs(beta) > DataType(0);
116 constexpr int dim = SpaceHelp::dim;
119 constexpr int num_loc_dofs = SpaceType::DofMappingType::dof_count;
121 constexpr int num_loc_verts = SpaceType::MeshType::template IndexSet<dim, 0>::Type::num_indices;
129 FEAT_PRAGMA_OMP(parallel
for)
130 for(
Index idx = 0; idx < coloring_size; ++idx)
133 LocalVectorType loc_vec(DataType(0));
134 LocalVectorType local_conv_dofs(DataType(0));
135 LocalVectorType local_prim_dofs(DataType(0));
140 IndexType cell = IndexType(coloring_map[idx]);
142 const IndexType* local_dofs = cell_to_dof + cell*num_loc_dofs;
143 SpaceHelp::set_coefficients(local_coeffs, local_dofs_w, nodes, cell);
144 LAFEM::template VectorGatherScatterHelper<SpaceType, DataType, IndexType>::gather_vector_dense(local_prim_dofs,
145 (
const VecValueType*)primal_data, IndexType(vec_size), local_dofs, DataType(1));
147 LAFEM::template VectorGatherScatterHelper<SpaceType, DataType, IndexType>::gather_vector_dense(local_conv_dofs,
148 (
const VecValueType*)conv_data, IndexType(vec_size), local_dofs, DataType(1));
150 VoxelAssembly::Kernel::burgers_velo_material_defect_assembly_kernel<SpaceHelp>(loc_vec, local_prim_dofs, local_conv_dofs, local_coeffs, cub_pt, cub_wg, num_cubs,
151 burgers_params, need_convection, visco_func);
153 LAFEM::template VectorGatherScatterHelper<SpaceType, DataType, IndexType>::scatter_vector_dense(loc_vec,
154 (VecValueType*)vector_data, IndexType(vec_size), local_dofs, alpha);
162 template<
typename Space_,
typename DT_,
typename IT_>
165 const DT_* conv_data,
168 const std::vector<int*>& coloring_maps,
169 const std::vector<Index>& coloring_map_sizes,
173 switch(material_type)
175 case MaterialType::carreau:
177 for(
Index col = 0; col <
Index(coloring_maps.size()); ++col)
180 FEAT::Intern::MatrixGatherScatterPolicy::useLocalSortHelper>(
181 matrix_data.data, conv_data, matrix_data.row_ptr, matrix_data.col_idx, matrix_data.num_rows, matrix_data.num_cols,
186 (
const int*) coloring_maps[col], coloring_map_sizes[col], dof_mapping.
cell_to_dof_sorter,
193 case MaterialType::carreauYasuda:
195 for(
Index col = 0; col <
Index(coloring_maps.size()); ++col)
198 FEAT::Intern::MatrixGatherScatterPolicy::useLocalSortHelper>(
199 matrix_data.data, conv_data, matrix_data.row_ptr, matrix_data.col_idx, matrix_data.num_rows, matrix_data.num_cols,
204 (
const int*) coloring_maps[col], coloring_map_sizes[col], dof_mapping.
cell_to_dof_sorter,
215 template<
typename Space_,
typename DT_,
typename IT_>
218 const DT_* conv_data,
219 const DT_* primal_data,
222 const std::vector<int*>& coloring_maps,
223 const std::vector<Index>& coloring_map_sizes,
227 switch(material_type)
229 case MaterialType::carreau:
231 for(
Index col = 0; col <
Index(coloring_maps.size()); ++col)
233 VoxelAssembly::Kernel::template full_burgers_vm_assembler_vector_bd_host<Space_, DT_, IT_, Intern::ViscoFunctor<DT_, MaterialType::carreau>>(vector_data, conv_data, primal_data,
238 (
const int*) coloring_maps[col], coloring_map_sizes[col], burgers_params,
244 case MaterialType::carreauYasuda:
246 for(
Index col = 0; col <
Index(coloring_maps.size()); ++col)
248 VoxelAssembly::Kernel::template full_burgers_vm_assembler_vector_bd_host<Space_, DT_, IT_, Intern::ViscoFunctor<DT_, MaterialType::carreauYasuda>>(vector_data, conv_data, primal_data,
253 (
const int*) coloring_maps[col], coloring_map_sizes[col], burgers_params,
280#ifdef FEAT_HAVE_HALFMATH
295#ifdef FEAT_HAVE_HALFMATH
312#ifdef FEAT_HAVE_HALFMATH
327#ifdef FEAT_HAVE_HALFMATH
Standard Lagrange-2 Finite-Element space class template.
Tiny Matrix class template.
Tiny Vector class template.
T_ sqrt(T_ x)
Returns the square-root of a value.
T_ abs(T_ x)
Returns the absolute value.
void assemble_burgers_velo_material_defect_host(const Space_ &space, DT_ *vector_data, const DT_ *conv_data, const DT_ *primal_data, const AssemblyCubatureData< DT_ > &cubature, const AssemblyMappingData< DT_, IT_ > &dof_mapping, const std::vector< int * > &coloring_maps, const std::vector< Index > &coloring_map_sizes, DT_ alpha, const AssemblyBurgersData< DT_ > &burgers_params, const AssemblyMaterialData< DT_ > &material_params, MaterialType material_type)
Host kernel wrapper for the defect burgers assembler.
void assemble_burgers_velo_material_csr_host(const Space_ &space, const CSRMatrixData< DT_, IT_ > &matrix_data, const DT_ *conv_data, const AssemblyCubatureData< DT_ > &cubature, const AssemblyMappingData< DT_, IT_ > &dof_mapping, const std::vector< int * > &coloring_maps, const std::vector< Index > &coloring_map_sizes, DT_ alpha, const AssemblyBurgersData< DT_ > &burgers_params, const AssemblyMaterialData< DT_ > &material_params, MaterialType material_type)
Host kernel wrapper for the full matrix burgers assembler.
Namespace for different voxel based assembly methods.
MaterialType
Enum for different material types.
__half Half
Half data type.
std::uint64_t Index
Index data type.
Data for burgers assembler.
A data field for a cubature rule.
const DT_ * cub_wg
The cubature weights.
int num_cubs
Number of cubtaure points.
const void * cub_pt
The cubature point data array.
A data field for all necessary values that define the dof mapping for assembly.
Index cell_num
The number of cells.
const IT_ * cell_to_dof
The cell to dof, where cell_to_dof[i],..., cell_to_dof[i+cell_dofs-1] are the dofs of one cell.
Index node_size
The number of nodes.
const IT_ * cell_to_dof_sorter
Array of sortingindices of cell_to_dof.
const void * nodes
An array of the nodes fitting to the cell_to_dof mapping.
Data for material burgers assembler.