DPsim
Loading...
Searching...
No Matches
GpuSparseAdapter.h
1/* Copyright 2017-2021 Institute for Automation of Complex Power Systems,
2 * EONERC, RWTH Aachen University
3 *
4 * This Source Code Form is subject to the terms of the Mozilla Public
5 * License, v. 2.0. If a copy of the MPL was not distributed with this
6 * file, You can obtain one at https://mozilla.org/MPL/2.0/.
7 *********************************************************************************/
8
9#pragma once
10
11#include <bitset>
12#include <iostream>
13#include <list>
14#include <memory>
15#include <unordered_map>
16#include <vector>
17
18#include <dpsim/Config.h>
19#include <dpsim/Definitions.h>
20#include <dpsim/DirectLinearSolver.h>
21
22#include <cusolverSp.h>
23#include <dpsim/Cuda/CudaUtility.h>
24
25namespace DPsim {
27protected:
28 // #### Attributes required for GPU ####
30 cusparseHandle_t mCusparsehandle = nullptr;
31 cusolverSpHandle_t mCusolverhandle = nullptr;
32
34 std::unique_ptr<cuda::CudaMatrix<double, int>> mSysMat = nullptr;
35 std::unique_ptr<Eigen::PermutationMatrix<Eigen::Dynamic>> mTransp = nullptr;
36
43
44 void iluPreconditioner();
45
46 void performFactorization(SparseMatrix &systemMatrix);
47
48private:
50 cuda::Vector<char> pBuffer = 0;
51 cusparseMatDescr_t descr_L = nullptr;
52 cusparseMatDescr_t descr_U = nullptr;
53 csrsv2Info_t info_L = nullptr;
54 csrsv2Info_t info_U = nullptr;
55
56 void checkCusparseStatus(cusparseStatus_t status,
57 std::string additionalInfo = "cuSparse Error:");
58
59public:
62
65
67 virtual ~GpuSparseAdapter();
68
70 virtual void preprocessing(SparseMatrix &systemMatrix,
71 std::vector<std::pair<UInt, UInt>>
72 &listVariableSystemMatrixEntries) override;
73
75 virtual void factorize(SparseMatrix &systemMatrix) override;
76
78 virtual void refactorize(SparseMatrix &systemMatrix) override;
79
81 virtual void partialRefactorize(SparseMatrix &systemMatrix,
82 std::vector<std::pair<UInt, UInt>> &
83 listVariableSystemMatrixEntries) override;
84
86 virtual Matrix solve(Matrix &rightSideVector) override;
87};
88} // namespace DPsim
DirectLinearSolver()=default
Constructor.
virtual ~GpuSparseAdapter()
Destructor.
virtual void refactorize(SparseMatrix &systemMatrix) override
refactorization without partial pivoting
cuda::Vector< double > mGpuLhsVec
LHS-Vector.
cuda::Vector< double > mGpuRhsVec
RHS-Vector.
virtual Matrix solve(Matrix &rightSideVector) override
solution function for a right hand side
DirectLinearSolver()=default
Constructor with logging.
virtual void preprocessing(SparseMatrix &systemMatrix, std::vector< std::pair< UInt, UInt > > &listVariableSystemMatrixEntries) override
preprocessing function pre-ordering and scaling the matrix
std::unique_ptr< cuda::CudaMatrix< double, int > > mSysMat
Systemmatrix on Device.
virtual void factorize(SparseMatrix &systemMatrix) override
factorization function with partial pivoting
cuda::Vector< double > mGpuIntermediateVec
Intermediate Vector.
cusparseHandle_t mCusparsehandle
Solver-Handle.
virtual void partialRefactorize(SparseMatrix &systemMatrix, std::vector< std::pair< UInt, UInt > > &listVariableSystemMatrixEntries) override
partial refactorization withouth partial pivoting