Nektar++
PreconCfsBRJ.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File PreconCfsBRJ.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 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: PreconCfsBRJ header
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_PRECONCFSBRJ
36 #define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_PRECONCFSBRJ
37 
39 
40 namespace Nektar
41 {
42 /**
43  * Block Relaxed(weighted) Jacobi iterative (BRJ) Preconditioner for CFS
44  *
45  */
46 class PreconCfsBRJ : public PreconCfsOp
47 {
48 public:
49  friend class MemoryManager<PreconCfsBRJ>;
50 
51  /// Creates an instance of this class
55  const LibUtilities::CommSharedPtr &vComm)
56  {
58  pFields, pSession, vComm);
59  return p;
60  }
61 
62  /// Name of the class
63  static std::string className;
64 
67  const LibUtilities::CommSharedPtr &vComm);
69 
70  virtual bool UpdatePreconMatCheck(const Array<OneD, const NekDouble> &res,
71  const NekDouble dtLambda);
72 
73 protected:
76 
85 
87 
88  virtual void v_InitObject();
89 
90 private:
91  virtual void v_DoPreconCfs(
93  const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
94  const bool &flag);
95 
96  virtual void v_BuildPreconCfs(
98  const Array<OneD, const Array<OneD, NekDouble>> &intmp,
99  const NekDouble time, const NekDouble lambda);
100 
101  template <typename DataType, typename TypeNekBlkMatSharedPtr>
102  void PreconBlkDiag(
104  const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray,
105  const TypeNekBlkMatSharedPtr &PreconMatVars,
106  const DataType &tmpDataType);
107 
108  template <typename DataType>
109  void MinusOffDiag2Rhs(
111  const int nvariables, const int nCoeffs,
112  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
113  Array<OneD, Array<OneD, NekDouble>> &outarray, bool flagUpdateDervFlux,
114  Array<OneD, Array<OneD, NekDouble>> &FwdFluxDeriv,
115  Array<OneD, Array<OneD, NekDouble>> &BwdFluxDeriv,
117  TensorOfArray3D<NekDouble> &wspTrace,
118  Array<OneD, Array<OneD, DataType>> &wspTraceDataType,
119  const TensorOfArray4D<DataType> &TraceJacArray,
120  const TensorOfArray4D<DataType> &TraceJacDerivArray,
121  const Array<OneD, const Array<OneD, DataType>> &TraceJacDerivSign,
122  const TensorOfArray5D<DataType> &TraceIPSymJacArray);
123 
124  template <typename TypeNekBlkMatSharedPtr>
128  const int &nscale = 1);
129 
131  const Array<OneD, unsigned int> nrow,
132  const Array<OneD, unsigned int> ncol)
133  {
134  mat =
136  SNekMatSharedPtr loc_matNvar;
137  for (int nelm = 0; nelm < nrow.size(); ++nelm)
138  {
139  int nrowsVars = nrow[nelm];
140  int ncolsVars = ncol[nelm];
141 
143  nrowsVars, ncolsVars, 0.0);
144  mat->SetBlock(nelm, nelm, loc_matNvar);
145  }
146  }
147 };
148 } // namespace Nektar
149 
150 #endif
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.
virtual void v_InitObject()
virtual void v_DoPreconCfs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacSingle
Definition: PreconCfsBRJ.h:79
TensorOfArray4D< NekSingle > m_TraceJacArraySingle
Definition: PreconCfsBRJ.h:80
Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > m_PreconMatVarsSingle
Definition: PreconCfsBRJ.h:77
static std::string className
Name of the class.
Definition: PreconCfsBRJ.h:63
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacDerivSingle
Definition: PreconCfsBRJ.h:81
Array< OneD, Array< OneD, NekSingle > > m_TraceJacDerivSignSingle
Definition: PreconCfsBRJ.h:83
void MinusOffDiag2Rhs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int nvariables, const int nCoeffs, const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, bool flagUpdateDervFlux, Array< OneD, Array< OneD, NekDouble >> &FwdFluxDeriv, Array< OneD, Array< OneD, NekDouble >> &BwdFluxDeriv, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &wspTrace, Array< OneD, Array< OneD, DataType >> &wspTraceDataType, const TensorOfArray4D< DataType > &TraceJacArray, const TensorOfArray4D< DataType > &TraceJacDerivArray, const Array< OneD, const Array< OneD, DataType >> &TraceJacDerivSign, const TensorOfArray5D< DataType > &TraceIPSymJacArray)
PrecType m_PreconMatStorage
Definition: PreconCfsBRJ.h:86
TensorOfArray4D< NekSingle > m_TraceJacDerivArraySingle
Definition: PreconCfsBRJ.h:82
PreconCfsBRJ(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
static PreconCfsOpSharedPtr create(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
Creates an instance of this class.
Definition: PreconCfsBRJ.h:52
void AllocatePreconBlkDiagCoeff(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, const int &nscale=1)
virtual void v_BuildPreconCfs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, const Array< OneD, NekDouble >> &intmp, const NekDouble time, const NekDouble lambda)
TensorOfArray5D< NekSingle > m_TraceIPSymJacArraySingle
Definition: PreconCfsBRJ.h:84
void AllocateNekBlkMatDig(SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
Definition: PreconCfsBRJ.h:130
virtual bool UpdatePreconMatCheck(const Array< OneD, const NekDouble > &res, const NekDouble dtLambda)
void PreconBlkDiag(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const TypeNekBlkMatSharedPtr &PreconMatVars, const DataType &tmpDataType)
SNekBlkMatSharedPtr m_PreconMatSingle
Definition: PreconCfsBRJ.h:78
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
std::shared_ptr< SNekBlkMat > SNekBlkMatSharedPtr
std::shared_ptr< SNekMat > SNekMatSharedPtr
PrecType
Definition: PreconCfs.h:46
std::shared_ptr< PreconCfsOp > PreconCfsOpSharedPtr
Definition: PreconCfsOp.h:152
double NekDouble