GpuSender.hpp
Go to the documentation of this file.
1/*
2 Copyright 2022-2023 SINTEF AS, 2025 Equinor ASA
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_GPUISTL_GPUSENDER_HPP
20#define OPM_GPUISTL_GPUSENDER_HPP
21
22#include <dune/istl/owneroverlapcopy.hh>
23
26
27#include <mpi.h>
28
29#include <memory>
30#include <mutex>
31
32namespace Opm::gpuistl
33{
34
42template <class field_type, class OwnerOverlapCopyCommunicationType>
44{
45public:
47
48 GPUSender(const OwnerOverlapCopyCommunicationType& cpuOwnerOverlapCopy)
49 : m_cpuOwnerOverlapCopy(cpuOwnerOverlapCopy)
50 {
51 }
52
53 virtual ~GPUSender() = default;
54
64 virtual void copyOwnerToAll(const X& source, X& dest) const = 0;
65 virtual void initIndexSet() const = 0;
66
74 void project(X& x) const
75 {
76 std::call_once(m_initializedIndices, [&]() { initIndexSet(); });
77 x.setZeroAtIndexSet(*m_indicesCopy);
78 }
79
89 void dot(const X& x, const X& y, field_type& output) const
90 {
91 std::call_once(m_initializedIndices, [&]() { initIndexSet(); });
92
93 const auto dotAtRank = x.dot(y, *m_indicesOwner);
94 output = m_cpuOwnerOverlapCopy.communicator().sum(dotAtRank);
95 }
96
103 field_type norm(const X& x) const
104 {
105 auto xDotX = field_type(0);
106 dot(x, x, xDotX);
107
108 // using std::sqrt;
109 return std::sqrt(xDotX);
110 }
111
112
117 const ::Dune::Communication<MPI_Comm>& communicator() const
118 {
119 return m_cpuOwnerOverlapCopy.communicator();
120 }
121
122protected:
123 // Used to call the initIndexSet. Note that this is kind of a
124 // premature optimization, in the sense that we could just initialize these indices
125 // always, but they are not always used.
126 mutable std::once_flag m_initializedIndices;
127 mutable std::unique_ptr<GpuVector<int>> m_indicesOwner;
128 mutable std::unique_ptr<GpuVector<int>> m_indicesCopy;
129 const OwnerOverlapCopyCommunicationType& m_cpuOwnerOverlapCopy;
130};
131
132} // namespace Opm::gpuistl
133
134#endif // OPM_GPUISTL_GPUSENDER_HPP
GPUSender is a wrapper class for classes which will implement copOwnerToAll This is implemented with ...
Definition: GpuSender.hpp:44
std::unique_ptr< GpuVector< int > > m_indicesCopy
Definition: GpuSender.hpp:128
const OwnerOverlapCopyCommunicationType & m_cpuOwnerOverlapCopy
Definition: GpuSender.hpp:129
void project(X &x) const
project will project x to the owned subspace
Definition: GpuSender.hpp:74
virtual void copyOwnerToAll(const X &source, X &dest) const =0
copyOwnerToAll will copy the data in source to all processes.
std::unique_ptr< GpuVector< int > > m_indicesOwner
Definition: GpuSender.hpp:127
virtual ~GPUSender()=default
virtual void initIndexSet() const =0
void dot(const X &x, const X &y, field_type &output) const
dot will carry out the dot product between x and y on the owned indices, then sum up the result acros...
Definition: GpuSender.hpp:89
field_type norm(const X &x) const
norm computes the l^2-norm of x across processes.
Definition: GpuSender.hpp:103
GPUSender(const OwnerOverlapCopyCommunicationType &cpuOwnerOverlapCopy)
Definition: GpuSender.hpp:48
std::once_flag m_initializedIndices
Definition: GpuSender.hpp:126
const ::Dune::Communication< MPI_Comm > & communicator() const
communicator returns the MPI communicator used by this GPUSender
Definition: GpuSender.hpp:117
Definition: AmgxInterface.hpp:38