Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Classes | Public Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::StdRegions::StdMatrixKey Class Reference

#include <StdMatrixKey.h>

Inheritance diagram for Nektar::StdRegions::StdMatrixKey:
Inheritance graph
[legend]
Collaboration diagram for Nektar::StdRegions::StdMatrixKey:
Collaboration graph
[legend]

Classes

struct  opLess
 Used to lookup the create function in NekManager. More...
 

Public Member Functions

 StdMatrixKey (const StdRegions::MatrixType matrixType, const LibUtilities::ShapeType shapeType, const StdRegions::StdExpansion &stdExpansion, const ConstFactorMap &factorMap=NullConstFactorMap, const VarCoeffMap &varCoeffMap=NullVarCoeffMap, LibUtilities::PointsType nodalType=LibUtilities::eNoPointsType)
 
 StdMatrixKey (const StdMatrixKey &rhs, const StdRegions::MatrixType matrixType)
 
 StdMatrixKey (const StdMatrixKey &rhs)
 
virtual ~StdMatrixKey ()
 
MatrixType GetMatrixType () const
 
LibUtilities::ShapeType GetShapeType () const
 
LibUtilities::PointsType GetNodalPointsType () const
 
int GetNcoeffs () const
 
const Array< OneD, const
LibUtilities::BasisSharedPtr > & 
GetBase () const
 
std::vector< std::size_t > GetVarCoeffHashes () const
 
const LibUtilities::BasisSharedPtr GetBasis (int dir) const
 
int GetNConstFactors () const
 
NekDouble GetConstFactor (const ConstFactorType &factor) const
 
bool ConstFactorExists (const ConstFactorType &factor) const
 
const ConstFactorMapGetConstFactors () const
 
int GetNVarCoeff () const
 
const Array< OneD, const
NekDouble > & 
GetVarCoeff (const StdRegions::VarCoeffType &coeff) const
 
const VarCoeffMap GetVarCoeffAsMap (const VarCoeffType &coeff) const
 
const VarCoeffMapGetVarCoeffs () const
 
bool HasVarCoeff (const StdRegions::VarCoeffType &coeff) const
 

Protected Attributes

LibUtilities::ShapeType m_shapeType
 
Array< OneD, const
LibUtilities::BasisSharedPtr
m_base
 
unsigned int m_ncoeffs
 
MatrixType m_matrixType
 
LibUtilities::PointsType m_nodalPointsType
 
ConstFactorMap m_factors
 
VarCoeffMap m_varcoeffs
 
std::vector< std::size_t > m_varcoeff_hashes
 

Private Member Functions

 StdMatrixKey ()
 

Friends

bool operator< (const StdMatrixKey &lhs, const StdMatrixKey &rhs)
 Used for finding value given the key in NekManager. More...
 
bool operator== (const StdMatrixKey &lhs, const StdMatrixKey &rhs)
 
bool opLess::operator() (const StdMatrixKey &lhs, const StdMatrixKey &rhs) const
 

Detailed Description

Definition at line 52 of file StdMatrixKey.h.

Constructor & Destructor Documentation

Nektar::StdRegions::StdMatrixKey::StdMatrixKey ( const StdRegions::MatrixType  matrixType,
const LibUtilities::ShapeType  shapeType,
const StdRegions::StdExpansion stdExpansion,
const ConstFactorMap factorMap = NullConstFactorMap,
const VarCoeffMap varCoeffMap = NullVarCoeffMap,
LibUtilities::PointsType  nodalType = LibUtilities::eNoPointsType 
)

Definition at line 44 of file StdMatrixKey.cpp.

References Nektar::StdRegions::StdExpansion::GetTotPoints(), and m_varcoeff_hashes.

49  :
50  m_shapeType(shapeType),
51  m_base(stdExpansion.GetBase()),
52  m_ncoeffs(stdExpansion.GetNcoeffs()),
53  m_matrixType(matrixType),
54  m_nodalPointsType(nodalType),
55  m_factors(factorMap),
56  m_varcoeffs(varCoeffMap),
57  m_varcoeff_hashes(varCoeffMap.size())
58  {
59  // Create hash
60  int i = 0;
61  for (VarCoeffMap::const_iterator x = varCoeffMap.begin(); x != varCoeffMap.end(); ++x)
62  {
63  m_varcoeff_hashes[i] = boost::hash_range(x->second.begin(), x->second.begin() + stdExpansion.GetTotPoints());
64  boost::hash_combine(m_varcoeff_hashes[i], (int)x->first);
65  i++;
66  }
67  }
std::vector< std::size_t > m_varcoeff_hashes
Definition: StdMatrixKey.h:189
LibUtilities::PointsType m_nodalPointsType
Definition: StdMatrixKey.h:184
Array< OneD, const LibUtilities::BasisSharedPtr > m_base
Definition: StdMatrixKey.h:180
LibUtilities::ShapeType m_shapeType
Definition: StdMatrixKey.h:179
Nektar::StdRegions::StdMatrixKey::StdMatrixKey ( const StdMatrixKey rhs,
const StdRegions::MatrixType  matrixType 
)

Definition at line 69 of file StdMatrixKey.cpp.

70  :
71  m_shapeType(rhs.m_shapeType),
72  m_base(rhs.m_base),
73  m_ncoeffs(rhs.m_ncoeffs),
74  m_matrixType(matrixType),
75  m_nodalPointsType(rhs.m_nodalPointsType),
76  m_factors(rhs.m_factors),
77  m_varcoeffs(rhs.m_varcoeffs),
78  m_varcoeff_hashes(rhs.m_varcoeff_hashes)
79  {
80  }
std::vector< std::size_t > m_varcoeff_hashes
Definition: StdMatrixKey.h:189
LibUtilities::PointsType m_nodalPointsType
Definition: StdMatrixKey.h:184
Array< OneD, const LibUtilities::BasisSharedPtr > m_base
Definition: StdMatrixKey.h:180
LibUtilities::ShapeType m_shapeType
Definition: StdMatrixKey.h:179
Nektar::StdRegions::StdMatrixKey::StdMatrixKey ( const StdMatrixKey rhs)

Definition at line 82 of file StdMatrixKey.cpp.

82  :
83  m_shapeType(rhs.m_shapeType),
84  m_base(rhs.m_base),
85  m_ncoeffs(rhs.m_ncoeffs),
86  m_matrixType(rhs.m_matrixType),
87  m_nodalPointsType(rhs.m_nodalPointsType),
88  m_factors(rhs.m_factors),
89  m_varcoeffs(rhs.m_varcoeffs),
90  m_varcoeff_hashes(rhs.m_varcoeff_hashes)
91  {
92  }
std::vector< std::size_t > m_varcoeff_hashes
Definition: StdMatrixKey.h:189
LibUtilities::PointsType m_nodalPointsType
Definition: StdMatrixKey.h:184
Array< OneD, const LibUtilities::BasisSharedPtr > m_base
Definition: StdMatrixKey.h:180
LibUtilities::ShapeType m_shapeType
Definition: StdMatrixKey.h:179
virtual Nektar::StdRegions::StdMatrixKey::~StdMatrixKey ( )
inlinevirtual

Definition at line 67 of file StdMatrixKey.h.

68  {
69  }
Nektar::StdRegions::StdMatrixKey::StdMatrixKey ( )
private

Member Function Documentation

bool Nektar::StdRegions::StdMatrixKey::ConstFactorExists ( const ConstFactorType factor) const
inline
const Array<OneD, const LibUtilities::BasisSharedPtr>& Nektar::StdRegions::StdMatrixKey::GetBase ( ) const
inline

Definition at line 102 of file StdMatrixKey.h.

References m_base.

Referenced by Nektar::StdRegions::operator<<().

103  {
104  return m_base;
105  }
Array< OneD, const LibUtilities::BasisSharedPtr > m_base
Definition: StdMatrixKey.h:180
const LibUtilities::BasisSharedPtr Nektar::StdRegions::StdMatrixKey::GetBasis ( int  dir) const
inline

Definition at line 112 of file StdMatrixKey.h.

References m_base.

113  {
114  return(m_base[dir]);
115  }
Array< OneD, const LibUtilities::BasisSharedPtr > m_base
Definition: StdMatrixKey.h:180
NekDouble Nektar::StdRegions::StdMatrixKey::GetConstFactor ( const ConstFactorType factor) const
inline

Definition at line 122 of file StdMatrixKey.h.

References ASSERTL1, Nektar::StdRegions::ConstFactorTypeMap, and m_factors.

Referenced by Nektar::LocalRegions::PyrExp::CreateMatrix(), Nektar::LocalRegions::NodalTriExp::CreateMatrix(), Nektar::LocalRegions::TetExp::CreateMatrix(), Nektar::LocalRegions::PrismExp::CreateMatrix(), Nektar::LocalRegions::QuadExp::CreateMatrix(), Nektar::LocalRegions::TriExp::CreateMatrix(), Nektar::LocalRegions::SegExp::CreateMatrix(), Nektar::LocalRegions::HexExp::CreateMatrix(), Nektar::StdRegions::StdExpansion::HelmholtzMatrixOp_MatFree_GenericImpl(), Nektar::StdRegions::StdExpansion::LinearAdvectionDiffusionReactionMatrixOp_MatFree(), Nektar::LocalRegions::Expansion1D::v_GenMatrix(), Nektar::LocalRegions::Expansion3D::v_GenMatrix(), Nektar::LocalRegions::Expansion2D::v_GenMatrix(), Nektar::StdRegions::StdQuadExp::v_GenMatrix(), Nektar::StdRegions::StdTriExp::v_GenMatrix(), Nektar::StdRegions::StdPrismExp::v_GenMatrix(), Nektar::StdRegions::StdTetExp::v_GenMatrix(), Nektar::StdRegions::StdSegExp::v_HelmholtzMatrixOp(), Nektar::LocalRegions::SegExp::v_HelmholtzMatrixOp(), Nektar::StdRegions::StdExpansion2D::v_HelmholtzMatrixOp_MatFree(), Nektar::StdRegions::StdExpansion3D::v_HelmholtzMatrixOp_MatFree(), Nektar::StdRegions::StdTriExp::v_SVVLaplacianFilter(), Nektar::StdRegions::StdPrismExp::v_SVVLaplacianFilter(), Nektar::StdRegions::StdQuadExp::v_SVVLaplacianFilter(), Nektar::StdRegions::StdTetExp::v_SVVLaplacianFilter(), and Nektar::StdRegions::StdHexExp::v_SVVLaplacianFilter().

123  {
124  ConstFactorMap::const_iterator x = m_factors.find(factor);
125  ASSERTL1(x != m_factors.end(),
126  "Constant factor not defined: "
127  + std::string(StdRegions::ConstFactorTypeMap[factor]));
128  return x->second;
129  }
const char *const ConstFactorTypeMap[]
Definition: StdRegions.hpp:241
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
const ConstFactorMap& Nektar::StdRegions::StdMatrixKey::GetConstFactors ( ) const
inline
MatrixType Nektar::StdRegions::StdMatrixKey::GetMatrixType ( ) const
inline

Definition at line 82 of file StdMatrixKey.h.

References m_matrixType.

Referenced by Nektar::StdRegions::StdExpansion::CreateGeneralMatrix(), Nektar::LocalRegions::PyrExp::CreateMatrix(), Nektar::LocalRegions::NodalTriExp::CreateMatrix(), Nektar::LocalRegions::TetExp::CreateMatrix(), Nektar::LocalRegions::PrismExp::CreateMatrix(), Nektar::LocalRegions::QuadExp::CreateMatrix(), Nektar::LocalRegions::TriExp::CreateMatrix(), Nektar::LocalRegions::SegExp::CreateMatrix(), Nektar::LocalRegions::HexExp::CreateMatrix(), Nektar::LocalRegions::PyrExp::CreateStaticCondMatrix(), Nektar::LocalRegions::NodalTriExp::CreateStaticCondMatrix(), Nektar::LocalRegions::TetExp::CreateStaticCondMatrix(), Nektar::LocalRegions::PrismExp::CreateStaticCondMatrix(), Nektar::LocalRegions::QuadExp::CreateStaticCondMatrix(), Nektar::LocalRegions::TriExp::CreateStaticCondMatrix(), Nektar::LocalRegions::SegExp::CreateStaticCondMatrix(), Nektar::LocalRegions::HexExp::CreateStaticCondMatrix(), Nektar::StdRegions::StdExpansion::GeneralMatrixOp(), Nektar::StdRegions::StdExpansion::GeneralMatrixOp_MatFree(), Nektar::LocalRegions::MatrixKey::opLess::operator()(), Nektar::StdRegions::operator<<(), Nektar::LocalRegions::Expansion1D::v_GenMatrix(), Nektar::StdRegions::StdPointExp::v_GenMatrix(), Nektar::LocalRegions::Expansion3D::v_GenMatrix(), Nektar::LocalRegions::PyrExp::v_GenMatrix(), Nektar::LocalRegions::Expansion2D::v_GenMatrix(), Nektar::StdRegions::StdNodalPrismExp::v_GenMatrix(), Nektar::StdRegions::StdNodalTetExp::v_GenMatrix(), Nektar::LocalRegions::TetExp::v_GenMatrix(), Nektar::LocalRegions::NodalTriExp::v_GenMatrix(), Nektar::LocalRegions::PrismExp::v_GenMatrix(), Nektar::StdRegions::StdNodalTriExp::v_GenMatrix(), Nektar::StdRegions::StdQuadExp::v_GenMatrix(), Nektar::LocalRegions::TriExp::v_GenMatrix(), Nektar::LocalRegions::QuadExp::v_GenMatrix(), Nektar::StdRegions::StdTriExp::v_GenMatrix(), Nektar::StdRegions::StdSegExp::v_GenMatrix(), Nektar::LocalRegions::SegExp::v_GenMatrix(), Nektar::StdRegions::StdPrismExp::v_GenMatrix(), Nektar::StdRegions::StdTetExp::v_GenMatrix(), and Nektar::LocalRegions::HexExp::v_GenMatrix().

83  {
84  return m_matrixType;
85  }
int Nektar::StdRegions::StdMatrixKey::GetNcoeffs ( void  ) const
inline

Definition at line 97 of file StdMatrixKey.h.

References m_ncoeffs.

Referenced by Nektar::StdRegions::operator<<().

98  {
99  return m_ncoeffs;
100  }
int Nektar::StdRegions::StdMatrixKey::GetNConstFactors ( ) const
inline

Definition at line 117 of file StdMatrixKey.h.

References m_factors.

118  {
119  return m_factors.size();
120  }
LibUtilities::PointsType Nektar::StdRegions::StdMatrixKey::GetNodalPointsType ( ) const
inline

Definition at line 92 of file StdMatrixKey.h.

References m_nodalPointsType.

Referenced by Nektar::StdRegions::StdExpansion::CreateGeneralMatrix().

93  {
94  return m_nodalPointsType;
95  }
LibUtilities::PointsType m_nodalPointsType
Definition: StdMatrixKey.h:184
int Nektar::StdRegions::StdMatrixKey::GetNVarCoeff ( ) const
inline
LibUtilities::ShapeType Nektar::StdRegions::StdMatrixKey::GetShapeType ( ) const
inline
const Array<OneD, const NekDouble>& Nektar::StdRegions::StdMatrixKey::GetVarCoeff ( const StdRegions::VarCoeffType coeff) const
inline

Definition at line 152 of file StdMatrixKey.h.

References ASSERTL1, m_varcoeffs, and Nektar::StdRegions::VarCoeffTypeMap.

Referenced by GetVarCoeffAsMap(), Nektar::StdRegions::StdExpansion::LaplacianMatrixOp_MatFree(), Nektar::StdRegions::StdExpansion::LinearAdvectionDiffusionReactionMatrixOp_MatFree(), Nektar::StdRegions::StdExpansion::MassMatrixOp_MatFree(), Nektar::LocalRegions::Expansion2D::v_GenMatrix(), Nektar::StdRegions::StdTriExp::v_SVVLaplacianFilter(), Nektar::StdRegions::StdQuadExp::v_SVVLaplacianFilter(), and Nektar::StdRegions::StdExpansion::WeakDerivMatrixOp_MatFree().

153  {
154  VarCoeffMap::const_iterator x = m_varcoeffs.find(coeff);
155  ASSERTL1(x != m_varcoeffs.end(),
156  "Variable coefficient not defined: "
157  + std::string(StdRegions::VarCoeffTypeMap[coeff]));
158  return x->second;
159  }
const char *const VarCoeffTypeMap[]
Definition: StdRegions.hpp:213
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
const VarCoeffMap Nektar::StdRegions::StdMatrixKey::GetVarCoeffAsMap ( const VarCoeffType coeff) const
inline

Definition at line 161 of file StdMatrixKey.h.

References GetVarCoeff().

Referenced by Nektar::LocalRegions::Expansion2D::v_GenMatrix().

162  {
163  VarCoeffMap m;
164  m[coeff] = GetVarCoeff(coeff);
165  return m;
166  }
const Array< OneD, const NekDouble > & GetVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:152
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:226
std::vector<std::size_t> Nektar::StdRegions::StdMatrixKey::GetVarCoeffHashes ( ) const
inline

Definition at line 107 of file StdMatrixKey.h.

References m_varcoeff_hashes.

Referenced by Nektar::StdRegions::operator<<().

108  {
109  return m_varcoeff_hashes;
110  }
std::vector< std::size_t > m_varcoeff_hashes
Definition: StdMatrixKey.h:189
const VarCoeffMap& Nektar::StdRegions::StdMatrixKey::GetVarCoeffs ( ) const
inline
bool Nektar::StdRegions::StdMatrixKey::HasVarCoeff ( const StdRegions::VarCoeffType coeff) const
inline

Friends And Related Function Documentation

bool operator< ( const StdMatrixKey lhs,
const StdMatrixKey rhs 
)
friend

Used for finding value given the key in NekManager.

Definition at line 100 of file StdMatrixKey.cpp.

101  {
102  if(lhs.m_matrixType < rhs.m_matrixType)
103  {
104  return true;
105  }
106 
107  if(lhs.m_matrixType > rhs.m_matrixType)
108  {
109  return false;
110  }
111 
112  if(lhs.m_ncoeffs < rhs.m_ncoeffs) // probably do not need this check since checking the m_base below?
113  {
114  return true;
115  }
116 
117  if(lhs.m_ncoeffs > rhs.m_ncoeffs)
118  {
119  return false;
120  }
121 
122  for(unsigned int i = 0; i < LibUtilities::ShapeTypeDimMap[lhs.m_shapeType]; ++i)
123  {
124  if(lhs.m_base[i].get() < rhs.m_base[i].get())
125  {
126  return true;
127  }
128 
129  if(lhs.m_base[i].get() > rhs.m_base[i].get())
130  {
131  return false;
132  }
133  }
134 
135  if(lhs.m_factors.size() < rhs.m_factors.size())
136  {
137  return true;
138  }
139  else if(lhs.m_factors.size() > rhs.m_factors.size())
140  {
141  return false;
142  }
143  else
144  {
145  ConstFactorMap::const_iterator x, y;
146  for(x = lhs.m_factors.begin(), y = rhs.m_factors.begin();
147  x != lhs.m_factors.end(); ++x, ++y)
148  {
149  if (x->second < y->second)
150  {
151  return true;
152  }
153  if (x->second > y->second)
154  {
155  return false;
156  }
157  }
158  }
159 
160  if(lhs.m_varcoeffs.size() < rhs.m_varcoeffs.size())
161  {
162  return true;
163  }
164 
165  if(lhs.m_varcoeffs.size() > rhs.m_varcoeffs.size())
166  {
167  return false;
168  }
169 
170  for (unsigned int i = 0; i < lhs.m_varcoeff_hashes.size(); ++i)
171  {
172  if(lhs.m_varcoeff_hashes[i] < rhs.m_varcoeff_hashes[i])
173  {
174  return true;
175  }
176  if(lhs.m_varcoeff_hashes[i] > rhs.m_varcoeff_hashes[i])
177  {
178  return false;
179  }
180  }
181 
182  if(lhs.m_nodalPointsType < rhs.m_nodalPointsType)
183  {
184  return true;
185  }
186 
187  if(lhs.m_nodalPointsType > rhs.m_nodalPointsType)
188  {
189  return false;
190  }
191 
192  return false;
193  }
StandardMatrixTag & lhs
const unsigned int ShapeTypeDimMap[SIZE_ShapeType]
Definition: ShapeType.hpp:81
bool operator== ( const StdMatrixKey lhs,
const StdMatrixKey rhs 
)
friend

Definition at line 195 of file StdMatrixKey.cpp.

196  {
197  if(lhs.m_matrixType != rhs.m_matrixType)
198  {
199  return false;
200  }
201 
202  if(lhs.m_ncoeffs != rhs.m_ncoeffs)
203  {
204  return false;
205  }
206 
207  for(unsigned int i = 0; i < LibUtilities::ShapeTypeDimMap[lhs.m_shapeType]; ++i)
208  {
209  if(lhs.m_base[i].get() != rhs.m_base[i].get())
210  {
211  return false;
212  }
213  }
214 
215  if(lhs.m_factors.size() != rhs.m_factors.size())
216  {
217  return false;
218  }
219  else
220  {
221  ConstFactorMap::const_iterator x, y;
222  for(x = lhs.m_factors.begin(), y = rhs.m_factors.begin();
223  x != lhs.m_factors.end(); ++x, ++y)
224  {
225  if (x->second != y->second)
226  {
227  return false;
228  }
229  }
230  }
231 
232  if(lhs.m_nodalPointsType != rhs.m_nodalPointsType)
233  {
234  return false;
235  }
236 
237  if(lhs.m_varcoeffs.size() != rhs.m_varcoeffs.size())
238  {
239  return false;
240  }
241 
242  for (unsigned int i = 0; i < lhs.m_varcoeff_hashes.size(); ++i)
243  {
244  if(lhs.m_varcoeff_hashes[i] != rhs.m_varcoeff_hashes[i])
245  {
246  return false;
247  }
248  }
249 
250  VarCoeffMap::const_iterator x;
251  for (x = lhs.m_varcoeffs.begin(); x != lhs.m_varcoeffs.end(); ++x)
252  {
253  VarCoeffMap::const_iterator y;
254  // Check var coeff is found
255  if ((y = rhs.m_varcoeffs.find(x->first)) == rhs.m_varcoeffs.end())
256  {
257  return false;
258  }
259 
260  if (x->second != y->second)
261  {
262  return false;
263  }
264  }
265  for (unsigned int i = 0; i < lhs.m_varcoeffs.size(); ++i)
266  {
267  if(lhs.m_varcoeff_hashes[i] != rhs.m_varcoeff_hashes[i])
268  {
269  return false;
270  }
271  }
272 
273  return true;
274  }
StandardMatrixTag & lhs
const unsigned int ShapeTypeDimMap[SIZE_ShapeType]
Definition: ShapeType.hpp:81
bool opLess::operator() ( const StdMatrixKey lhs,
const StdMatrixKey rhs 
) const
friend

Member Data Documentation

Array<OneD, const LibUtilities::BasisSharedPtr> Nektar::StdRegions::StdMatrixKey::m_base
protected
ConstFactorMap Nektar::StdRegions::StdMatrixKey::m_factors
protected
MatrixType Nektar::StdRegions::StdMatrixKey::m_matrixType
protected
unsigned int Nektar::StdRegions::StdMatrixKey::m_ncoeffs
protected
LibUtilities::PointsType Nektar::StdRegions::StdMatrixKey::m_nodalPointsType
protected
LibUtilities::ShapeType Nektar::StdRegions::StdMatrixKey::m_shapeType
protected
std::vector<std::size_t> Nektar::StdRegions::StdMatrixKey::m_varcoeff_hashes
protected
VarCoeffMap Nektar::StdRegions::StdMatrixKey::m_varcoeffs
protected