preCICE
Loading...
Searching...
No Matches
KokkosPUMKernels.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Kokkos_Core.hpp>
7
8template <typename MemorySpace>
10 MatrixOffsetView<MemorySpace> dst, int nCluster);
11
12// returns true, if successful, currently not tuned as it is not performance critical
13template <typename MemorySpace>
14bool compute_weights(const int nCluster,
15 const int avgOutClusterSize,
16 const offset_1d_type nWeights,
17 const int nMeshVertices,
18 const int dim,
23 const CompactPolynomialC2 &w,
24 VectorView<MemorySpace> normalizedWeights);
25
26// the input assembly
27template <typename EvalFunctionType, typename MemorySpace>
29 int nCluster, // Number of local systems
30 int dim, // Dimension of points
31 int avgClusterSize,
32 int maxInClusterSize,
33 EvalFunctionType f,
34 const VectorOffsetView<MemorySpace> &inOffsets, // vertex offsets (length N+1)
35 const GlobalIDView<MemorySpace> &globalInIDs,
36 const MeshView<MemorySpace> &inCoords, // meshes
37 const MatrixOffsetView<MemorySpace> &matrixOffsets,
39
40template <typename EvalFunctionType, typename MemorySpace>
41void do_batched_assembly(int nCluster, // Number of local systems
42 int dim, // Dimension of points
43 int avgClusterSize,
44 EvalFunctionType f,
45 const VectorOffsetView<MemorySpace> &inOffsets, // vertex offsets (length N+1)
46 const GlobalIDView<MemorySpace> &globalInIDs,
47 const MeshView<MemorySpace> &inCoords, // meshes
48 const VectorOffsetView<MemorySpace> &targetOffsets,
49 const GlobalIDView<MemorySpace> &globalTargetIDs,
50 const MeshView<MemorySpace> &targetCoords,
51 const MatrixOffsetView<MemorySpace> &matrixOffsets,
53
54template <typename MemorySpace>
55void do_batched_qr(int nCluster,
56 int dim,
57 int avgClusterSize,
58 int maxClusterSize,
60 GlobalIDView<MemorySpace> globalInIDs,
65
66// currently unused
67template <typename MemorySpace>
68void do_qr_solve(int nCluster,
69 int dim,
70 int maxInClusterSize,
72 GlobalIDView<MemorySpace> globalInIDs,
78 const VectorView<MemorySpace> weights,
80 GlobalIDView<MemorySpace> globalOutIDs,
82 MeshView<MemorySpace> outMesh);
83
84template <typename MemorySpace>
85void do_batched_lu(
86 int nCluster,
87 int avgClusterSize,
88 const MatrixOffsetView<MemorySpace> &matrixOffsets,
90
91template <bool polynomial, bool evaluation_op_available, typename EvalFunctionType, typename MemorySpace>
93 int nCluster,
94 int dim,
95 int avgInClusterSize,
96 int maxInClusterSize,
97 int maxOutClusterSize,
98 EvalFunctionType f,
99 const VectorOffsetView<MemorySpace> &rhsOffsets,
100 const GlobalIDView<MemorySpace> &globalRhsIDs,
102 const MatrixOffsetView<MemorySpace> &matrixOffsets,
103 const VectorView<MemorySpace> &matrices,
104 const VectorView<MemorySpace> &normalizedWeights,
105 const MatrixOffsetView<MemorySpace> &evalOffsets,
106 const VectorView<MemorySpace> &evalMat,
107 const VectorOffsetView<MemorySpace> &outOffsets,
108 const GlobalIDView<MemorySpace> &globalOutIDs,
110 const MeshView<MemorySpace> &inMesh,
111 const MeshView<MemorySpace> &outMesh,
112 const VectorView<MemorySpace> &qrMatrix,
113 const VectorView<MemorySpace> &qrTau,
114 const PivotView<MemorySpace> &qrP);
115
116} // namespace precice::mapping::kernel
117
void do_batched_qr(int nCluster, int dim, int avgClusterSize, int maxClusterSize, VectorOffsetView< MemorySpace > inOffsets, GlobalIDView< MemorySpace > globalInIDs, MeshView< MemorySpace > inMesh, VectorView< MemorySpace > qrMatrix, VectorView< MemorySpace > qrTau, PivotView< MemorySpace > qrP)
void do_batched_assembly(int nCluster, int dim, int avgClusterSize, EvalFunctionType f, const VectorOffsetView< MemorySpace > &inOffsets, const GlobalIDView< MemorySpace > &globalInIDs, const MeshView< MemorySpace > &inCoords, const VectorOffsetView< MemorySpace > &targetOffsets, const GlobalIDView< MemorySpace > &globalTargetIDs, const MeshView< MemorySpace > &targetCoords, const MatrixOffsetView< MemorySpace > &matrixOffsets, VectorView< MemorySpace > matrices)
bool compute_weights(const int nCluster, const int avgOutClusterSize, const offset_1d_type nWeights, const int nMeshVertices, const int dim, VectorOffsetView< MemorySpace > offsets, MeshView< MemorySpace > centers, GlobalIDView< MemorySpace > globalIDs, MeshView< MemorySpace > mesh, const CompactPolynomialC2 &w, VectorView< MemorySpace > normalizedWeights)
void do_input_assembly(int nCluster, int dim, int avgClusterSize, int maxInClusterSize, EvalFunctionType f, const VectorOffsetView< MemorySpace > &inOffsets, const GlobalIDView< MemorySpace > &globalInIDs, const MeshView< MemorySpace > &inCoords, const MatrixOffsetView< MemorySpace > &matrixOffsets, VectorView< MemorySpace > matrices)
void do_qr_solve(int nCluster, int dim, int maxInClusterSize, VectorOffsetView< MemorySpace > inOffsets, GlobalIDView< MemorySpace > globalInIDs, VectorView< MemorySpace > inData, MeshView< MemorySpace > inMesh, VectorView< MemorySpace > qrMatrix, VectorView< MemorySpace > qrTau, PivotView< MemorySpace > qrP, const VectorView< MemorySpace > weights, VectorOffsetView< MemorySpace > outOffsets, GlobalIDView< MemorySpace > globalOutIDs, VectorView< MemorySpace > outData, MeshView< MemorySpace > outMesh)
void do_batched_solve(int nCluster, int dim, int avgInClusterSize, int maxInClusterSize, int maxOutClusterSize, EvalFunctionType f, const VectorOffsetView< MemorySpace > &rhsOffsets, const GlobalIDView< MemorySpace > &globalRhsIDs, VectorView< MemorySpace > rhs, const MatrixOffsetView< MemorySpace > &matrixOffsets, const VectorView< MemorySpace > &matrices, const VectorView< MemorySpace > &normalizedWeights, const MatrixOffsetView< MemorySpace > &evalOffsets, const VectorView< MemorySpace > &evalMat, const VectorOffsetView< MemorySpace > &outOffsets, const GlobalIDView< MemorySpace > &globalOutIDs, VectorView< MemorySpace > out, const MeshView< MemorySpace > &inMesh, const MeshView< MemorySpace > &outMesh, const VectorView< MemorySpace > &qrMatrix, const VectorView< MemorySpace > &qrTau, const PivotView< MemorySpace > &qrP)
void compute_offsets(const VectorOffsetView< MemorySpace > src1, const VectorOffsetView< MemorySpace > src2, MatrixOffsetView< MemorySpace > dst, int nCluster)
void do_batched_lu(int nCluster, int avgClusterSize, const MatrixOffsetView< MemorySpace > &matrixOffsets, VectorView< MemorySpace > matrices)
Kokkos::View< int *, MemorySpace > PivotView
Kokkos::View< offset_2d_type *, MemorySpace > MatrixOffsetView
Kokkos::View< offset_1d_type *, MemorySpace > VectorOffsetView
Kokkos::View< VertexID *, MemorySpace > GlobalIDView
ExecutionSpace::size_type offset_1d_type
Kokkos::View< double *, MemorySpace > VectorView
Kokkos::View< double **, Kokkos::LayoutRight, MemorySpace > MeshView