6#include <kernel/voxel_assembly/burgers_assembler.hpp>
7#include <kernel/lafem/matrix_gather_scatter_helper.hpp>
8#include <kernel/lafem/vector_gather_scatter_helper.hpp>
9#include <kernel/util/math.hpp>
10#include <kernel/lafem/vector_mirror.hpp>
11#include <kernel/global/vector.hpp>
15 namespace VoxelAssembly
20 template<
typename Space_,
typename DT_,
typename IT_, FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
21 void full_burgers_assembler_matrix1_bcsr_host(DT_* matrix_data,
const DT_* conv_data,
22 const IT_* matrix_row_ptr,
const IT_* matrix_col_idx,
Index matrix_num_rows,
Index matrix_num_cols,
24 const DT_* cub_wg,
int num_cubs, DT_ alpha,
25 const IT_* cell_to_dof, [[maybe_unused]]
Index cell_num,
27 const int* coloring_map,
Index coloring_size,
const IT_* cell_to_dof_sorter,
31 typedef Space_ SpaceType;
33 typedef IT_ IndexType;
35 const DataType& beta{burgers_params.beta};
36 const DataType& sd_delta{burgers_params.sd_delta};
37 const DataType& sd_v_norm{burgers_params.sd_v_norm};
39 const DataType tol_eps =
Math::sqrt(Math::eps<DataType>());
42 const bool need_streamline = (
Math::abs(sd_delta) > DataType(0)) && (sd_v_norm > tol_eps);
43 const bool need_convection =
Math::abs(beta) > DataType(0);
47 constexpr int dim = SpaceHelp::dim;
50 constexpr int num_loc_dofs = SpaceType::DofMappingType::dof_count;
52 constexpr int num_loc_verts = SpaceType::MeshType::template IndexSet<dim, 0>::Type::num_indices;
60 FEAT_PRAGMA_OMP(parallel
for)
61 for(
Index idx = 0; idx < coloring_size; ++idx)
66 LocalMatrixType loc_mat;
67 LocalVectorType local_conv_dofs(DataType(0));
69 IndexType cell = IndexType(coloring_map[idx]);
72 const IndexType* local_dofs = cell_to_dof + cell*num_loc_dofs;
73 const IndexType* local_dof_sorter = cell_to_dof_sorter + cell*num_loc_dofs;
75 SpaceHelp::set_coefficients(local_coeffs, local_dofs_w, nodes, cell);
77 if(need_convection || need_streamline)
79 LAFEM::template VectorGatherScatterHelper<SpaceType, DataType, IndexType>::gather_vector_dense(local_conv_dofs,
80 (
const VecValueType*)conv_data, IndexType(matrix_num_rows), local_dofs,DataType(1));
83 VoxelAssembly::Kernel::burgers_mat_assembly_kernel<SpaceHelp>(loc_mat, local_conv_dofs, local_coeffs, cub_pt, cub_wg, num_cubs, burgers_params,
84 need_streamline, need_convection, tol_eps);
87 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);
92 template<
typename Space_,
typename DT_,
typename IT_>
93 void full_burgers_assembler_vector_bd_host(DT_* vector_data,
94 const DT_* conv_data,
const DT_* primal_data,
Index vec_size,
96 const DT_* cub_wg,
int num_cubs, DT_ alpha,
97 const IT_* cell_to_dof, [[maybe_unused]]
Index cell_num,
99 const int* coloring_map,
Index coloring_size,
103 typedef Space_ SpaceType;
104 typedef DT_ DataType;
105 typedef IT_ IndexType;
107 const DataType& beta{burgers_params.beta};
109 const bool need_convection =
Math::abs(beta) > DataType(0);
113 constexpr int dim = SpaceHelp::dim;
116 constexpr int num_loc_dofs = SpaceType::DofMappingType::dof_count;
118 constexpr int num_loc_verts = SpaceType::MeshType::template IndexSet<dim, 0>::Type::num_indices;
126 FEAT_PRAGMA_OMP(parallel
for)
127 for(
Index idx = 0; idx < coloring_size; ++idx)
130 LocalVectorType loc_vec(DataType(0));
131 LocalVectorType local_conv_dofs(DataType(0));
132 LocalVectorType local_prim_dofs(DataType(0));
137 IndexType cell = IndexType(coloring_map[idx]);
139 const IndexType* local_dofs = cell_to_dof + cell*num_loc_dofs;
140 SpaceHelp::set_coefficients(local_coeffs, local_dofs_w, nodes, cell);
141 LAFEM::template VectorGatherScatterHelper<SpaceType, DataType, IndexType>::gather_vector_dense(local_prim_dofs,
142 (
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));
151 VoxelAssembly::Kernel::burgers_defect_assembly_kernel<SpaceHelp>(loc_vec, local_prim_dofs, local_conv_dofs, local_coeffs, cub_pt, cub_wg, num_cubs,
152 burgers_params, need_convection);
154 LAFEM::template VectorGatherScatterHelper<SpaceType, DataType, IndexType>::scatter_vector_dense(loc_vec,
155 (VecValueType*)vector_data, IndexType(vec_size), local_dofs, alpha);
173 template<
typename DT_,
int dim_>
178 FEAT_PRAGMA_OMP(parallel
for reduction(max : max_val))
179 for(
int i = 0; i < int(vec_size); ++i)
182 max_val =
Math::max(convect[i].norm_euclid(), max_val);
193 template<
typename Space_,
typename DT_,
typename IT_>
196 const DT_* conv_data,
199 const std::vector<int*>& coloring_maps,
200 const std::vector<Index>& coloring_map_sizes,
203 for(
Index col = 0; col <
Index(coloring_maps.size()); ++col)
205 VoxelAssembly::Kernel::template full_burgers_assembler_matrix1_bcsr_host<Space_, DT_, IT_, FEAT::Intern::MatrixGatherScatterPolicy::useLocalSortHelper>(
206 matrix_data.data, conv_data, matrix_data.row_ptr, matrix_data.col_idx, matrix_data.num_rows, matrix_data.num_cols,
211 (
const int*) coloring_maps[col], coloring_map_sizes[col], dof_mapping.
cell_to_dof_sorter,
218 template<
typename Space_,
typename DT_,
typename IT_>
221 const DT_* conv_data,
222 const DT_* primal_data,
225 const std::vector<int*>& coloring_maps,
226 const std::vector<Index>& coloring_map_sizes,
229 for(
Index col = 0; col <
Index(coloring_maps.size()); ++col)
231 VoxelAssembly::Kernel::template full_burgers_assembler_vector_bd_host<Space_, DT_, IT_>(vector_data, conv_data, primal_data,
236 (
const int*) coloring_maps[col], coloring_map_sizes[col], burgers_params
241 template<
typename DT_,
typename IT_,
int dim_>
247 DT_ glob_res = DT_(0);
249 VoxelAssembly::Kernel::template set_sd_v_norm_host<DT_, dim_>(vec_data, &glob_res, convect.template size<LAFEM::Perspective::native>());
254 template<
typename DT_,
typename IT_,
int dim_>
258 const auto* gate = convect.get_gate();
261 local_norm = gate->max(local_norm);
283#ifdef FEAT_HAVE_HALFMATH
298#ifdef FEAT_HAVE_HALFMATH
309#ifdef FEAT_HAVE_HALFMATH
318#ifdef FEAT_HAVE_HALFMATH
332#ifdef FEAT_HAVE_HALFMATH
347#ifdef FEAT_HAVE_HALFMATH
358#ifdef FEAT_HAVE_HALFMATH
367#ifdef FEAT_HAVE_HALFMATH
Global vector wrapper class template.
Blocked Dense data vector class template.
Handles vector prolongation, restriction and serialization.
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.
T_ max(T_ a, T_ b)
Returns the maximum of two values.
void assemble_burgers_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)
Host kernel wrapper for the defect burgers assembler.
void assemble_burgers_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)
Host kernel wrapper for the full matrix burgers assembler.
DT_ get_sd_v_norm_host(const LAFEM::DenseVectorBlocked< DT_, IT_, dim_ > &convect)
Host kernel wrapper for the local sd_v_norm calculation.
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.
Namespace for different voxel based assembly methods.
__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.