vector_operations.hpp
Go to the documentation of this file.
1/*
2 Copyright 2022-2023 SINTEF AS
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef OPM_CUISTL_VECTOR_OPERATIONS_HPP
20#define OPM_CUISTL_VECTOR_OPERATIONS_HPP
21#include <cstddef>
22#include <cublas_v2.h>
23namespace Opm::cuistl::detail
24{
25
32template <class T>
33void setVectorValue(T* deviceData, size_t numberOfElements, const T& value);
34
41template <class T>
42void setZeroAtIndexSet(T* deviceData, size_t numberOfElements, const int* indices);
43
57template <class T>
58T innerProductAtIndices(cublasHandle_t cublasHandle, const T* deviceA, const T* deviceB, T* buffer, size_t numberOfElements, const int* indices);
59
60template <class T>
61void prepareSendBuf(const T* deviceA, T* buffer, size_t numberOfElements, const int* indices);
62template <class T>
63void syncFromRecvBuf(T* deviceA, T* buffer, size_t numberOfElements, const int* indices);
64
77template <class T>
78void weightedDiagMV(const T* squareBlockVector,
79 const size_t numberOfRows,
80 const size_t blocksize,
81 T relaxationFactor,
82 const T* srcVec,
83 T* dstVec);
84} // namespace Opm::cuistl::detail
85#endif
Definition: cublas_safe_call.hpp:32
void prepareSendBuf(const T *deviceA, T *buffer, size_t numberOfElements, const int *indices)
void syncFromRecvBuf(T *deviceA, T *buffer, size_t numberOfElements, const int *indices)
void setVectorValue(T *deviceData, size_t numberOfElements, const T &value)
setVectorValue sets every element of deviceData to value
T innerProductAtIndices(cublasHandle_t cublasHandle, const T *deviceA, const T *deviceB, T *buffer, size_t numberOfElements, const int *indices)
innerProductAtIndices computes the inner product between deviceA[indices] and deviceB[indices]
void weightedDiagMV(const T *squareBlockVector, const size_t numberOfRows, const size_t blocksize, T relaxationFactor, const T *srcVec, T *dstVec)
Compue the weighted matrix vector product where the matrix is diagonal, the diagonal is a vector,...
void setZeroAtIndexSet(T *deviceData, size_t numberOfElements, const int *indices)
setZeroAtIndexSet sets deviceData to zero in the indices of contained in indices