FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
defo_assembler.cpp
1// FEAT3: Finite Element Analysis Toolbox, Version 3
2// Copyright (C) 2010 by Stefan Turek & the FEAT group
3// FEAT3 is released under the GNU General Public License version 3,
4// see the file 'copyright.txt' in the top level directory for details.
5
6#include <kernel/voxel_assembly/defo_assembler.hpp>
7#include <kernel/lafem/matrix_gather_scatter_helper.hpp>
8
9namespace FEAT
10{
11 namespace VoxelAssembly
12 {
13 namespace Kernel
14 {
15
16 template<typename Space_, typename DT_, typename IT_, FEAT::Intern::MatrixGatherScatterPolicy pol_ = FEAT::Intern::MatrixGatherScatterPolicy::useLocalOps>
17 void defo_assembler_matrix1_csr_host(DT_* matrix_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,
22 const Tiny::Vector<DT_, Space_::world_dim>* nodes, [[maybe_unused]] Index node_size,
23 const int* coloring_map, Index coloring_size,
24 const IT_* cell_to_dof_sorter, DT_ nu)
25 {
26 // define types
27 typedef Space_ SpaceType;
28 typedef DT_ DataType;
29 typedef IT_ IndexType;
30
32
33 constexpr int dim = SpaceHelp::dim;
34
35 // define local sizes
36 constexpr int num_loc_dofs = SpaceType::DofMappingType::dof_count;
37 // get number of nodes per element
38 constexpr int num_loc_verts = SpaceType::MeshType::template IndexSet<dim, 0>::Type::num_indices;
39 //our local datatypes
40 // typedef Tiny::Vector<DataType, dim> VecValueType;
41 typedef Tiny::Matrix<DataType, dim, dim> MatValueType;
42 // typedef Tiny::Vector<VecValueType, num_loc_dofs> LocalVectorType;
44
45
46
47 FEAT_PRAGMA_OMP(parallel for)
48 for(Index idx = 0; idx < coloring_size; ++idx)
49 {
50 // define local coefficients
52 // and local matrix
53 LocalMatrixType loc_mat;
54 // now do work for this cell
55 IndexType cell = IndexType(coloring_map[idx]);
56 // std::cout << "Starting with cell " << cell << "\n";
57 const IndexSetWrapper<IndexType> local_dofs_w{cell_to_dof, IndexType(num_loc_dofs)};
58 const IndexType* local_dofs = cell_to_dof + cell*num_loc_dofs;
59 const IndexType* local_dof_sorter = cell_to_dof_sorter + cell*num_loc_dofs;
60
61 SpaceHelp::set_coefficients(local_coeffs, local_dofs_w, nodes, cell);
62
63 // To minimize code repetition, we call the same kernel from cuda and non cuda build...
64 VoxelAssembly::Kernel::defo_assembly_kernel<SpaceHelp, LocalMatrixType>(loc_mat, local_coeffs, cub_pt, cub_wg, num_cubs, nu);
65
66 // scatter
67 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);
68 }
69 }
70 }
71
72 namespace Arch
73 {
74 template<typename Space_, typename DT_, typename IT_>
75 void assemble_defo_csr_host([[maybe_unused]] const Space_& space,
76 const CSRMatrixData<DT_, IT_>& matrix_data,
77 const AssemblyCubatureData<DT_>& cubature,
78 const AssemblyMappingData<DT_, IT_>& dof_mapping,
79 const std::vector<int*>& coloring_maps,
80 const std::vector<Index>& coloring_map_sizes,
81 DT_ alpha, DT_ nu)
82 {
83 for(Index col = 0; col < Index(coloring_maps.size()); ++col)
84 {
85 VoxelAssembly::Kernel::template defo_assembler_matrix1_csr_host<Space_, DT_, IT_, FEAT::Intern::MatrixGatherScatterPolicy::useLocalSortHelper>(
86 matrix_data.data, matrix_data.row_ptr, matrix_data.col_idx, matrix_data.num_rows, matrix_data.num_cols,
87 (const typename Tiny::Vector<DT_, Space_::world_dim>*) cubature.cub_pt,
88 cubature.cub_wg, cubature.num_cubs, alpha, dof_mapping.cell_to_dof, dof_mapping.cell_num,
89 (const typename Tiny::Vector<DT_, Space_::world_dim>*) dof_mapping.nodes, dof_mapping.node_size,
90 (const int*) coloring_maps[col], coloring_map_sizes[col], dof_mapping.cell_to_dof_sorter, nu
91 );
92 }
93 }
94 }
95 }
96}
97
98//instantiate the templates
99using namespace FEAT;
100using namespace FEAT::VoxelAssembly;
101
102
103/*--------------------Defo Assembler Q2Quad-------------------------------------------------*/
104template void Arch::assemble_defo_csr_host(const Q2StandardQuad&, const CSRMatrixData<double, std::uint32_t>&, const AssemblyCubatureData<double>&, const AssemblyMappingData<double, std::uint32_t>&,
105 const std::vector<int*>&, const std::vector<Index>&, double, double);
106template void Arch::assemble_defo_csr_host(const Q2StandardQuad&, const CSRMatrixData<float, std::uint32_t>&, const AssemblyCubatureData<float>&, const AssemblyMappingData<float, std::uint32_t>&,
107 const std::vector<int*>&, const std::vector<Index>&, float, float);
108template void Arch::assemble_defo_csr_host(const Q2StandardQuad&, const CSRMatrixData<double, std::uint64_t>&, const AssemblyCubatureData<double>&, const AssemblyMappingData<double, std::uint64_t>&,
109 const std::vector<int*>&, const std::vector<Index>&, double, double);
110template void Arch::assemble_defo_csr_host(const Q2StandardQuad&, const CSRMatrixData<float, std::uint64_t>&, const AssemblyCubatureData<float>&, const AssemblyMappingData<float, std::uint64_t>&,
111 const std::vector<int*>&, const std::vector<Index>&, float, float);
112
113
114/*---------------------Defo Assembler Q2Hexa------------------------------------------------------*/
115template void Arch::assemble_defo_csr_host(const Q2StandardHexa&, const CSRMatrixData<double, std::uint32_t>&, const AssemblyCubatureData<double>&, const AssemblyMappingData<double, std::uint32_t>&,
116 const std::vector<int*>&, const std::vector<Index>&, double, double);
117template void Arch::assemble_defo_csr_host(const Q2StandardHexa&, const CSRMatrixData<float, std::uint32_t>&, const AssemblyCubatureData<float>&, const AssemblyMappingData<float, std::uint32_t>&,
118 const std::vector<int*>&, const std::vector<Index>&, float, float);
119template void Arch::assemble_defo_csr_host(const Q2StandardHexa&, const CSRMatrixData<double, std::uint64_t>&, const AssemblyCubatureData<double>&, const AssemblyMappingData<double, std::uint64_t>&,
120 const std::vector<int*>&, const std::vector<Index>&, double, double);
121template void Arch::assemble_defo_csr_host(const Q2StandardHexa&, const CSRMatrixData<float, std::uint64_t>&, const AssemblyCubatureData<float>&, const AssemblyMappingData<float, std::uint64_t>&,
122 const std::vector<int*>&, const std::vector<Index>&, float, float);
Standard Lagrange-2 Finite-Element space class template.
Definition: element.hpp:39
Tiny Matrix class template.
Tiny Vector class template.
Namespace for different voxel based assembly methods.
FEAT namespace.
Definition: adjactor.hpp:12
std::uint64_t Index
Index data type.
const void * cub_pt
The cubature point data array.
A data field for all necessary values that define the dof mapping for assembly.
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.
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.