Nektar++
Expansion.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Expansion.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: Header file for Expansion routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef EXPANSION_H
36#define EXPANSION_H
37
43#include <map>
44#include <memory>
45#include <vector>
46
48{
49
50class Expansion;
51class MatrixKey;
52
54
56{
64};
65
66typedef std::shared_ptr<Expansion> ExpansionSharedPtr;
67typedef std::weak_ptr<Expansion> ExpansionWeakPtr;
68typedef std::vector<ExpansionSharedPtr> ExpansionVector;
69typedef std::map<MetricType, Array<OneD, NekDouble>> MetricMap;
70
71class Expansion : virtual public StdRegions::StdExpansion
72{
73public:
75 SpatialDomains::GeometrySharedPtr pGeom); // default constructor.
76 LOCAL_REGIONS_EXPORT Expansion(const Expansion &pSrc); // copy constructor.
78
79 LOCAL_REGIONS_EXPORT void SetTraceExp(const int traceid,
82
85
87 const LocalRegions::MatrixKey &mkey);
88
90 const StdRegions::MatrixType mtype,
94
96
98
100 CreateIndexMap(const IndexMapKey &ikey);
101
104
106 const;
107
110 const StdRegions::MatrixType matrixType);
111
114
116 const NekDouble *data, const std::vector<unsigned int> &nummodes,
117 const int nmodes_offset, NekDouble *coeffs,
118 std::vector<LibUtilities::BasisType> &fromType);
119
121 const int edge, const std::shared_ptr<Expansion> &EdgeExp,
124 Array<OneD, NekDouble> &outarray);
126 const int edge, const std::shared_ptr<Expansion> &EdgeExp,
128 Array<OneD, NekDouble> &outarray);
130 const int face, const std::shared_ptr<Expansion> &FaceExp,
132 Array<OneD, NekDouble> &outarray);
134 const int dir, const Array<OneD, const NekDouble> &inarray,
137 Array<OneD, NekDouble> &outarray);
140
143 Array<OneD, Array<OneD, NekDouble>> &d0factors,
144 Array<OneD, Array<OneD, NekDouble>> &d1factors);
145
147 {
148 return m_indexMapManager[ikey];
149 }
150
152 const int dir, const Array<OneD, const NekDouble> &inarray,
154 {
155 v_AlignVectorToCollapsedDir(dir, inarray, outarray);
156 }
157
159
161
162 inline int GetLeftAdjacentElementTrace() const;
163
164 inline int GetRightAdjacentElementTrace() const;
165
166 inline void SetAdjacentElementExp(int traceid, ExpansionSharedPtr &e);
167
169 {
170 return v_GetTraceOrient(trace);
171 }
172
175 Array<OneD, NekDouble> &outarray)
176 {
177 v_SetCoeffsToOrientation(dir, inarray, outarray);
178 }
179
180 /// Divided by the metric jacobi and quadrature weights
182 const Array<OneD, const NekDouble> &inarray,
183 Array<OneD, NekDouble> &outarray)
184 {
185 v_DivideByQuadratureMetric(inarray, outarray);
186 }
187
188 /**
189 * @brief Extract the metric factors to compute the contravariant
190 * fluxes along edge \a edge and stores them into \a outarray
191 * following the local edge orientation (i.e. anticlockwise
192 * convention).
193 */
194 inline void GetTraceQFactors(const int trace,
195 Array<OneD, NekDouble> &outarray)
196 {
197 v_GetTraceQFactors(trace, outarray);
198 }
199
200 inline void GetTracePhysVals(
201 const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp,
202 const Array<OneD, const NekDouble> &inarray,
203 Array<OneD, NekDouble> &outarray,
205 {
206 v_GetTracePhysVals(trace, TraceExp, inarray, outarray, orient);
207 }
208
209 inline void GetTracePhysMap(const int edge, Array<OneD, int> &outarray)
210 {
211 v_GetTracePhysMap(edge, outarray);
212 }
213
215 Array<OneD, int> &idmap, const int nq0,
216 const int nq1)
217 {
218 v_ReOrientTracePhysMap(orient, idmap, nq0, nq1);
219 }
220
222
223 inline void ComputeTraceNormal(const int id)
224 {
226 }
227
229 {
230 return v_GetPhysNormals();
231 }
232
234 {
235 v_SetPhysNormals(normal);
236 }
237
238 inline void SetUpPhysNormals(const int trace)
239 {
240 v_SetUpPhysNormals(trace);
241 }
242
244 const int traceid, const Array<OneD, const NekDouble> &primCoeffs,
245 DNekMatSharedPtr &inoutmat)
246 {
247 v_AddRobinMassMatrix(traceid, primCoeffs, inoutmat);
248 }
249
250 inline void TraceNormLen(const int traceid, NekDouble &h, NekDouble &p)
251 {
252 v_TraceNormLen(traceid, h, p);
253 }
254
256 const int traceid, const Array<OneD, const NekDouble> &primCoeffs,
257 const Array<OneD, NekDouble> &incoeffs, Array<OneD, NekDouble> &coeffs)
258 {
259 v_AddRobinTraceContribution(traceid, primCoeffs, incoeffs, coeffs);
260 }
261
263 GetElmtBndNormDirElmtLen(const int nbnd) const;
264
267
268protected:
271
272 std::map<int, ExpansionWeakPtr> m_traceExp;
276 std::map<int, NormalVector> m_traceNormals;
281
282 /// the element length in each element boundary(Vertex, edge
283 /// or face) normal direction calculated based on the local
284 /// m_metricinfo times the standard element length (which is
285 /// 2.0)
286 std::map<int, Array<OneD, NekDouble>> m_elmtBndNormDirElmtLen;
287
291 const Array<OneD, const NekDouble> &direction,
293
294 Array<OneD, NekDouble> GetMF(const int dir, const int shapedim,
295 const StdRegions::VarCoeffMap &varcoeffs);
296
297 Array<OneD, NekDouble> GetMFDiv(const int dir,
298 const StdRegions::VarCoeffMap &varcoeffs);
299
300 Array<OneD, NekDouble> GetMFMag(const int dir,
301 const StdRegions::VarCoeffMap &varcoeffs);
302
304 const Array<OneD, const NekDouble> &inarray,
305 Array<OneD, NekDouble> &outarray) override;
306
307 virtual void v_DivideByQuadratureMetric(
308 const Array<OneD, const NekDouble> &inarray,
309 Array<OneD, NekDouble> &outarray);
310
312 {
313 }
314
315 int v_GetCoordim() const override
316 {
317 return m_geom->GetCoordim();
318 }
319
320 void v_GetCoords(Array<OneD, NekDouble> &coords_1,
321 Array<OneD, NekDouble> &coords_2,
322 Array<OneD, NekDouble> &coords_3) override;
323
325 const LocalRegions::MatrixKey &mkey);
326
327 virtual void v_DropLocMatrix(const LocalRegions::MatrixKey &mkey);
328
330 const DNekScalMatSharedPtr &r_bnd,
331 const StdRegions::MatrixType matrixType);
332
334 const DNekScalMatSharedPtr &r_bnd);
335
336 virtual void v_ExtractDataToCoeffs(
337 const NekDouble *data, const std::vector<unsigned int> &nummodes,
338 const int nmodes_offset, NekDouble *coeffs,
339 std::vector<LibUtilities::BasisType> &fromType);
340
341 virtual void v_AddEdgeNormBoundaryInt(
342 const int edge, const std::shared_ptr<Expansion> &EdgeExp,
345 Array<OneD, NekDouble> &outarray);
346 virtual void v_AddEdgeNormBoundaryInt(
347 const int edge, const std::shared_ptr<Expansion> &EdgeExp,
349 Array<OneD, NekDouble> &outarray);
350 virtual void v_AddFaceNormBoundaryInt(
351 const int face, const std::shared_ptr<Expansion> &FaceExp,
353 Array<OneD, NekDouble> &outarray);
354 virtual void v_DGDeriv(const int dir,
355 const Array<OneD, const NekDouble> &inarray,
358 Array<OneD, NekDouble> &outarray);
359 virtual NekDouble v_VectorFlux(
360 const Array<OneD, Array<OneD, NekDouble>> &vec);
361
362 virtual void v_NormalTraceDerivFactors(
364 Array<OneD, Array<OneD, NekDouble>> &d0factors,
365 Array<OneD, Array<OneD, NekDouble>> &d1factors);
366
367 virtual void v_AlignVectorToCollapsedDir(
368 const int dir, const Array<OneD, const NekDouble> &inarray,
369 Array<OneD, Array<OneD, NekDouble>> &outarray);
370
371 virtual StdRegions::Orientation v_GetTraceOrient(int trace);
372
375 Array<OneD, NekDouble> &outarray) override;
376
377 virtual void v_GetTraceQFactors(const int trace,
378 Array<OneD, NekDouble> &outarray);
379
380 virtual void v_GetTracePhysVals(
381 const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp,
382 const Array<OneD, const NekDouble> &inarray,
384
385 virtual void v_GetTracePhysMap(const int edge, Array<OneD, int> &outarray);
386
387 virtual void v_ReOrientTracePhysMap(const StdRegions::Orientation orient,
388 Array<OneD, int> &idmap, const int nq0,
389 const int nq1 = -1);
390
391 virtual void v_ComputeTraceNormal(const int id);
392
394
396
397 virtual void v_SetUpPhysNormals(const int id);
398
399 virtual void v_AddRobinMassMatrix(
400 const int face, const Array<OneD, const NekDouble> &primCoeffs,
401 DNekMatSharedPtr &inoutmat);
402
403 virtual void v_AddRobinTraceContribution(
404 const int traceid, const Array<OneD, const NekDouble> &primCoeffs,
405 const Array<OneD, NekDouble> &incoeffs, Array<OneD, NekDouble> &coeffs);
406
407 virtual void v_TraceNormLen(const int traceid, NekDouble &h, NekDouble &p);
408
409 virtual void v_GenTraceExp(const int traceid, ExpansionSharedPtr &exp);
410
411private:
412};
413
415{
416 ASSERTL1(traceid < GetNtraces(), "Trace is out of range.");
417
418 ExpansionSharedPtr returnval;
419
420 if (m_traceExp.count(traceid))
421 {
422 // Use stored value
423 returnval = m_traceExp[traceid].lock();
424 }
425 else
426 {
427 // Generate trace exp
428 v_GenTraceExp(traceid, returnval);
429 }
430
431 return returnval;
432}
433
434inline void Expansion::SetTraceExp(const int traceid, ExpansionSharedPtr &exp)
435{
436 ASSERTL1(traceid < GetNtraces(), "Trace out of range.");
437
438 m_traceExp[traceid] = exp;
439}
440
442{
443 ASSERTL1(m_elementLeft.lock().get(), "Left adjacent element not set.");
444 return m_elementLeft.lock();
445}
446
448{
449 ASSERTL1(m_elementLeft.lock().get(), "Right adjacent element not set.");
450
451 return m_elementRight.lock();
452}
453
455{
456 return m_elementTraceLeft;
457}
458
460{
461 return m_elementTraceRight;
462}
463
464inline void Expansion::SetAdjacentElementExp(int traceid,
466{
467 if (m_elementLeft.lock().get())
468 {
469 m_elementRight = exp;
470 m_elementTraceRight = traceid;
471 }
472 else
473 {
474 m_elementLeft = exp;
475 m_elementTraceLeft = traceid;
476 }
477}
478
479} // namespace Nektar::LocalRegions
480
481#endif
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define LOCAL_REGIONS_EXPORT
std::map< int, NormalVector > m_traceNormals
Definition: Expansion.h:276
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
Definition: Expansion.h:286
virtual void v_SetPhysNormals(Array< OneD, const NekDouble > &normal)
Definition: Expansion.cpp:879
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:167
void SetTraceExp(const int traceid, ExpansionSharedPtr &f)
Definition: Expansion.h:434
virtual void v_TraceNormLen(const int traceid, NekDouble &h, NekDouble &p)
Definition: Expansion.cpp:907
IndexMapValuesSharedPtr CreateIndexMap(const IndexMapKey &ikey)
Definition: Expansion.cpp:181
void AddRobinTraceContribution(const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
Definition: Expansion.h:255
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:273
void GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Definition: Expansion.h:209
void DropLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:90
const SpatialDomains::GeomFactorsSharedPtr & GetMetricInfo() const
Definition: Expansion.cpp:246
ExpansionSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion.h:441
DNekMatSharedPtr BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:101
virtual void v_DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:434
virtual void v_GenTraceExp(const int traceid, ExpansionSharedPtr &exp)
Definition: Expansion.cpp:861
void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: Expansion.cpp:818
virtual void v_AddRobinTraceContribution(const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
Definition: Expansion.cpp:898
ExpansionWeakPtr m_elementRight
Definition: Expansion.h:278
void AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Expansion.h:151
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
Definition: Expansion.cpp:603
virtual void v_GetTraceQFactors(const int trace, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:826
virtual void v_ComputeTraceNormal(const int id)
Definition: Expansion.cpp:868
void GetTraceQFactors(const int trace, Array< OneD, NekDouble > &outarray)
Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into o...
Definition: Expansion.h:194
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
Definition: Expansion.cpp:272
void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: Expansion.cpp:530
Array< OneD, NekDouble > GetMFMag(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:709
void SetUpPhysNormals(const int trace)
Definition: Expansion.h:238
void SetAdjacentElementExp(int traceid, ExpansionSharedPtr &e)
Definition: Expansion.h:464
virtual DNekMatSharedPtr v_BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:739
const Array< OneD, const NekDouble > & GetPhysNormals(void)
Definition: Expansion.h:228
virtual DNekScalMatSharedPtr v_GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:265
int GetLeftAdjacentElementTrace() const
Definition: Expansion.h:454
void DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Divided by the metric jacobi and quadrature weights.
Definition: Expansion.h:181
void GetTracePhysVals(const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=StdRegions::eNoOrientation)
Definition: Expansion.h:200
virtual void v_AddFaceNormBoundaryInt(const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:775
virtual void v_DGDeriv(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:784
virtual DNekMatSharedPtr v_BuildTransformationMatrix(const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
Definition: Expansion.cpp:731
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen(const int nbnd) const
Definition: Expansion.cpp:914
virtual void v_AddRobinMassMatrix(const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
Definition: Expansion.cpp:890
std::map< int, ExpansionWeakPtr > m_traceExp
Definition: Expansion.h:272
virtual void v_NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)
Definition: Expansion.cpp:803
void AddEdgeNormBoundaryInt(const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:114
ExpansionWeakPtr m_elementLeft
Definition: Expansion.h:277
void AddFaceNormBoundaryInt(const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:129
virtual void v_DropLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:413
void AddRobinMassMatrix(const int traceid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
Definition: Expansion.h:243
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &vec)
Definition: Expansion.cpp:794
int v_GetCoordim() const override
Definition: Expansion.h:315
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:274
virtual StdRegions::Orientation v_GetTraceOrient(int trace)
Definition: Expansion.cpp:813
void DGDeriv(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:136
void StdDerivBaseOnTraceMat(Array< OneD, DNekMatSharedPtr > &DerivMat)
Definition: Expansion.cpp:474
void ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
Definition: Expansion.cpp:106
void v_MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: Expansion.cpp:419
ExpansionSharedPtr GetTraceExp(const int traceid)
Definition: Expansion.h:414
void SetPhysNormals(Array< OneD, const NekDouble > &normal)
Definition: Expansion.h:233
virtual void v_ReOrientTracePhysMap(const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1=-1)
Definition: Expansion.cpp:852
void ReOrientTracePhysMap(const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
Definition: Expansion.h:214
void NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)
Definition: Expansion.cpp:150
virtual void v_ComputeLaplacianMetric()
Definition: Expansion.h:311
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals()
Definition: Expansion.cpp:873
virtual void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
Definition: Expansion.cpp:746
ExpansionSharedPtr GetRightAdjacentElementExp() const
Definition: Expansion.h:447
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Expansion.cpp:923
StdRegions::Orientation GetTraceOrient(int trace)
Definition: Expansion.h:168
virtual void v_AddEdgeNormBoundaryInt(const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:756
IndexMapValuesSharedPtr GetIndexMap(const IndexMapKey &ikey)
Definition: Expansion.h:146
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:84
virtual void v_SetUpPhysNormals(const int id)
Definition: Expansion.cpp:885
virtual void v_GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Definition: Expansion.cpp:845
Array< OneD, NekDouble > GetMFDiv(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:686
int GetRightAdjacentElementTrace() const
Definition: Expansion.h:459
void TraceNormLen(const int traceid, NekDouble &h, NekDouble &p)
Definition: Expansion.h:250
const NormalVector & GetTraceNormal(const int id)
Definition: Expansion.cpp:251
Array< OneD, NekDouble > GetMF(const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:633
virtual void v_GetTracePhysVals(const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
Definition: Expansion.cpp:834
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:43
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLess > m_indexMapManager
Definition: Expansion.h:270
void ComputeTraceNormal(const int id)
Definition: Expansion.h:223
void SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: Expansion.h:173
DNekMatSharedPtr BuildTransformationMatrix(const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
Definition: Expansion.cpp:95
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &vec)
Definition: Expansion.cpp:145
The base class for all shapes.
Definition: StdExpansion.h:65
int GetNtraces() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:351
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:126
std::weak_ptr< Expansion > ExpansionWeakPtr
Definition: Expansion.h:67
Array< OneD, Array< OneD, NekDouble > > NormalVector
Definition: Expansion.h:53
std::vector< ExpansionSharedPtr > ExpansionVector
Definition: Expansion.h:68
std::map< MetricType, Array< OneD, NekDouble > > MetricMap
Definition: Expansion.h:69
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:60
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:51
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:431
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:376
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:375
StdRegions::ConstFactorMap factors
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble