Project
Loading...
Searching...
No Matches
MinResSolve.h
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
15
16#ifndef ALICEO2_FWDALIGN_MINRESSOLVE_H
17#define ALICEO2_FWDALIGN_MINRESSOLVE_H
18
19#include <TObject.h>
20#include <TVectorD.h>
21#include <TString.h>
22
23namespace o2
24{
25namespace fwdalign
26{
27
28class MatrixSq;
29class MatrixSparse;
30class SymBDMatrix;
31
36class MinResSolve : public TObject
37{
38
39 public:
40 enum { kPreconBD = 1,
44 enum { kSolMinRes,
47
48 public:
51
53 MinResSolve(const MatrixSq* mat, const TVectorD* rhs);
54
56 MinResSolve(const MatrixSq* mat, const double* rhs);
57
60
62 ~MinResSolve() override;
63
66
68 Bool_t SolveMinRes(Double_t* VecSol, Int_t precon = 0, int itnlim = 2000, double rtol = 1e-12);
69
71 Bool_t SolveMinRes(TVectorD& VecSol, Int_t precon = 0, int itnlim = 2000, double rtol = 1e-12);
72
74 Bool_t SolveFGMRES(Double_t* VecSol, Int_t precon = 0, int itnlim = 2000, double rtol = 1e-12, int nkrylov = 60);
75
77 Bool_t SolveFGMRES(TVectorD& VecSol, Int_t precon = 0, int itnlim = 2000, double rtol = 1e-12, int nkrylov = 60);
78
80 Bool_t InitAuxMinRes();
81
83 Bool_t InitAuxFGMRES(int nkrylov);
84
86 void ApplyPrecon(const TVectorD& vecRHS, TVectorD& vecOut) const;
87
89 void ApplyPrecon(const double* vecRHS, double* vecOut) const;
90
92 Int_t BuildPrecon(Int_t val = 0);
93
94 Int_t GetPrecon() const { return fPrecon; }
95
97 void ClearAux();
98
100 Int_t BuildPreconBD(Int_t hwidth);
101
103 Int_t BuildPreconILUK(Int_t lofM);
104
106 Int_t BuildPreconILUKDense(Int_t lofM);
107
109 Int_t PreconILUKsymb(Int_t lofM);
110
112 Int_t PreconILUKsymbDense(Int_t lofM);
113
114 protected:
115 Int_t fSize;
116 Int_t fPrecon;
118 Double_t* fRHS;
119
120 Double_t* fPVecY;
121 Double_t* fPVecR1; // aux. space
122 Double_t* fPVecR2; // aux. space
123 Double_t* fPVecV; // aux. space
124 Double_t* fPVecW; // aux. space
125 Double_t* fPVecW1; // aux. space
126 Double_t* fPVecW2; // aux. space
127 Double_t** fPvv; // aux. space
128 Double_t** fPvz; // aux. space
129 Double_t** fPhh; // aux. space
130 Double_t* fDiagLU; // aux space
131 MatrixSparse* fMatL; // aux. space
132 MatrixSparse* fMatU; // aux. space
133 SymBDMatrix* fMatBD; // aux. space
134
136};
137
138} // namespace fwdalign
139} // namespace o2
140
141#endif
for solving large system of linear equations
Definition MinResSolve.h:37
MinResSolve()
default constructor
Int_t BuildPreconBD(Int_t hwidth)
build Band-Diagonal preconditioner
MatrixSq * fMatrix
matrix defining the equations
void ApplyPrecon(const TVectorD &vecRHS, TVectorD &vecOut) const
apply precond.
Int_t BuildPreconILUKDense(Int_t lofM)
ILUK preconditioner.
MinResSolve & operator=(const MinResSolve &rhs)
assignment op.
Bool_t InitAuxFGMRES(int nkrylov)
init auxiliary space for fgmres
Int_t fPrecon
preconditioner type
Int_t PreconILUKsymb(Int_t lofM)
ILUK preconditioner.
Double_t * fPVecY
aux. space
Int_t BuildPreconILUK(Int_t lofM)
ILUK preconditioner.
ClassDefOverride(MinResSolve, 0)
void ClearAux()
clear aux. space
Int_t BuildPrecon(Int_t val=0)
preconditioner building
~MinResSolve() override
destructor
Double_t * fRHS
right hand side
Int_t fSize
dimension of the input matrix
Int_t PreconILUKsymbDense(Int_t lofM)
ILUK preconditioner.
Bool_t InitAuxMinRes()
init auxiliary space for MinRes
Bool_t SolveMinRes(Double_t *VecSol, Int_t precon=0, int itnlim=2000, double rtol=1e-12)
MINRES method (for symmetric matrices)
Bool_t SolveFGMRES(Double_t *VecSol, Int_t precon=0, int itnlim=2000, double rtol=1e-12, int nkrylov=60)
FGMRES method (for general symmetric matrices)
Symmetric Matrix Class.
Definition SymBDMatrix.h:33
GLenum src
Definition glcorearb.h:1767
GLuint GLfloat * val
Definition glcorearb.h:1582
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...