Nektar++
StdNodalTetExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: StdNodalTetExp.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 tetrahedral 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::StdTetData::getNumberOfCoefficients(
46 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
47 3, Ba, Bb, Bc),
48 StdExpansion3D(LibUtilities::StdTetData::getNumberOfCoefficients(
49 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
50 Ba, Bb, Bc),
51 StdTetExp(Ba, Bb, Bc), m_nodalPointsKey(Ba.GetNumModes(), Ntype)
52{
53 ASSERTL0(Ba.GetNumModes() <= Bb.GetNumModes(),
54 "order in 'a' direction is higher than order "
55 "in 'b' direction");
56 ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
57 "order in 'a' direction is higher than order "
58 "in 'c' direction");
59 ASSERTL0(Bb.GetNumModes() <= Bc.GetNumModes(),
60 "order in 'b' direction is higher than order "
61 "in 'c' direction");
62}
63
65{
66 return true;
67}
68
69//-------------------------------
70// Nodal basis specific routines
71//-------------------------------
72
74 Array<OneD, NekDouble> &outarray)
75{
79 DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
80
82 NekVector<NekDouble> modal(m_ncoeffs, outarray, eWrapper);
83 modal = (*inv_vdm) * nodal;
84}
85
86// Operate with transpose of NodalToModal transformation
88 const Array<OneD, const NekDouble> &inarray,
89 Array<OneD, NekDouble> &outarray)
90{
94 DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
95
96 NekVector<NekDouble> nodal(m_ncoeffs, inarray, eCopy);
97 NekVector<NekDouble> modal(m_ncoeffs, outarray, eWrapper);
98 modal = Transpose(*inv_vdm) * nodal;
99}
100
102 Array<OneD, NekDouble> &outarray)
103{
106 DNekMatSharedPtr vdm = GetStdMatrix(Nkey);
107
108 // Multiply out matrix
109 NekVector<NekDouble> modal(m_ncoeffs, inarray, eWrapper);
110 NekVector<NekDouble> nodal(m_ncoeffs, outarray, eWrapper);
111 nodal = (*vdm) * modal;
112}
113
117{
119}
120
122{
123 int i, j;
127
129 GetNodalPoints(r, s, t);
130
131 // Store the values of m_phys in a temporary array
132 int nqtot = GetTotPoints();
133 Array<OneD, NekDouble> tmp_phys(nqtot);
134
135 for (i = 0; i < m_ncoeffs; ++i)
136 {
137 // fill physical space with mode i
138 StdTetExp::v_FillMode(i, tmp_phys);
139
140 // interpolate mode i to the Nodal points 'j' and
141 // store in outarray
142 for (j = 0; j < m_ncoeffs; ++j)
143 {
144 c[0] = r[j];
145 c[1] = s[j];
146 c[2] = t[j];
147 (*Mat)(j, i) = StdExpansion3D::v_PhysEvaluate(c, tmp_phys);
148 }
149 }
150
151 return Mat;
152}
153
154//---------------------------------------
155// Transforms
156//---------------------------------------
157
159 Array<OneD, NekDouble> &outarray)
160{
161 v_BwdTrans_SumFac(inarray, outarray);
162}
163
165 const Array<OneD, const NekDouble> &inarray,
166 Array<OneD, NekDouble> &outarray)
167{
169 NodalToModal(inarray, tmp);
170 StdTetExp::v_BwdTrans_SumFac(tmp, outarray);
171}
172
174 Array<OneD, NekDouble> &outarray)
175{
176 v_IProductWRTBase(inarray, outarray);
177
178 // get Mass matrix inverse
181 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
182
183 // copy inarray in case inarray == outarray
184 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
186
187 out = (*matsys) * in;
188}
189
190//---------------------------------------
191// Inner product functions
192//---------------------------------------
193
195 const Array<OneD, const NekDouble> &inarray,
196 Array<OneD, NekDouble> &outarray)
197{
198 v_IProductWRTBase_SumFac(inarray, outarray);
199}
200
202 const Array<OneD, const NekDouble> &inarray,
203 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
204{
205 StdTetExp::v_IProductWRTBase_SumFac(inarray, outarray, multiplybyweights);
206 NodalToModalTranspose(outarray, outarray);
207}
208
210 const int dir, const Array<OneD, const NekDouble> &inarray,
211 Array<OneD, NekDouble> &outarray)
212{
213 v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
214}
215
217 const int dir, const Array<OneD, const NekDouble> &inarray,
218 Array<OneD, NekDouble> &outarray)
219{
220 StdTetExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
221 NodalToModalTranspose(outarray, outarray);
222}
223
224//---------------------------------------
225// Evaluation functions
226//---------------------------------------
227
228void StdNodalTetExp::v_FillMode(const int mode,
229 Array<OneD, NekDouble> &outarray)
230{
231 ASSERTL2(mode >= m_ncoeffs,
232 "calling argument mode is larger than total expansion order");
233
234 Vmath::Zero(m_ncoeffs, outarray, 1);
235 outarray[mode] = 1.0;
236 v_BwdTrans(outarray, outarray);
237}
238
239//---------------------------------------
240// Mapping functions
241//---------------------------------------
242
243int StdNodalTetExp::v_GetVertexMap(const int localVertexId,
244 [[maybe_unused]] bool useCoeffPacking)
245{
246 ASSERTL0(localVertexId >= 0 && localVertexId <= 3,
247 "Local Vertex ID must be between 0 and 3");
248 return localVertexId;
249}
250
252{
253 unsigned int i;
254 const unsigned int nBndryCoeff = NumBndryCoeffs();
255
256 if (outarray.size() != nBndryCoeff)
257 {
258 outarray = Array<OneD, unsigned int>(nBndryCoeff);
259 }
260
261 for (i = 0; i < nBndryCoeff; i++)
262 {
263 outarray[i] = i;
264 }
265}
266
268{
269 unsigned int i;
270 const unsigned int nBndryCoeff = NumBndryCoeffs();
271
272 if (outarray.size() != m_ncoeffs - nBndryCoeff)
273 {
274 outarray = Array<OneD, unsigned int>(m_ncoeffs - nBndryCoeff);
275 }
276
277 for (i = nBndryCoeff; i < m_ncoeffs; i++)
278 {
279 outarray[i - nBndryCoeff] = i;
280 }
281}
282
283//---------------------------------------
284// Wrapper functions
285//---------------------------------------
286
288{
290
291 switch (mkey.GetMatrixType())
292 {
293 case eNBasisTrans:
294 Mat = GenNBasisTransMatrix();
295 break;
296 default:
298 break;
299 }
300
301 return Mat;
302}
303
305{
306 return StdNodalTetExp::v_GenMatrix(mkey);
307}
308} // 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
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:83
void NodalToModalTranspose(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
int v_GetVertexMap(const int localVertexId, bool useCoeffPacking=false) override
void NodalToModal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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.
void ModalToNodal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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 GetNodalPoints(Array< OneD, const NekDouble > &x, Array< OneD, const NekDouble > &y, Array< OneD, const NekDouble > &z)
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTDerivBase_SumFac(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
StdNodalTetExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, const LibUtilities::PointsType Ntype)
LibUtilities::PointsKey m_nodalPointsKey
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:808
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:617
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdTetExp.cpp:493
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:293
LibUtilities::ShapeType DetShapeType() const
Definition: StdTetExp.h:56
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