Nektar++
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 <boost/core/ignore_unused.hpp>
36 
37 #include <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44 namespace StdRegions
45 {
46 StdNodalTriExp::StdNodalTriExp(const LibUtilities::BasisKey &Ba,
47  const LibUtilities::BasisKey &Bb,
49  : StdExpansion(LibUtilities::StdTriData::getNumberOfCoefficients(
50  Ba.GetNumModes(), Bb.GetNumModes()),
51  2, Ba, Bb),
52  StdExpansion2D(LibUtilities::StdTriData::getNumberOfCoefficients(
53  Ba.GetNumModes(), Bb.GetNumModes()),
54  Ba, Bb),
55  StdTriExp(Ba, Bb), m_nodalPointsKey(Ba.GetNumModes(), Ntype)
56 {
57  ASSERTL0(m_base[0]->GetNumModes() == m_base[1]->GetNumModes(),
58  "Nodal basis initiated with different orders in the a "
59  "and b directions");
60 }
61 
64  m_nodalPointsKey(T.m_nodalPointsKey)
65 {
66 }
67 
69 {
70  return true;
71 }
72 
73 //-------------------------------
74 // Nodal basis specific routines
75 //-------------------------------
76 
78  Array<OneD, NekDouble> &outarray)
79 {
83  DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
84 
85  NekVector<NekDouble> nodal(m_ncoeffs, inarray, eWrapper);
86  NekVector<NekDouble> modal(m_ncoeffs, outarray, eWrapper);
87  modal = (*inv_vdm) * nodal;
88 }
89 
90 // Operate with transpose of NodalToModal transformation
92  const Array<OneD, const NekDouble> &inarray,
93  Array<OneD, NekDouble> &outarray)
94 {
98  DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
99 
100  NekVector<NekDouble> nodal(m_ncoeffs, inarray, eCopy);
101  NekVector<NekDouble> modal(m_ncoeffs, outarray, eWrapper);
102  modal = Transpose(*inv_vdm) * nodal;
103 }
104 
106  Array<OneD, NekDouble> &outarray)
107 {
110  DNekMatSharedPtr vdm = GetStdMatrix(Nkey);
111 
112  // Multiply out matrix
113  NekVector<NekDouble> modal(m_ncoeffs, inarray, eWrapper);
114  NekVector<NekDouble> nodal(m_ncoeffs, outarray, eWrapper);
115  nodal = (*vdm) * modal;
116 }
117 
120 {
121  LibUtilities::PointsManager()[m_nodalPointsKey]->GetPoints(x, y);
122 }
123 
125 {
126  int i, j;
129  DNekMatSharedPtr Mat;
130 
132  GetNodalPoints(r, s);
133 
134  // Store the values of m_phys in a temporary array
135  int nqtot = GetTotPoints();
136  Array<OneD, NekDouble> phys(nqtot);
137 
138  for (i = 0; i < m_ncoeffs; ++i)
139  {
140  // fill physical space with mode i
141  StdTriExp::v_FillMode(i, phys);
142 
143  // interpolate mode i to the Nodal points 'j' and
144  // store in outarray
145  for (j = 0; j < m_ncoeffs; ++j)
146  {
147  c[0] = r[j];
148  c[1] = s[j];
149  (*Mat)(j, i) = StdExpansion2D::v_PhysEvaluate(c, phys);
150  }
151  }
152  return Mat;
153 }
154 
155 //---------------------------------------
156 // Transforms
157 //---------------------------------------
158 
160  Array<OneD, NekDouble> &outarray)
161 {
162  v_BwdTrans_SumFac(inarray, outarray);
163 }
164 
166  const Array<OneD, const NekDouble> &inarray,
167  Array<OneD, NekDouble> &outarray)
168 {
170  NodalToModal(inarray, tmp);
171  StdTriExp::v_BwdTrans_SumFac(tmp, outarray);
172 }
173 
175  Array<OneD, NekDouble> &outarray)
176 {
177  v_IProductWRTBase(inarray, outarray);
178 
179  // get Mass matrix inverse
182  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
183 
184  // copy inarray in case inarray == outarray
185  NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
186  NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
187 
188  out = (*matsys) * in;
189 }
190 
191 //---------------------------------------
192 // Inner product functions
193 //---------------------------------------
194 
196  const Array<OneD, const NekDouble> &inarray,
197  Array<OneD, NekDouble> &outarray)
198 {
199  v_IProductWRTBase_SumFac(inarray, outarray);
200 }
201 
203  const Array<OneD, const NekDouble> &inarray,
204  Array<OneD, NekDouble> &outarray, bool multiplybyweights)
205 {
206  StdTriExp::v_IProductWRTBase_SumFac(inarray, outarray, multiplybyweights);
207  NodalToModalTranspose(outarray, outarray);
208 }
209 
211  const int dir, const Array<OneD, const NekDouble> &inarray,
212  Array<OneD, NekDouble> &outarray)
213 {
214  v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
215 }
216 
218  const int dir, const Array<OneD, const NekDouble> &inarray,
219  Array<OneD, NekDouble> &outarray)
220 {
221  StdTriExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
222  NodalToModalTranspose(outarray, outarray);
223 }
224 
225 //---------------------------------------
226 // Evaluation functions
227 //---------------------------------------
228 
229 void StdNodalTriExp::v_FillMode(const int mode,
230  Array<OneD, NekDouble> &outarray)
231 {
232  ASSERTL2(mode >= m_ncoeffs,
233  "calling argument mode is larger than total expansion order");
234 
235  Vmath::Zero(m_ncoeffs, outarray, 1);
236  outarray[mode] = 1.0;
237  v_BwdTrans(outarray, outarray);
238 }
239 
240 //---------------------------
241 // Helper functions
242 //---------------------------
243 
245 {
246  return 3 + (GetBasisNumModes(0) - 2) + 2 * (GetBasisNumModes(1) - 2);
247 }
248 
249 //--------------------------
250 // Mappings
251 //--------------------------
252 
253 int StdNodalTriExp::v_GetVertexMap(const int localVertexId,
254  bool useCoeffPacking)
255 {
256  boost::ignore_unused(useCoeffPacking);
257  ASSERTL0(localVertexId >= 0 && localVertexId <= 2,
258  "Local Vertex ID must be between 0 and 2");
259  return localVertexId;
260 }
261 
263  Array<OneD, unsigned int> &maparray,
264  Array<OneD, int> &signarray,
265  Orientation edgeOrient, int P,
266  int Q)
267 {
268  boost::ignore_unused(Q);
269 
270  ASSERTL0(eid >= 0 && eid <= 2, "Local Edge ID must be between 0 and 2");
271 
272  const int nEdgeCoeffs = GetTraceNcoeffs(eid);
273 
274  ASSERTL0(P == -1 || P == nEdgeCoeffs,
275  "Nodal triangle not set up to deal with variable "
276  "polynomial order.");
277 
278  if (maparray.size() != nEdgeCoeffs)
279  {
280  maparray = Array<OneD, unsigned int>(nEdgeCoeffs);
281  }
282 
283  if (signarray.size() != nEdgeCoeffs)
284  {
285  signarray = Array<OneD, int>(nEdgeCoeffs, 1);
286  }
287  else
288  {
289  fill(signarray.get(), signarray.get() + nEdgeCoeffs, 1);
290  }
291 
292  Orientation orient = edgeOrient;
293  if (eid == 2)
294  {
295  orient = orient == eForwards ? eBackwards : eForwards;
296  }
297 
298  maparray[0] = eid;
299  maparray[nEdgeCoeffs - 1] = eid == 2 ? 0 : eid + 1;
300  for (int i = 2; i < nEdgeCoeffs; i++)
301  {
302  maparray[i - 1] = eid * (nEdgeCoeffs - 2) + 1 + i;
303  }
304 
305  if (orient == eBackwards)
306  {
307  reverse(maparray.get(), maparray.get() + nEdgeCoeffs);
308  }
309 }
310 
312  const int eid, Array<OneD, unsigned int> &maparray,
313  Array<OneD, int> &signarray, const Orientation edgeOrient)
314 {
315  ASSERTL0(eid >= 0 && eid <= 2, "Local Edge ID must be between 0 and 2");
316 
317  const int nEdgeIntCoeffs = GetTraceNcoeffs(eid) - 2;
318 
319  if (maparray.size() != nEdgeIntCoeffs)
320  {
321  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
322  }
323 
324  if (signarray.size() != nEdgeIntCoeffs)
325  {
326  signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
327  }
328  else
329  {
330  fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
331  }
332 
333  Orientation orient = edgeOrient;
334  if (eid == 2)
335  {
336  orient = orient == eForwards ? eBackwards : eForwards;
337  }
338 
339  for (int i = 0; i < nEdgeIntCoeffs; i++)
340  {
341  maparray[i] = eid * nEdgeIntCoeffs + 3 + i;
342  }
343 
344  if (orient == eBackwards)
345  {
346  reverse(maparray.get(), maparray.get() + nEdgeIntCoeffs);
347  }
348 }
349 
351 {
352  unsigned int i;
353  if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
354  {
356  }
357 
358  for (i = NumBndryCoeffs(); i < GetNcoeffs(); i++)
359  {
360  outarray[i - NumBndryCoeffs()] = i;
361  }
362 }
363 
365 {
366  unsigned int i;
367  if (outarray.size() != NumBndryCoeffs())
368  {
370  }
371 
372  for (i = 0; i < NumBndryCoeffs(); i++)
373  {
374  outarray[i] = i;
375  }
376 }
377 
378 //---------------------------------------
379 // Wrapper functions
380 //---------------------------------------
381 
383 {
384  DNekMatSharedPtr Mat;
385 
386  switch (mkey.GetMatrixType())
387  {
388  case eNBasisTrans:
389  Mat = GenNBasisTransMatrix();
390  break;
391  default:
393  break;
394  }
395 
396  return Mat;
397 }
398 
400 {
401  return StdNodalTriExp::v_GenMatrix(mkey);
402 }
403 
404 //---------------------------------------
405 // Operator evaluation functions
406 //---------------------------------------
407 
409  Array<OneD, NekDouble> &outarray,
410  const StdMatrixKey &mkey)
411 {
412  StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
413 }
414 
416  const Array<OneD, const NekDouble> &inarray,
417  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
418 {
420  mkey);
421 }
422 
424  const int k1, const int k2, const Array<OneD, const NekDouble> &inarray,
425  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
426 
427 {
428  StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
429 }
430 
432  const int i, const Array<OneD, const NekDouble> &inarray,
433  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
434 {
435  StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
436 }
437 
439  const Array<OneD, const NekDouble> &inarray,
440  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
441 {
443  mkey);
444 }
445 
446 //---------------------------------------
447 // Private helper functions
448 //---------------------------------------
449 
450 } // namespace StdRegions
451 } // 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:50
PointsType GetPointsType() const
Definition: Points.h:109
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 GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
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)
Definition: StdExpansion.h:609
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.
Definition: StdExpansion.h:373
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:267
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.
Definition: StdExpansion.h:175
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
Definition: StdMatrixKey.h:85
void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
void NodalToModal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) 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
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTriExp.cpp:246
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
Definition: StdTriExp.cpp:628
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTriExp.cpp:523
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdTriExp.cpp:450
PointsManagerT & PointsManager(void)
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:400
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:344
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:492