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 
38 #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():
47  StdTriExp(),
48  m_nodalPointsKey()
49  {
50  }
51 
53  const LibUtilities::BasisKey &Ba,
54  const LibUtilities::BasisKey &Bb,
56  StdExpansion (LibUtilities::StdTriData::getNumberOfCoefficients(
57  Ba.GetNumModes(),
58  Bb.GetNumModes()),
59  2,Ba,Bb),
60  StdExpansion2D(LibUtilities::StdTriData::getNumberOfCoefficients(
61  Ba.GetNumModes(),
62  Bb.GetNumModes()),
63  Ba,Bb),
64  StdTriExp (Ba,Bb),
65  m_nodalPointsKey(Ba.GetNumModes(),Ntype)
66  {
67  ASSERTL0(m_base[0]->GetNumModes() == m_base[1]->GetNumModes(),
68  "Nodal basis initiated with different orders in the a "
69  "and b directions");
70  }
71 
73  StdExpansion(T),
74  StdExpansion2D(T),
75  StdTriExp(T),
76  m_nodalPointsKey(T.m_nodalPointsKey)
77  {
78  }
79 
81  {
82  }
83 
85  {
86  return true;
87  }
88 
89  //-------------------------------
90  // Nodal basis specific routines
91  //-------------------------------
92 
94  const Array<OneD, const NekDouble>& inarray,
95  Array<OneD, NekDouble>& outarray)
96  {
100  DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
101 
102  NekVector<NekDouble> nodal(m_ncoeffs,inarray,eWrapper);
103  NekVector<NekDouble> modal(m_ncoeffs,outarray,eWrapper);
104  modal = (*inv_vdm) * nodal;
105  }
106 
107  // Operate with transpose of NodalToModal transformation
109  const Array<OneD, const NekDouble>& inarray,
110  Array<OneD, NekDouble>& outarray)
111  {
115  DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
116 
117  NekVector<NekDouble> nodal(m_ncoeffs,inarray,eCopy);
118  NekVector<NekDouble> modal(m_ncoeffs,outarray,eWrapper);
119  modal = Transpose(*inv_vdm) * nodal;
120  }
121 
123  const Array<OneD, const NekDouble>& inarray,
124  Array<OneD, NekDouble>& outarray)
125  {
126  StdMatrixKey Nkey(eNBasisTrans, DetShapeType(), *this,
129  DNekMatSharedPtr vdm = GetStdMatrix(Nkey);
130 
131  // Multiply out matrix
132  NekVector<NekDouble> modal(m_ncoeffs,inarray,eWrapper);
133  NekVector<NekDouble> nodal(m_ncoeffs,outarray,eWrapper);
134  nodal = (*vdm)*modal;
135  }
136 
140  {
142  }
143 
145  {
146  int i,j;
149  DNekMatSharedPtr Mat;
150 
153  GetNodalPoints(r,s);
154 
155  //Store the values of m_phys in a temporary array
156  int nqtot = GetTotPoints();
157  Array<OneD,NekDouble> phys(nqtot);
158 
159  for(i = 0; i < m_ncoeffs; ++i)
160  {
161  // fill physical space with mode i
162  StdTriExp::v_FillMode(i,phys);
163 
164  // interpolate mode i to the Nodal points 'j' and
165  // store in outarray
166  for(j = 0; j < m_ncoeffs; ++j)
167  {
168  c[0] = r[j];
169  c[1] = s[j];
170  (*Mat)(j,i) = StdTriExp::v_PhysEvaluate(c,phys);
171  }
172  }
173  return Mat;
174  }
175 
176 
177  //---------------------------------------
178  // Transforms
179  //---------------------------------------
180 
182  const Array<OneD, const NekDouble>& inarray,
183  Array<OneD, NekDouble>& outarray)
184  {
185  v_BwdTrans_SumFac(inarray,outarray);
186  }
187 
189  const Array<OneD, const NekDouble>& inarray,
190  Array<OneD, NekDouble>& outarray)
191  {
193  NodalToModal(inarray,tmp);
194  StdTriExp::v_BwdTrans_SumFac(tmp,outarray);
195  }
196 
198  const Array<OneD, const NekDouble>& inarray,
199  Array<OneD, NekDouble>& outarray)
200  {
201  v_IProductWRTBase(inarray,outarray);
202 
203  // get Mass matrix inverse
204  StdMatrixKey masskey(eInvMass, DetShapeType(), *this,
207  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
208 
209  // copy inarray in case inarray == outarray
210  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
212 
213  out = (*matsys)*in;
214  }
215 
216 
217  //---------------------------------------
218  // Inner product functions
219  //---------------------------------------
220 
222  const Array<OneD, const NekDouble>& inarray,
223  Array<OneD, NekDouble>& outarray)
224  {
225  v_IProductWRTBase_SumFac(inarray,outarray);
226  }
227 
229  const Array<OneD, const NekDouble>& inarray,
230  Array<OneD, NekDouble>& outarray,
231  bool multiplybyweights)
232  {
233  StdTriExp::v_IProductWRTBase_SumFac(inarray,outarray,multiplybyweights);
234  NodalToModalTranspose(outarray,outarray);
235  }
236 
238  const int dir,
239  const Array<OneD, const NekDouble>& inarray,
240  Array<OneD, NekDouble>& outarray)
241  {
242  v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
243  }
244 
246  const int dir,
247  const Array<OneD, const NekDouble>& inarray,
248  Array<OneD, NekDouble>& outarray)
249  {
250  StdTriExp::v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
251  NodalToModalTranspose(outarray,outarray);
252  }
253 
254  //---------------------------------------
255  // Evaluation functions
256  //---------------------------------------
257 
259  const int mode,
260  Array<OneD, NekDouble> &outarray)
261  {
262  ASSERTL2(mode >= m_ncoeffs,
263  "calling argument mode is larger than total expansion order");
264 
265  Vmath::Zero(m_ncoeffs, outarray, 1);
266  outarray[mode] = 1.0;
267  v_BwdTrans(outarray,outarray);
268  }
269 
270 
271  //---------------------------
272  // Helper functions
273  //---------------------------
274 
276  {
277  return 3 + (GetBasisNumModes(0)-2) + 2*(GetBasisNumModes(1)-2);
278  }
279 
280  //--------------------------
281  // Mappings
282  //--------------------------
283 
284 
285  int StdNodalTriExp::v_GetVertexMap(const int localVertexId,
286  bool useCoeffPacking)
287  {
288  boost::ignore_unused(useCoeffPacking);
289  ASSERTL0(localVertexId >= 0 && localVertexId <= 2,
290  "Local Vertex ID must be between 0 and 2");
291  return localVertexId;
292  }
293 
295  const int eid,
296  Array<OneD, unsigned int> &maparray,
297  Array<OneD, int> &signarray,
298  Orientation edgeOrient,
299  int P, int Q)
300  {
301  boost::ignore_unused(Q);
302 
303  ASSERTL0(eid >= 0 && eid <= 2,
304  "Local Edge ID must be between 0 and 2");
305 
306  const int nEdgeCoeffs = GetTraceNcoeffs(eid);
307 
308  ASSERTL0(P == -1 || P == nEdgeCoeffs,
309  "Nodal triangle not set up to deal with variable "
310  "polynomial order.");
311 
312  if (maparray.size() != nEdgeCoeffs)
313  {
314  maparray = Array<OneD, unsigned int>(nEdgeCoeffs);
315  }
316 
317  if (signarray.size() != nEdgeCoeffs)
318  {
319  signarray = Array<OneD, int>(nEdgeCoeffs,1);
320  }
321  else
322  {
323  fill(signarray.get(), signarray.get()+nEdgeCoeffs, 1);
324  }
325 
326  Orientation orient = edgeOrient;
327  if (eid == 2)
328  {
329  orient = orient == eForwards ? eBackwards : eForwards;
330  }
331 
332  maparray[0] = eid;
333  maparray[nEdgeCoeffs-1] = eid == 2 ? 0 : eid+1;
334  for (int i = 2; i < nEdgeCoeffs; i++)
335  {
336  maparray[i-1] = eid*(nEdgeCoeffs-2)+1+i;
337  }
338 
339  if (orient == eBackwards)
340  {
341  reverse(maparray.get(), maparray.get()+nEdgeCoeffs);
342  }
343  }
344 
346  const int eid,
347  Array<OneD, unsigned int> &maparray,
348  Array<OneD, int> &signarray,
349  const Orientation edgeOrient)
350  {
351  ASSERTL0(eid >= 0 && eid <= 2,
352  "Local Edge ID must be between 0 and 2");
353 
354  const int nEdgeIntCoeffs = GetTraceNcoeffs(eid)-2;
355 
356  if (maparray.size() != nEdgeIntCoeffs)
357  {
358  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
359  }
360 
361  if (signarray.size() != nEdgeIntCoeffs)
362  {
363  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
364  }
365  else
366  {
367  fill(signarray.get(), signarray.get()+nEdgeIntCoeffs, 1);
368  }
369 
370  Orientation orient = edgeOrient;
371  if (eid == 2)
372  {
373  orient = orient == eForwards ? eBackwards : eForwards;
374  }
375 
376  for (int i = 0; i < nEdgeIntCoeffs; i++)
377  {
378  maparray[i] = eid*nEdgeIntCoeffs+3+i;
379  }
380 
381  if (orient == eBackwards)
382  {
383  reverse(maparray.get(), maparray.get()+nEdgeIntCoeffs);
384  }
385  }
386 
388  Array<OneD, unsigned int>& outarray)
389  {
390  unsigned int i;
391  if (outarray.size() != GetNcoeffs()-NumBndryCoeffs())
392  {
393  outarray = Array<OneD, unsigned int>(
395  }
396 
397  for (i = NumBndryCoeffs(); i < GetNcoeffs(); i++)
398  {
399  outarray[i-NumBndryCoeffs()] = i;
400  }
401  }
402 
404  Array<OneD, unsigned int>& outarray)
405  {
406  unsigned int i;
407  if (outarray.size()!=NumBndryCoeffs())
408  {
410  }
411 
412  for (i = 0; i < NumBndryCoeffs(); i++)
413  {
414  outarray[i] = i;
415  }
416  }
417 
418 
419  //---------------------------------------
420  // Wrapper functions
421  //---------------------------------------
422 
424  {
425  DNekMatSharedPtr Mat;
426 
427  switch(mkey.GetMatrixType())
428  {
429  case eNBasisTrans:
430  Mat = GenNBasisTransMatrix();
431  break;
432  default:
434  break;
435  }
436 
437  return Mat;
438  }
439 
441  const StdMatrixKey &mkey)
442  {
443  return StdNodalTriExp::v_GenMatrix(mkey);
444  }
445 
446 
447  //---------------------------------------
448  // Operator evaluation functions
449  //---------------------------------------
450 
452  const Array<OneD, const NekDouble> &inarray,
453  Array<OneD, NekDouble> &outarray,
454  const StdMatrixKey &mkey)
455  {
456  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
457  }
458 
460  const Array<OneD, const NekDouble> &inarray,
461  Array<OneD, NekDouble> &outarray,
462  const StdMatrixKey &mkey)
463  {
465  inarray,outarray,mkey);
466  }
467 
469  const int k1,
470  const int k2,
471  const Array<OneD, const NekDouble> &inarray,
472  Array<OneD, NekDouble> &outarray,
473  const StdMatrixKey &mkey)
474 
475  {
477  k1,k2,inarray,outarray,mkey);
478  }
479 
481  const int i,
482  const Array<OneD, const NekDouble> &inarray,
483  Array<OneD, NekDouble> &outarray,
484  const StdMatrixKey &mkey)
485  {
486  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
487  }
488 
490  const Array<OneD, const NekDouble> &inarray,
491  Array<OneD, NekDouble> &outarray,
492  const StdMatrixKey &mkey)
493  {
495  inarray,outarray,mkey);
496  }
497 
498 
499  //---------------------------------------
500  // Private helper functions
501  //---------------------------------------
502 
503  } // end of namespace
504 } // end of namespace
505 
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:274
Describes the specification for a Basis.
Definition: Basis.h:50
PointsType GetPointsType() const
Definition: Points.h:112
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)
This function evaluates the expansion at a single (arbitrary) point of the domain.
The base class for all shapes.
Definition: StdExpansion.h:63
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
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:617
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:376
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:265
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:171
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:81
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void NodalToModal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
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)
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Transform a given function from physical quadrature space to coefficient space.
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::PointsKey m_nodalPointsKey
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
virtual void v_GetTraceToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient=eForwards, int P=-1, int Q=-1)
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculate the inner product of inarray with respect to the basis B=base0[p]*base1[pq] and put into ou...
virtual void v_GetTraceInteriorToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation edgeOrient=eForwards)
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Backward tranform for triangular elements.
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:259
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdTriExp.cpp:485
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:703
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:597
PointsManagerT & PointsManager(void)
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:315
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:273
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
P
Definition: main.py:133