Nektar++
NekLinSysIterGMRES.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NekLinSysIterGMRES.h
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: NekLinSysIterGMRES header
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#ifndef NEKTAR_LIB_UTILITIES_LINEAR_ALGEBRA_NEK_LINSYS_ITERAT_GMRES_H
37#define NEKTAR_LIB_UTILITIES_LINEAR_ALGEBRA_NEK_LINSYS_ITERAT_GMRES_H
38
40
42{
43/// A global linear system.
45
47{
48public:
49 /// Support creation through MemoryManager.
51
54 const LibUtilities::CommSharedPtr &vRowComm, const int nDimen,
55 const NekSysKey &pKey)
56 {
59 pSession, vRowComm, nDimen, pKey);
60 p->InitObject();
61 return p;
62 }
63
64 static std::string className;
65
68 const LibUtilities::CommSharedPtr &vRowComm, const int nDimen,
69 const NekSysKey &pKey = NekSysKey());
71
73 {
75 }
76
77protected:
78 // This is maximum gmres restart iteration
80
81 // This is maximum bandwidth of Hessenburg matrix
82 // if use truncted Gmres(m)
84
85 // This is the maximum number of solution vectors that can be stored
86 // For example, in gmres, it is the max number of Krylov space
87 // search directions can be stored
88 // It determines the max storage usage
90
92
96
97 void v_InitObject() override;
98
99 int v_SolveSystem(const int nGlobal,
100 const Array<OneD, const NekDouble> &pInput,
101 Array<OneD, NekDouble> &pOutput, const int nDir) override;
102
103private:
104 /// Actual iterative solve-GMRES
105 int DoGMRES(const int pNumRows, const Array<OneD, const NekDouble> &pInput,
106 Array<OneD, NekDouble> &pOutput, const int pNumDir);
107
108 /// Actual iterative gmres solver for one restart
109 NekDouble DoGmresRestart(const bool restarted, const bool truncted,
110 const int nGlobal,
111 const Array<OneD, const NekDouble> &pInput,
112 Array<OneD, NekDouble> &pOutput, const int nDir);
113
114 // Arnoldi process
115 void DoArnoldi(const int starttem, const int endtem, const int nGlobal,
116 const int nDir, Array<OneD, NekDouble> &w,
117 // V[nd] current search direction
118 Array<OneD, NekDouble> &Vsingle1,
119 // V[nd+1] new search direction
120 Array<OneD, NekDouble> &Vsingle2,
121 // One line of Hessenburg matrix
122 Array<OneD, NekDouble> &hsingle);
123
124 // QR fatorization through Givens rotation
125 void DoGivensRotation(const int starttem, const int endtem,
126 const int nGlobal, const int nDir,
128 Array<OneD, NekDouble> &hsingle,
130
131 // Backward calculation to calculate coeficients
132 // of least square problem
133 // To notice, Hessenburg's columnns and rows are reverse
134 void DoBackward(const int number, Array<OneD, Array<OneD, NekDouble>> &A,
137
138 static std::string lookupIds[];
139 static std::string def;
140
141 // Hessenburg matrix
143 // Hesseburg matrix after rotation
145 // Total search directions
147};
148} // namespace Nektar::LibUtilities
149
150#endif
#define LIB_UTILITIES_EXPORT
int v_SolveSystem(const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir) override
NekLinSysIterGMRES(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey=NekSysKey())
int DoGMRES(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int pNumDir)
Actual iterative solve-GMRES.
void DoBackward(const int number, Array< OneD, Array< OneD, NekDouble > > &A, const Array< OneD, const NekDouble > &b, Array< OneD, NekDouble > &y)
Array< OneD, Array< OneD, NekDouble > > m_V_total
Array< OneD, Array< OneD, NekDouble > > m_Upper
Array< OneD, Array< OneD, NekDouble > > m_hes
void DoGivensRotation(const int starttem, const int endtem, const int nGlobal, const int nDir, Array< OneD, NekDouble > &c, Array< OneD, NekDouble > &s, Array< OneD, NekDouble > &hsingle, Array< OneD, NekDouble > &eta)
NekDouble DoGmresRestart(const bool restarted, const bool truncted, const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir)
Actual iterative gmres solver for one restart.
void DoArnoldi(const int starttem, const int endtem, const int nGlobal, const int nDir, Array< OneD, NekDouble > &w, Array< OneD, NekDouble > &Vsingle1, Array< OneD, NekDouble > &Vsingle2, Array< OneD, NekDouble > &hsingle)
static NekLinSysIterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< NekLinSysIter > NekLinSysIterSharedPtr
Definition: NekLinSysIter.h:47
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::vector< double > w(NPUPPER)
double NekDouble