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 <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
37
38namespace Nektar::StdRegions
39{
40
42 const LibUtilities::BasisKey &Bb,
43 const LibUtilities::BasisKey &Bc,
45 : StdExpansion(LibUtilities::StdPrismData::getNumberOfCoefficients(
46 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
47 3, Ba, Bb, Bc),
48 StdExpansion3D(LibUtilities::StdPrismData::getNumberOfCoefficients(
49 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
50 Ba, Bb, Bc),
51 StdPrismExp(Ba, Bb, Bc), m_nodalPointsKey(Ba.GetNumModes(), Ntype)
52{
53 ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
54 "order in 'a' direction is higher than order "
55 "in 'c' direction");
56}
57
59{
60 return true;
61}
62
63//-------------------------------
64// Nodal basis specific routines
65//-------------------------------
66
68 Array<OneD, NekDouble> &outarray)
69{
73 DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
74
76 NekVector<NekDouble> modal(m_ncoeffs, outarray, eWrapper);
77 modal = (*inv_vdm) * nodal;
78}
79
80// Operate with transpose of NodalToModal transformation
82 const Array<OneD, const NekDouble> &inarray,
83 Array<OneD, NekDouble> &outarray)
84{
88 DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
89
90 NekVector<NekDouble> nodal(m_ncoeffs, inarray, eCopy);
91 NekVector<NekDouble> modal(m_ncoeffs, outarray, eWrapper);
92 modal = Transpose(*inv_vdm) * nodal;
93}
94
96 Array<OneD, NekDouble> &outarray)
97{
100 DNekMatSharedPtr vdm = GetStdMatrix(Nkey);
101
102 // Multiply out matrix
103 NekVector<NekDouble> modal(m_ncoeffs, inarray, eWrapper);
104 NekVector<NekDouble> nodal(m_ncoeffs, outarray, eWrapper);
105 nodal = (*vdm) * modal;
106}
107
111{
112 // Get 3D nodal distribution.
114}
115
117{
118 int i, j;
122
124 GetNodalPoints(r, s, t);
125
126 // Store the values of m_phys in a temporary array
127 int nqtot = GetTotPoints();
128 Array<OneD, NekDouble> tmp_phys(nqtot);
129
130 for (i = 0; i < m_ncoeffs; ++i)
131 {
132 // fill physical space with mode i
133 StdPrismExp::v_FillMode(i, tmp_phys);
134
135 // interpolate mode i to the Nodal points 'j' and
136 // store in outarray
137 for (j = 0; j < m_ncoeffs; ++j)
138 {
139 c[0] = r[j];
140 c[1] = s[j];
141 c[2] = t[j];
142 (*Mat)(j, i) = StdExpansion3D::v_PhysEvaluate(c, tmp_phys);
143 }
144 }
145
146 return Mat;
147}
148
149//---------------------------------------
150// Transforms
151//---------------------------------------
152
154 Array<OneD, NekDouble> &outarray)
155{
156 v_BwdTrans_SumFac(inarray, outarray);
157}
158
160 const Array<OneD, const NekDouble> &inarray,
161 Array<OneD, NekDouble> &outarray)
162{
164 NodalToModal(inarray, tmp);
165 StdPrismExp::v_BwdTrans_SumFac(tmp, outarray);
166}
167
169 Array<OneD, NekDouble> &outarray)
170{
171 v_IProductWRTBase(inarray, outarray);
172
173 // get Mass matrix inverse
176 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
177
178 // copy inarray in case inarray == outarray
179 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
181
182 out = (*matsys) * in;
183}
184
185//---------------------------------------
186// Inner product functions
187//---------------------------------------
188
190 const Array<OneD, const NekDouble> &inarray,
191 Array<OneD, NekDouble> &outarray)
192{
193 v_IProductWRTBase_SumFac(inarray, outarray);
194}
195
197 const Array<OneD, const NekDouble> &inarray,
198 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
199{
200 StdPrismExp::v_IProductWRTBase_SumFac(inarray, outarray, multiplybyweights);
201 NodalToModalTranspose(outarray, outarray);
202}
203
205 const int dir, const Array<OneD, const NekDouble> &inarray,
206 Array<OneD, NekDouble> &outarray)
207{
208 v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
209}
210
212 const int dir, const Array<OneD, const NekDouble> &inarray,
213 Array<OneD, NekDouble> &outarray)
214{
215 StdPrismExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
216 NodalToModalTranspose(outarray, outarray);
217}
218
219//---------------------------------------
220// Evaluation functions
221//---------------------------------------
222
224 Array<OneD, NekDouble> &outarray)
225{
226 ASSERTL2(mode >= m_ncoeffs,
227 "calling argument mode is larger than total expansion order");
228
229 Vmath::Zero(m_ncoeffs, outarray, 1);
230 outarray[mode] = 1.0;
231 v_BwdTrans(outarray, outarray);
232}
233
234//---------------------------------------
235// Mapping functions
236//---------------------------------------
237
238int StdNodalPrismExp::v_GetVertexMap(const int localVertexId,
239 [[maybe_unused]] bool useCoeffPacking)
240{
241 ASSERTL0(false, "Needs setting up");
242 return localVertexId;
243}
244
246{
247 unsigned int i;
248 const unsigned int nBndryCoeff = NumBndryCoeffs();
249
250 if (outarray.size() != nBndryCoeff)
251 {
252 outarray = Array<OneD, unsigned int>(nBndryCoeff);
253 }
254
255 for (i = 0; i < nBndryCoeff; i++)
256 {
257 outarray[i] = i;
258 }
259}
260
262{
263 unsigned int i;
264 const unsigned int nBndryCoeff = NumBndryCoeffs();
265
266 if (outarray.size() != m_ncoeffs - nBndryCoeff)
267 {
268 outarray = Array<OneD, unsigned int>(m_ncoeffs - nBndryCoeff);
269 }
270
271 for (i = nBndryCoeff; i < m_ncoeffs; i++)
272 {
273 outarray[i - nBndryCoeff] = i;
274 }
275}
276
277//---------------------------------------
278// Wrapper functions
279//---------------------------------------
280
282{
284
285 switch (mkey.GetMatrixType())
286 {
287 case eNBasisTrans:
288 Mat = GenNBasisTransMatrix();
289 break;
290 default:
292 break;
293 }
294
295 return Mat;
296}
297
299{
301}
302} // namespace Nektar::StdRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
Describes the specification for a Basis.
Definition: Basis.h:45
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:74
PointsType GetPointsType() const
Definition: Points.h:90
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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:65
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
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:367
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:83
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
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
Transform a given function from physical quadrature space to coefficient space.
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)
StdNodalPrismExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, const LibUtilities::PointsType Ntype)
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:45
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
PointsManagerT & PointsManager(void)
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:431
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:376
std::vector< double > z(NPUPPER)
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.hpp:273