votca 2025.1-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// Standard includes
25#include <cassert>
26
27// Local VOTCA includes
28#include "eigen.h"
29
30#ifdef USE_CUDA
31#include "cudapipeline.h"
32#endif
33
67
68namespace votca {
69namespace xtp {
70
72 public:
74 static Index UsingGPUs();
75 static Index AvailableGPUs();
76 static void SetNoGPUs(Index number);
77
78 // 3c multiply
79 void setOperators(const std::vector<Eigen::MatrixXd>& tensor,
80 const Eigen::MatrixXd& rightoperator);
81 void MultiplyRight(Eigen::MatrixXd& matrix, Index OpenmpThread);
82
83 // 3c
84 void setOperators(const Eigen::MatrixXd& leftoperator,
85 const Eigen::MatrixXd& rightoperator);
86 void MultiplyLeftRight(Eigen::MatrixXd& matrix, Index OpenmpThread);
87
88 // RPA
89 void createTemporaries(Index rows, Index cols);
90 void PushMatrix(const Eigen::MatrixXd& mat, Index OpenmpThread);
91 void A_TDA(const Eigen::VectorXd& vec, Index OpenmpThread);
92
93 // Hd + Hqp + Hd2
94 void createTemporaries(const Eigen::VectorXd& vec,
95 const Eigen::MatrixXd& input, Index rows1, Index rows2,
96 Index cols);
97 void PrepareMatrix1(Eigen::MatrixXd& mat, Index OpenmpThread);
98 void SetTempZero(Index OpenmpThread);
99 void PrepareMatrix2(const Eigen::Block<const Eigen::MatrixXd>& mat, bool Hd2,
100 Index OpenmpThread);
101 void Addvec(const Eigen::VectorXd& row, Index OpenmpThread);
102 void MultiplyRow(Index row, Index OpenmpThread);
103
104 // Hx
105 void createAdditionalTemporaries(Index rows, Index cols);
106 void PushMatrix1(const Eigen::MatrixXd& mat, Index OpenmpThread);
107 void MultiplyBlocks(const Eigen::Block<const Eigen::MatrixXd>& mat, Index i1,
108 Index i2, Index OpenmpThread);
109
110 Eigen::MatrixXd getReductionVar();
111
112 private:
113 template <class T>
115 public:
116 DefaultReference() = default;
117 DefaultReference(T object) : p(&object) {};
118
119 DefaultReference& operator=(const T& object) {
120 p = &object;
121 return *this;
122 }
123
124 const T& operator()() {
125 assert(p != nullptr && "Dangling reference!");
126 return *p;
127 }
128
129 private:
130 const T* p = nullptr;
131 };
132
136
137 struct CPU_data {
138
139 Eigen::MatrixXd& reduce() { return reduce_mat; }
140 void InitializeReduce(Index rows, Index cols) {
141 reduce_mat = Eigen::MatrixXd::Zero(rows, cols);
142 }
143
144 void InitializeVec(Index size) { temp_vec = Eigen::VectorXd::Zero(size); }
145
147 Eigen::MatrixXd temp_mat;
148 Eigen::VectorXd temp_vec;
149 Eigen::MatrixXd reduce_mat;
150 };
151
152 std::vector<CPU_data> cpus_;
153
156
158
159 Index getParentThreadId(Index OpenmpThreadId) const;
160
161 Index getLocalThreadId(Index ParentThreadId) const;
162
163 Index getNumberThreads() const;
164
165#ifdef USE_CUDA
166 bool isGPUthread(Index ParentThreadId) const;
167
168 struct GPU_data {
169
170 explicit GPU_data(Index i)
171 : Id(i), pipeline(std::make_unique<CudaPipeline>(int(i))) {
172 ;
173 }
174
175 Index Id;
176 std::unique_ptr<CudaPipeline> pipeline;
177 std::vector<std::unique_ptr<CudaMatrix>> temp;
178
179 CudaMatrix& Mat(Index i) { return *temp[i]; }
180 CudaPipeline& pipe() { return *pipeline; }
181 void activateGPU() { checkCuda(cudaSetDevice(pipeline->getDeviceId())); }
182
183 void push_back(const Eigen::MatrixXd& m) {
184 temp.push_back(std::make_unique<CudaMatrix>(m, pipeline->get_stream()));
185 }
186 void push_back(Index rows, Index cols) {
187 temp.push_back(
188 std::make_unique<CudaMatrix>(rows, cols, pipeline->get_stream()));
189 }
190
191 void resize(Index id, Index rows, Index cols) {
192 temp[id] =
193 std::make_unique<CudaMatrix>(rows, cols, pipeline->get_stream());
194 }
195 };
196
197 std::vector<GPU_data> gpus_;
198 static bool isInVector(Index Id, const std::vector<GPU_data>& vec);
199#endif
200};
201
202} // namespace xtp
203} // namespace votca
204
205#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.
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
DefaultReference< Eigen::MatrixXd > ref_mat
void InitializeReduce(Index rows, Index cols)