Nektar++
StdNodalPrismExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: StdNodalPrismExp.cpp
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: Nodal prismatic routines built upon StdExpansion3D
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <boost/core/ignore_unused.hpp>
36
37#include <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
39
40namespace Nektar
41{
42namespace StdRegions
43{
44
46{
47}
48
50 const LibUtilities::BasisKey &Bb,
51 const LibUtilities::BasisKey &Bc,
53 : StdExpansion(LibUtilities::StdPrismData::getNumberOfCoefficients(
54 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
55 3, Ba, Bb, Bc),
56 StdExpansion3D(LibUtilities::StdPrismData::getNumberOfCoefficients(
57 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
58 Ba, Bb, Bc),
59 StdPrismExp(Ba, Bb, Bc), m_nodalPointsKey(Ba.GetNumModes(), Ntype)
60{
61 ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
62 "order in 'a' direction is higher than order "
63 "in 'c' direction");
64}
65
68 m_nodalPointsKey(T.m_nodalPointsKey)
69{
70}
71
72// Destructor
74{
75}
76
78{
79 return true;
80}
81
82//-------------------------------
83// Nodal basis specific routines
84//-------------------------------
85
87 Array<OneD, NekDouble> &outarray)
88{
92 DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
93
95 NekVector<NekDouble> modal(m_ncoeffs, outarray, eWrapper);
96 modal = (*inv_vdm) * nodal;
97}
98
99// Operate with transpose of NodalToModal transformation
101 const Array<OneD, const NekDouble> &inarray,
102 Array<OneD, NekDouble> &outarray)
103{
107 DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
108
109 NekVector<NekDouble> nodal(m_ncoeffs, inarray, eCopy);
110 NekVector<NekDouble> modal(m_ncoeffs, outarray, eWrapper);
111 modal = Transpose(*inv_vdm) * nodal;
112}
113
115 Array<OneD, NekDouble> &outarray)
116{
119 DNekMatSharedPtr vdm = GetStdMatrix(Nkey);
120
121 // Multiply out matrix
122 NekVector<NekDouble> modal(m_ncoeffs, inarray, eWrapper);
123 NekVector<NekDouble> nodal(m_ncoeffs, outarray, eWrapper);
124 nodal = (*vdm) * modal;
125}
126
130{
131 // Get 3D nodal distribution.
133}
134
136{
137 int i, j;
141
143 GetNodalPoints(r, s, t);
144
145 // Store the values of m_phys in a temporary array
146 int nqtot = GetTotPoints();
147 Array<OneD, NekDouble> tmp_phys(nqtot);
148
149 for (i = 0; i < m_ncoeffs; ++i)
150 {
151 // fill physical space with mode i
152 StdPrismExp::v_FillMode(i, tmp_phys);
153
154 // interpolate mode i to the Nodal points 'j' and
155 // store in outarray
156 for (j = 0; j < m_ncoeffs; ++j)
157 {
158 c[0] = r[j];
159 c[1] = s[j];
160 c[2] = t[j];
161 (*Mat)(j, i) = StdExpansion3D::v_PhysEvaluate(c, tmp_phys);
162 }
163 }
164
165 return Mat;
166}
167
168//---------------------------------------
169// Transforms
170//---------------------------------------
171
173 Array<OneD, NekDouble> &outarray)
174{
175 v_BwdTrans_SumFac(inarray, outarray);
176}
177
179 const Array<OneD, const NekDouble> &inarray,
180 Array<OneD, NekDouble> &outarray)
181{
183 NodalToModal(inarray, tmp);
184 StdPrismExp::v_BwdTrans_SumFac(tmp, outarray);
185}
186
188 Array<OneD, NekDouble> &outarray)
189{
190 v_IProductWRTBase(inarray, outarray);
191
192 // get Mass matrix inverse
195 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
196
197 // copy inarray in case inarray == outarray
198 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
200
201 out = (*matsys) * in;
202}
203
204//---------------------------------------
205// Inner product functions
206//---------------------------------------
207
209 const Array<OneD, const NekDouble> &inarray,
210 Array<OneD, NekDouble> &outarray)
211{
212 v_IProductWRTBase_SumFac(inarray, outarray);
213}
214
216 const Array<OneD, const NekDouble> &inarray,
217 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
218{
219 StdPrismExp::v_IProductWRTBase_SumFac(inarray, outarray, multiplybyweights);
220 NodalToModalTranspose(outarray, outarray);
221}
222
224 const int dir, const Array<OneD, const NekDouble> &inarray,
225 Array<OneD, NekDouble> &outarray)
226{
227 v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
228}
229
231 const int dir, const Array<OneD, const NekDouble> &inarray,
232 Array<OneD, NekDouble> &outarray)
233{
234 StdPrismExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
235 NodalToModalTranspose(outarray, outarray);
236}
237
238//---------------------------------------
239// Evaluation functions
240//---------------------------------------
241
243 Array<OneD, NekDouble> &outarray)
244{
245 ASSERTL2(mode >= m_ncoeffs,
246 "calling argument mode is larger than total expansion order");
247
248 Vmath::Zero(m_ncoeffs, outarray, 1);
249 outarray[mode] = 1.0;
250 v_BwdTrans(outarray, outarray);
251}
252
253//---------------------------------------
254// Mapping functions
255//---------------------------------------
256
257/*
258void StdNodalTriExp::v_GetFaceToElementMap(
259 const int fid,
260 const FaceOrientation faceOrient,
261 Array<OneD, unsigned int> &maparray,
262 Array<OneD, int> &signarray,
263 int nummodesA,
264 int nummodesB)
265{
266 int P, Q, i, j, k, idx = 0, nFaceCoeffs = 0;
267
268 ASSERTL0(fid >= 0 && fid <= 3,
269 "Local face ID must be between 0 and 3");
270
271 if (nummodesA == -1)
272 {
273 switch(fid)
274 {
275 case 0:
276 nummodesA = m_base[0]->GetNumModes();
277 nummodesB = m_base[1]->GetNumModes();
278 break;
279 case 1:
280 nummodesA = m_base[0]->GetNumModes();
281 nummodesB = m_base[2]->GetNumModes();
282 break;
283 case 2:
284 case 3:
285 nummodesA = m_base[1]->GetNumModes();
286 nummodesB = m_base[2]->GetNumModes();
287 break;
288 }
289 }
290
291 P = nummodesA;
292 Q = nummodesB;
293 nFaceCoeffs = Q + ((P-1)*(1 + 2*(Q-1) - (P-1)))/2;
294
295 if (maparray.size() != nFaceCoeffs)
296 {
297 maparray = Array<OneD, unsigned int>(nFaceCoeffs);
298 }
299
300 if (signarray.size() != nFaceCoeffs)
301 {
302 signarray = Array<OneD, int>(nFaceCoeffs,1);
303 }
304 else
305 {
306 fill(signarray.get(), signarray.get()+nFaceCoeffs, 1);
307 }
308
309 switch(fid)
310 {
311 case 0:
312 // Add vertices.
313 maparray[idx++] = 0;
314 maparray[idx++] = 1;
315 maparray[idx++] = 2;
316
317 // Add edges.
318 for (i = 2; i < P; ++i)
319 {
320 maparray[idx++] = ;
321 }
322 }
323}
324*/
325
326int StdNodalPrismExp::v_GetVertexMap(const int localVertexId,
327 bool useCoeffPacking)
328{
329 boost::ignore_unused(useCoeffPacking);
330 ASSERTL0(false, "Needs setting up");
331 return localVertexId;
332}
333
335{
336 unsigned int i;
337 const unsigned int nBndryCoeff = NumBndryCoeffs();
338
339 if (outarray.size() != nBndryCoeff)
340 {
341 outarray = Array<OneD, unsigned int>(nBndryCoeff);
342 }
343
344 for (i = 0; i < nBndryCoeff; i++)
345 {
346 outarray[i] = i;
347 }
348}
349
351{
352 unsigned int i;
353 const unsigned int nBndryCoeff = NumBndryCoeffs();
354
355 if (outarray.size() != m_ncoeffs - nBndryCoeff)
356 {
357 outarray = Array<OneD, unsigned int>(m_ncoeffs - nBndryCoeff);
358 }
359
360 for (i = nBndryCoeff; i < m_ncoeffs; i++)
361 {
362 outarray[i - nBndryCoeff] = i;
363 }
364}
365
366//---------------------------------------
367// Wrapper functions
368//---------------------------------------
369
371{
373
374 switch (mkey.GetMatrixType())
375 {
376 case eNBasisTrans:
377 Mat = GenNBasisTransMatrix();
378 break;
379 default:
381 break;
382 }
383
384 return Mat;
385}
386
388{
390}
391} // namespace StdRegions
392} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
Describes the specification for a Basis.
Definition: Basis.h:47
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:76
PointsType GetPointsType() const
Definition: Points.h:95
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
This function evaluates the expansion at a single (arbitrary) point of the domain.
The base class for all shapes.
Definition: StdExpansion.h:71
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:87
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
LibUtilities::PointsKey m_nodalPointsKey
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Inner product of inarray over region with respect to the object's default expansion basis; output in ...
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool mult=true) override
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
int v_GetVertexMap(const int localVertexId, bool useCoeffPacking=false) override
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
void GetNodalPoints(Array< OneD, const NekDouble > &x, Array< OneD, const NekDouble > &y, Array< OneD, const NekDouble > &z)
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculates the inner product of a given function f with the different modes of the expansion.
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void NodalToModal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void NodalToModalTranspose(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void ModalToNodal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Class representing a prismatic element in reference space.
Definition: StdPrismExp.h:49
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
PointsManagerT & PointsManager(void)
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:409
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:353
std::vector< double > z(NPUPPER)
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487