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 46 of file StdMatrixKey.cpp.

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

51  :
52  m_shapeType(shapeType),
53  m_base(stdExpansion.GetBase()),
54  m_ncoeffs(stdExpansion.GetNcoeffs()),
55  m_matrixType(matrixType),
56  m_nodalPointsType(nodalType),
57  m_factors(factorMap),
58  m_varcoeffs(varCoeffMap),
59  m_varcoeff_hashes(varCoeffMap.size())
60  {
61  // Create hash
62  int i = 0;
63  for (VarCoeffMap::const_iterator x = varCoeffMap.begin(); x != varCoeffMap.end(); ++x)
64  {
65  m_varcoeff_hashes[i] = boost::hash_range(x->second.begin(), x->second.begin() + stdExpansion.GetTotPoints());
66  boost::hash_combine(m_varcoeff_hashes[i], (int)x->first);
67  i++;
68  }
69  }
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 71 of file StdMatrixKey.cpp.

72  :
73  m_shapeType(rhs.m_shapeType),
74  m_base(rhs.m_base),
75  m_ncoeffs(rhs.m_ncoeffs),
76  m_matrixType(matrixType),
77  m_nodalPointsType(rhs.m_nodalPointsType),
78  m_factors(rhs.m_factors),
79  m_varcoeffs(rhs.m_varcoeffs),
80  m_varcoeff_hashes(rhs.m_varcoeff_hashes)
81  {
82  }
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 84 of file StdMatrixKey.cpp.

84  :
85  m_shapeType(rhs.m_shapeType),
86  m_base(rhs.m_base),
87  m_ncoeffs(rhs.m_ncoeffs),
88  m_matrixType(rhs.m_matrixType),
89  m_nodalPointsType(rhs.m_nodalPointsType),
90  m_factors(rhs.m_factors),
91  m_varcoeffs(rhs.m_varcoeffs),
92  m_varcoeff_hashes(rhs.m_varcoeff_hashes)
93  {
94  }
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:218
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:218
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 102 of file StdMatrixKey.cpp.

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

Definition at line 197 of file StdMatrixKey.cpp.

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