GpuObliviousMPISender.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_GPUOBLIVIOUSMPISENDER_HPP
20#define OPM_GPUISTL_GPUOBLIVIOUSMPISENDER_HPP
21
22#include <dune/istl/owneroverlapcopy.hh>
23
27
28#include <mpi.h>
29
30#include <memory>
31
32namespace Opm::gpuistl
33{
34
42template <class field_type, int block_size, class OwnerOverlapCopyCommunicationType>
43class GPUObliviousMPISender : public GPUSender<field_type, OwnerOverlapCopyCommunicationType>
44{
45public:
47
48 explicit GPUObliviousMPISender(const OwnerOverlapCopyCommunicationType& cpuOwnerOverlapCopy)
49 : GPUSender<field_type, OwnerOverlapCopyCommunicationType>(cpuOwnerOverlapCopy)
50 {
51 }
52
53 void copyOwnerToAll(const X& source, X& dest) const override
54 {
55 // TODO: [perf] Can we reduce copying from the GPU here?
56 // TODO: [perf] Maybe create a global buffer instead?
57 auto sourceAsDuneVector = source.template asDuneBlockVector<block_size>();
58 auto destAsDuneVector = dest.template asDuneBlockVector<block_size>();
59 this->m_cpuOwnerOverlapCopy.copyOwnerToAll(sourceAsDuneVector, destAsDuneVector);
60 dest.copyFromHost(destAsDuneVector);
61 }
62
63private:
64 void initIndexSet() const override
65 {
66 // We need indices that we we will use in the project, dot and norm calls.
67 // TODO: [premature perf] Can this be run once per instance? Or do we need to rebuild every time?
68 const auto& pis = this->m_cpuOwnerOverlapCopy.indexSet();
69 std::vector<int> indicesCopyOnCPU;
70 std::vector<int> indicesOwnerCPU;
71 for (const auto& index : pis) {
72 if (index.local().attribute() == Dune::OwnerOverlapCopyAttributeSet::copy) {
73 for (int component = 0; component < block_size; ++component) {
74 indicesCopyOnCPU.push_back(index.local().local() * block_size + component);
75 }
76 }
77
78 if (index.local().attribute() == Dune::OwnerOverlapCopyAttributeSet::owner) {
79 for (int component = 0; component < block_size; ++component) {
80 indicesOwnerCPU.push_back(index.local().local() * block_size + component);
81 }
82 }
83 }
84
85 this->m_indicesCopy = std::make_unique<GpuVector<int>>(indicesCopyOnCPU);
86 this->m_indicesOwner = std::make_unique<GpuVector<int>>(indicesOwnerCPU);
87 }
88};
89
90} // namespace Opm::gpuistl
91
92#endif // OPM_GPUISTL_GPUOBLIVIOUSMPISENDER_HPP
Derived class of GPUSender that handles MPI calls that should NOT use GPU direct communicatoin The im...
Definition: GpuObliviousMPISender.hpp:44
void copyOwnerToAll(const X &source, X &dest) const override
copyOwnerToAll will copy the data in source to all processes.
Definition: GpuObliviousMPISender.hpp:53
GPUObliviousMPISender(const OwnerOverlapCopyCommunicationType &cpuOwnerOverlapCopy)
Definition: GpuObliviousMPISender.hpp:48
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
std::unique_ptr< GpuVector< int > > m_indicesOwner
Definition: GpuSender.hpp:127
Definition: AmgxInterface.hpp:38