Nektar++
Loading...
Searching...
No Matches
StdNodalTriExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: StdNodalTriExp.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 triangle routines built upon StdExpansion2D
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
37
38using namespace std;
39
40namespace Nektar::StdRegions
41{
43 const LibUtilities::BasisKey &Bb,
45 : StdExpansion(LibUtilities::StdTriData::getNumberOfCoefficients(
46 Ba.GetNumModes(), Bb.GetNumModes()),
47 2, Ba, Bb),
48 StdExpansion2D(LibUtilities::StdTriData::getNumberOfCoefficients(
49 Ba.GetNumModes(), Bb.GetNumModes()),
50 Ba, Bb),
51 StdTriExp(Ba, Bb), m_nodalPointsKey(Ba.GetNumModes(), Ntype)
52{
53 ASSERTL0(m_base[0]->GetNumModes() == m_base[1]->GetNumModes(),
54 "Nodal basis initiated with different orders in the a "
55 "and b directions");
56}
57
59{
60 return true;
61}
62
63//-------------------------------
64// Nodal basis specific routines
65//-------------------------------
66
79
80// Operate with transpose of NodalToModal transformation
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
113
115{
116 int i, j;
120
122 GetNodalPoints(r, s);
123
124 // Store the values of m_phys in a temporary array
125 int nqtot = GetTotPoints();
126 Array<OneD, NekDouble> phys(nqtot);
127
128 for (i = 0; i < m_ncoeffs; ++i)
129 {
130 // fill physical space with mode i
131 StdTriExp::v_FillMode(i, phys);
132
133 // interpolate mode i to the Nodal points 'j' and
134 // store in outarray
135 for (j = 0; j < m_ncoeffs; ++j)
136 {
137 c[0] = r[j];
138 c[1] = s[j];
139 (*Mat)(j, i) = StdExpansion2D::v_PhysEvaluate(c, phys);
140 }
141 }
142 return Mat;
143}
144
145//---------------------------------------
146// Transforms
147//---------------------------------------
148
150 Array<OneD, NekDouble> &outarray)
151{
152 v_BwdTrans_SumFac(inarray, outarray);
153}
154
156 const Array<OneD, const NekDouble> &inarray,
157 Array<OneD, NekDouble> &outarray)
158{
160 v_NodalToModal(inarray, tmp);
161 StdTriExp::v_BwdTrans_SumFac(tmp, outarray);
162}
163
165 Array<OneD, NekDouble> &outarray)
166{
167 v_IProductWRTBase(inarray, outarray);
168
169 // get Mass matrix inverse
172 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
173
174 // copy inarray in case inarray == outarray
175 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
177
178 out = (*matsys) * in;
179}
180
181//---------------------------------------
182// Inner product functions
183//---------------------------------------
184
186 const Array<OneD, const NekDouble> &inarray,
187 Array<OneD, NekDouble> &outarray)
188{
189 v_IProductWRTBase_SumFac(inarray, outarray);
190}
191
193 const Array<OneD, const NekDouble> &inarray,
194 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
195{
196 StdTriExp::v_IProductWRTBase_SumFac(inarray, outarray, multiplybyweights);
197 NodalToModalTranspose(outarray, outarray);
198}
199
201 const int dir, const Array<OneD, const NekDouble> &inarray,
202 Array<OneD, NekDouble> &outarray)
203{
204 v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
205}
206
208 const int dir, const Array<OneD, const NekDouble> &inarray,
209 Array<OneD, NekDouble> &outarray)
210{
211 StdTriExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
212 NodalToModalTranspose(outarray, outarray);
213}
214
215//---------------------------------------
216// Evaluation functions
217//---------------------------------------
218
219void StdNodalTriExp::v_FillMode(const int mode,
220 Array<OneD, NekDouble> &outarray)
221{
222 ASSERTL2(mode >= m_ncoeffs,
223 "calling argument mode is larger than total expansion order");
224
225 Vmath::Zero(m_ncoeffs, outarray, 1);
226 outarray[mode] = 1.0;
227 v_BwdTrans(outarray, outarray);
228}
229
230//---------------------------
231// Helper functions
232//---------------------------
233
235{
236 return 3 + (GetBasisNumModes(0) - 2) + 2 * (GetBasisNumModes(1) - 2);
237}
238
239//--------------------------
240// Mappings
241//--------------------------
242
243int StdNodalTriExp::v_GetVertexMap(const int localVertexId,
244 [[maybe_unused]] bool useCoeffPacking)
245{
246 ASSERTL0(localVertexId >= 0 && localVertexId <= 2,
247 "Local Vertex ID must be between 0 and 2");
248 return localVertexId;
249}
250
253 Array<OneD, int> &signarray,
254 Orientation edgeOrient, int P,
255 [[maybe_unused]] int Q)
256{
257 ASSERTL0(eid >= 0 && eid <= 2, "Local Edge ID must be between 0 and 2");
258
259 const int nEdgeCoeffs = GetTraceNcoeffs(eid);
260
261 ASSERTL0(P == -1 || P == nEdgeCoeffs,
262 "Nodal triangle not set up to deal with variable "
263 "polynomial order.");
264
265 if (maparray.size() != nEdgeCoeffs)
266 {
267 maparray = Array<OneD, unsigned int>(nEdgeCoeffs);
268 }
269
270 if (signarray.size() != nEdgeCoeffs)
271 {
272 signarray = Array<OneD, int>(nEdgeCoeffs, 1);
273 }
274 else
275 {
276 fill(signarray.data(), signarray.data() + nEdgeCoeffs, 1);
277 }
278
279 Orientation orient = edgeOrient;
280 if (eid == 2)
281 {
282 orient = orient == eForwards ? eBackwards : eForwards;
283 }
284
285 maparray[0] = eid;
286 maparray[nEdgeCoeffs - 1] = eid == 2 ? 0 : eid + 1;
287 for (int i = 2; i < nEdgeCoeffs; i++)
288 {
289 maparray[i - 1] = eid * (nEdgeCoeffs - 2) + 1 + i;
290 }
291
292 if (orient == eBackwards)
293 {
294 reverse(maparray.data(), maparray.data() + nEdgeCoeffs);
295 }
296}
297
299 const int eid, Array<OneD, unsigned int> &maparray,
300 Array<OneD, int> &signarray, const Orientation edgeOrient)
301{
302 ASSERTL0(eid >= 0 && eid <= 2, "Local Edge ID must be between 0 and 2");
303
304 const int nEdgeIntCoeffs = GetTraceNcoeffs(eid) - 2;
305
306 if (maparray.size() != nEdgeIntCoeffs)
307 {
308 maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
309 }
310
311 if (signarray.size() != nEdgeIntCoeffs)
312 {
313 signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
314 }
315 else
316 {
317 fill(signarray.data(), signarray.data() + nEdgeIntCoeffs, 1);
318 }
319
320 Orientation orient = edgeOrient;
321 if (eid == 2)
322 {
323 orient = orient == eForwards ? eBackwards : eForwards;
324 }
325
326 for (int i = 0; i < nEdgeIntCoeffs; i++)
327 {
328 maparray[i] = eid * nEdgeIntCoeffs + 3 + i;
329 }
330
331 if (orient == eBackwards)
332 {
333 reverse(maparray.data(), maparray.data() + nEdgeIntCoeffs);
334 }
335}
336
338{
339 unsigned int i;
340 if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
341 {
343 }
344
345 for (i = NumBndryCoeffs(); i < GetNcoeffs(); i++)
346 {
347 outarray[i - NumBndryCoeffs()] = i;
348 }
349}
350
352{
353 unsigned int i;
354 if (outarray.size() != NumBndryCoeffs())
355 {
357 }
358
359 for (i = 0; i < NumBndryCoeffs(); i++)
360 {
361 outarray[i] = i;
362 }
363}
364
365//---------------------------------------
366// Wrapper functions
367//---------------------------------------
368
370{
372
373 switch (mkey.GetMatrixType())
374 {
375 case eNBasisTrans:
376 Mat = GenNBasisTransMatrix();
377 break;
378 default:
380 break;
381 }
382
383 return Mat;
384}
385
390
391//---------------------------------------
392// Operator evaluation functions
393//---------------------------------------
394
396 Array<OneD, NekDouble> &outarray,
397 const StdMatrixKey &mkey)
398{
399 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
400}
401
403 const Array<OneD, const NekDouble> &inarray,
404 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
405{
407 mkey);
408}
409
411 const int k1, const int k2, const Array<OneD, const NekDouble> &inarray,
412 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
413
414{
415 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
416}
417
419 const int i, const Array<OneD, const NekDouble> &inarray,
420 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
421{
422 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
423}
424
426 const Array<OneD, const NekDouble> &inarray,
427 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
428{
430 mkey);
431}
432
433//---------------------------------------
434// Private helper functions
435//---------------------------------------
436
437} // namespace Nektar::StdRegions
#define ASSERTL0(condition, msg)
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Describes the specification for a Basis.
Definition Basis.h:45
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.
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
MatrixType GetMatrixType() const
void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
void v_NodalToModal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
void NodalToModalTranspose(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetNodalPoints(Array< OneD, const NekDouble > &x, Array< OneD, const NekDouble > &y)
void ModalToNodal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) 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.
void v_GetTraceInteriorToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation edgeOrient=eForwards) override
void v_GetTraceToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient=eForwards, int P=-1, int Q=-1) override
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.
int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
LibUtilities::PointsKey m_nodalPointsKey
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) 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
void v_FillMode(const int mode, 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 multiplybyweights=true) override
PointsManagerT & PointsManager(void)
@ P
Monomial polynomials .
Definition BasisType.h:62
static ConstFactorMap NullConstFactorMap
static VarCoeffMap NullVarCoeffMap
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekMat > DNekMatSharedPtr
void Zero(int n, T *x, const int incx)
Zero vector.
Definition Vmath.hpp:273
STL namespace.