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 namespace Nektar
41 {
42 namespace LibUtilities
43 {
44 /// A global linear system.
45 class NekLinSysIterGMRES;
46 
48 {
49 public:
50  /// Support creation through MemoryManager.
51  friend class MemoryManager<NekLinSysIterGMRES>;
52 
55  const LibUtilities::CommSharedPtr &vComm, const int nDimen,
56  const NekSysKey &pKey)
57  {
60  vComm, nDimen, pKey);
61  p->InitObject();
62  return p;
63  }
64  static std::string className;
65 
68  const LibUtilities::CommSharedPtr &vComm, const int nDimen,
69  const NekSysKey &pKey = NekSysKey());
71 
73  {
75  }
76 
77 protected:
78  // This is maximum gmres restart iteration
80 
81  // This is maximum bandwidth of Hessenburg matrix
82  // if use truncted Gmres(m)
84 
85  bool m_NekLinSysLeftPrecon = false;
87 
88  bool m_DifferenceFlag0 = false;
89  bool m_DifferenceFlag1 = false;
90 
91  virtual void v_InitObject();
92 
93  virtual int v_SolveSystem(const int nGlobal,
94  const Array<OneD, const NekDouble> &pInput,
95  Array<OneD, NekDouble> &pOutput, const int nDir,
96  const NekDouble tol, const NekDouble factor);
97 
98 private:
99  /// Actual iterative solve-GMRES
100  int DoGMRES(const int pNumRows, const Array<OneD, const NekDouble> &pInput,
101  Array<OneD, NekDouble> &pOutput, const int pNumDir);
102  /// Actual iterative gmres solver for one restart
103  NekDouble DoGmresRestart(const bool restarted, const bool truncted,
104  const int nGlobal,
105  const Array<OneD, const NekDouble> &pInput,
106  Array<OneD, NekDouble> &pOutput, const int nDir);
107 
108  // Arnoldi process
109  void DoArnoldi(const int starttem, const int endtem, const int nGlobal,
110  const int nDir,
111  // V_total(:,1:nd) total search directions
112  Array<OneD, Array<OneD, NekDouble>> &V_local,
113  // V[nd] current search direction
114  Array<OneD, NekDouble> &Vsingle1,
115  // V[nd+1] new search direction
116  Array<OneD, NekDouble> &Vsingle2,
117  // One line of Hessenburg matrix
118  Array<OneD, NekDouble> &hsingle);
119  // QR fatorization through Givens rotation
120  void DoGivensRotation(const int starttem, const int endtem,
121  const int nGlobal, const int nDir,
123  Array<OneD, NekDouble> &hsingle,
125  // Backward calculation to calculate coeficients
126  // of least square problem
127  // To notice, Hessenburg's columnns and rows are reverse
128  void DoBackward(const int number, Array<OneD, Array<OneD, NekDouble>> &A,
131  static std::string lookupIds[];
132  static std::string def;
133 };
134 } // namespace LibUtilities
135 } // namespace Nektar
136 
137 #endif
#define LIB_UTILITIES_EXPORT
void DoArnoldi(const int starttem, const int endtem, const int nGlobal, const int nDir, Array< OneD, Array< OneD, NekDouble >> &V_local, Array< OneD, NekDouble > &Vsingle1, Array< OneD, NekDouble > &Vsingle2, Array< OneD, NekDouble > &hsingle)
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)
NekLinSysIterGMRES(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, const int nDimen, const NekSysKey &pKey=NekSysKey())
virtual int v_SolveSystem(const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir, const NekDouble tol, const NekDouble factor)
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)
static NekLinSysIterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, const int nDimen, const NekSysKey &pKey)
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.
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:46
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:54
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble