votca 2024.2-dev
Loading...
Searching...
No Matches
openmp_cuda.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2020 The VOTCA Development Team
3 * (http://www.votca.org)
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License")
6 *
7 * You may not use this file except in compliance with the License.
8 * You may obtain a copy of the License at
9 *
10 * http://www.apache.org/licenses/LICENSE-2.0
11 *
12 * Unless required by applicable law or agreed to in writing, software
13 * distributed under the License is distributed on an "AS IS" BASIS,
14 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
15 * See the License for the specific language governing permissions and
16 * limitations under the License.
17 *
18 */
19
20#pragma once
21#ifndef VOTCA_XTP_OPENMP_CUDA_H
22#define VOTCA_XTP_OPENMP_CUDA_H
23
24// Local VOTCA includes
25#include "eigen.h"
26
27#ifdef USE_CUDA
28#include "cudapipeline.h"
29#endif
30
65namespace votca {
66namespace xtp {
67
69 public:
71 static Index UsingGPUs();
72 static Index AvailableGPUs();
73 static void SetNoGPUs(Index number);
74
75 // 3c multiply
76 void setOperators(const std::vector<Eigen::MatrixXd>& tensor,
77 const Eigen::MatrixXd& rightoperator);
78 void MultiplyRight(Eigen::MatrixXd& matrix, Index OpenmpThread);
79
80 // 3c
81 void setOperators(const Eigen::MatrixXd& leftoperator,
82 const Eigen::MatrixXd& rightoperator);
83 void MultiplyLeftRight(Eigen::MatrixXd& matrix, Index OpenmpThread);
84
85 // RPA
86 void createTemporaries(Index rows, Index cols);
87 void PushMatrix(const Eigen::MatrixXd& mat, Index OpenmpThread);
88 void A_TDA(const Eigen::VectorXd& vec, Index OpenmpThread);
89
90 // Hd + Hqp + Hd2
91 void createTemporaries(const Eigen::VectorXd& vec,
92 const Eigen::MatrixXd& input, Index rows1, Index rows2,
93 Index cols);
94 void PrepareMatrix1(Eigen::MatrixXd& mat, Index OpenmpThread);
95 void SetTempZero(Index OpenmpThread);
96 void PrepareMatrix2(const Eigen::Block<const Eigen::MatrixXd>& mat, bool Hd2,
97 Index OpenmpThread);
98 void Addvec(const Eigen::VectorXd& row, Index OpenmpThread);
99 void MultiplyRow(Index row, Index OpenmpThread);
100
101 // Hx
102 void createAdditionalTemporaries(Index rows, Index cols);
103 void PushMatrix1(const Eigen::MatrixXd& mat, Index OpenmpThread);
104 void MultiplyBlocks(const Eigen::Block<const Eigen::MatrixXd>& mat, Index i1,
105 Index i2, Index OpenmpThread);
106
107 Eigen::MatrixXd getReductionVar();
108
109 private:
110 template <class T>
112 public:
113 DefaultReference() = default;
114 DefaultReference(T object) : p(&object) {};
115
116 DefaultReference& operator=(const T& object) {
117 p = &object;
118 return *this;
119 }
120
121 const T& operator()() {
122 assert(p != nullptr && "Dangling reference!");
123 return *p;
124 }
125
126 private:
127 const T* p = nullptr;
128 };
129
133
134 struct CPU_data {
135
136 Eigen::MatrixXd& reduce() { return reduce_mat; }
137 void InitializeReduce(Index rows, Index cols) {
138 reduce_mat = Eigen::MatrixXd::Zero(rows, cols);
139 }
140
141 void InitializeVec(Index size) { temp_vec = Eigen::VectorXd::Zero(size); }
142
144 Eigen::MatrixXd temp_mat;
145 Eigen::VectorXd temp_vec;
146 Eigen::MatrixXd reduce_mat;
147 };
148
149 std::vector<CPU_data> cpus_;
150
153
155
156 Index getParentThreadId(Index OpenmpThreadId) const;
157
158 Index getLocalThreadId(Index ParentThreadId) const;
159
160 Index getNumberThreads() const;
161
162#ifdef USE_CUDA
163 bool isGPUthread(Index ParentThreadId) const;
164
165 struct GPU_data {
166
167 explicit GPU_data(Index i)
168 : Id(i), pipeline(std::make_unique<CudaPipeline>(int(i))) {
169 ;
170 }
171
172 Index Id;
173 std::unique_ptr<CudaPipeline> pipeline;
174 std::vector<std::unique_ptr<CudaMatrix>> temp;
175
176 CudaMatrix& Mat(Index i) { return *temp[i]; }
177 CudaPipeline& pipe() { return *pipeline; }
178 void activateGPU() { checkCuda(cudaSetDevice(pipeline->getDeviceId())); }
179
180 void push_back(const Eigen::MatrixXd& m) {
181 temp.push_back(std::make_unique<CudaMatrix>(m, pipeline->get_stream()));
182 }
183 void push_back(Index rows, Index cols) {
184 temp.push_back(
185 std::make_unique<CudaMatrix>(rows, cols, pipeline->get_stream()));
186 }
187
188 void resize(Index id, Index rows, Index cols) {
189 temp[id] =
190 std::make_unique<CudaMatrix>(rows, cols, pipeline->get_stream());
191 }
192 };
193
194 std::vector<GPU_data> gpus_;
195 static bool isInVector(Index Id, const std::vector<GPU_data>& vec);
196#endif
197};
198
199} // namespace xtp
200} // namespace votca
201
202#endif // VOTCA_XTP_OPENMP_CUDA_H
DefaultReference & operator=(const T &object)
void A_TDA(const Eigen::VectorXd &vec, Index OpenmpThread)
void Addvec(const Eigen::VectorXd &row, Index OpenmpThread)
Index getParentThreadId(Index OpenmpThreadId) const
void PushMatrix1(const Eigen::MatrixXd &mat, Index OpenmpThread)
void MultiplyLeftRight(Eigen::MatrixXd &matrix, Index OpenmpThread)
std::vector< CPU_data > cpus_
Eigen::MatrixXd getReductionVar()
void SetTempZero(Index OpenmpThread)
void MultiplyRight(Eigen::MatrixXd &matrix, Index OpenmpThread)
void createTemporaries(Index rows, Index cols)
void setOperators(const std::vector< Eigen::MatrixXd > &tensor, const Eigen::MatrixXd &rightoperator)
void PushMatrix(const Eigen::MatrixXd &mat, Index OpenmpThread)
void MultiplyBlocks(const Eigen::Block< const Eigen::MatrixXd > &mat, Index i1, Index i2, Index OpenmpThread)
static Index AvailableGPUs()
void MultiplyRow(Index row, Index OpenmpThread)
void createAdditionalTemporaries(Index rows, Index cols)
void PrepareMatrix2(const Eigen::Block< const Eigen::MatrixXd > &mat, bool Hd2, Index OpenmpThread)
static void SetNoGPUs(Index number)
Index getLocalThreadId(Index ParentThreadId) const
static Index UsingGPUs()
DefaultReference< const Eigen::MatrixXd > rOP_
DefaultReference< const Eigen::MatrixXd > lOP_
static Index number_of_gpus
Index getNumberThreads() const
DefaultReference< const Eigen::VectorXd > vec_
void PrepareMatrix1(Eigen::MatrixXd &mat, Index OpenmpThread)
STL namespace.
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
DefaultReference< Eigen::MatrixXd > ref_mat
void InitializeReduce(Index rows, Index cols)