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),
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 
285  const int eid,
286  const Orientation edgeOrient,
287  Array<OneD, unsigned int> &maparray,
288  Array<OneD, int> &signarray,
289  int P)
290  {
291  ASSERTL0(eid >= 0 && eid <= 2,
292  "Local Edge ID must be between 0 and 2");
293 
294  ASSERTL0(P == -1, "Nodal triangle not set up to deal with variable"
295  "polynomial order.");
296  const int nEdgeCoeffs = GetEdgeNcoeffs(eid);
297 
298  if (maparray.num_elements() != nEdgeCoeffs)
299  {
300  maparray = Array<OneD, unsigned int>(nEdgeCoeffs);
301  }
302 
303  if (signarray.num_elements() != nEdgeCoeffs)
304  {
305  signarray = Array<OneD, int>(nEdgeCoeffs,1);
306  }
307  else
308  {
309  fill(signarray.get(), signarray.get()+nEdgeCoeffs, 1);
310  }
311 
312  Orientation orient = edgeOrient;
313  if (eid == 2)
314  {
315  orient = orient == eForwards ? eBackwards : eForwards;
316  }
317 
318  maparray[0] = eid;
319  maparray[1] = eid == 2 ? 0 : eid+1;
320  for (int i = 2; i < nEdgeCoeffs; i++)
321  {
322  maparray[i] = eid*(nEdgeCoeffs-2)+1+i;
323  }
324 
325  if (orient == eBackwards)
326  {
327  reverse(maparray.get(), maparray.get()+nEdgeCoeffs);
328  }
329  }
330 
331  int StdNodalTriExp::v_GetVertexMap(const int localVertexId,
332  bool useCoeffPacking)
333  {
334  boost::ignore_unused(useCoeffPacking);
335  ASSERTL0(localVertexId >= 0 && localVertexId <= 2,
336  "Local Vertex ID must be between 0 and 2");
337  return localVertexId;
338  }
339 
341  const int eid,
342  const Orientation edgeOrient,
343  Array<OneD, unsigned int> &maparray,
344  Array<OneD, int> &signarray)
345  {
346  ASSERTL0(eid >= 0 && eid <= 2,
347  "Local Edge ID must be between 0 and 2");
348 
349  const int nEdgeIntCoeffs = GetEdgeNcoeffs(eid)-2;
350 
351  if (maparray.num_elements() != nEdgeIntCoeffs)
352  {
353  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
354  }
355 
356  if (signarray.num_elements() != nEdgeIntCoeffs)
357  {
358  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
359  }
360  else
361  {
362  fill(signarray.get(), signarray.get()+nEdgeIntCoeffs, 1);
363  }
364 
365  Orientation orient = edgeOrient;
366  if (eid == 2)
367  {
368  orient = orient == eForwards ? eBackwards : eForwards;
369  }
370 
371  for (int i = 0; i < nEdgeIntCoeffs; i++)
372  {
373  maparray[i] = eid*nEdgeIntCoeffs+3+i;
374  }
375 
376  if (orient == eBackwards)
377  {
378  reverse(maparray.get(), maparray.get()+nEdgeIntCoeffs);
379  }
380  }
381 
383  Array<OneD, unsigned int>& outarray)
384  {
385  unsigned int i;
386  if (outarray.num_elements() != GetNcoeffs()-NumBndryCoeffs())
387  {
388  outarray = Array<OneD, unsigned int>(
390  }
391 
392  for (i = NumBndryCoeffs(); i < GetNcoeffs(); i++)
393  {
394  outarray[i-NumBndryCoeffs()] = i;
395  }
396  }
397 
399  Array<OneD, unsigned int>& outarray)
400  {
401  unsigned int i;
402  if (outarray.num_elements()!=NumBndryCoeffs())
403  {
405  }
406 
407  for (i = 0; i < NumBndryCoeffs(); i++)
408  {
409  outarray[i] = i;
410  }
411  }
412 
413 
414  //---------------------------------------
415  // Wrapper functions
416  //---------------------------------------
417 
419  {
420  DNekMatSharedPtr Mat;
421 
422  switch(mkey.GetMatrixType())
423  {
424  case eNBasisTrans:
425  Mat = GenNBasisTransMatrix();
426  break;
427  default:
429  break;
430  }
431 
432  return Mat;
433  }
434 
436  const StdMatrixKey &mkey)
437  {
438  return StdNodalTriExp::v_GenMatrix(mkey);
439  }
440 
441 
442  //---------------------------------------
443  // Operator evaluation functions
444  //---------------------------------------
445 
447  const Array<OneD, const NekDouble> &inarray,
448  Array<OneD, NekDouble> &outarray,
449  const StdMatrixKey &mkey)
450  {
451  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
452  }
453 
455  const Array<OneD, const NekDouble> &inarray,
456  Array<OneD, NekDouble> &outarray,
457  const StdMatrixKey &mkey)
458  {
460  inarray,outarray,mkey);
461  }
462 
464  const int k1,
465  const int k2,
466  const Array<OneD, const NekDouble> &inarray,
467  Array<OneD, NekDouble> &outarray,
468  const StdMatrixKey &mkey)
469 
470  {
472  k1,k2,inarray,outarray,mkey);
473  }
474 
476  const int i,
477  const Array<OneD, const NekDouble> &inarray,
478  Array<OneD, NekDouble> &outarray,
479  const StdMatrixKey &mkey)
480  {
481  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
482  }
483 
485  const Array<OneD, const NekDouble> &inarray,
486  Array<OneD, NekDouble> &outarray,
487  const StdMatrixKey &mkey)
488  {
490  inarray,outarray,mkey);
491  }
492 
493 
494  //---------------------------------------
495  // Private helper functions
496  //---------------------------------------
497 
498  } // end of namespace
499 } // end of namespace
500 
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
PointsType GetPointsType() const
Definition: Points.h:112
void ModalToNodal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &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...
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:469
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:177
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
STL namespace.
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:598
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:81
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:714
virtual void v_WeakDerivMatrixOp(const int i, 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)
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Backward tranform for triangular elements.
LibUtilities::PointsKey m_nodalPointsKey
The base class for all shapes.
Definition: StdExpansion.h:68
int GetEdgeNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge.
Definition: StdExpansion.h:286
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
virtual void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P)
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_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void GetNodalPoints(Array< OneD, const NekDouble > &x, Array< OneD, const NekDouble > &y)
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.
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
PointsManagerT & PointsManager(void)
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdTriExp.cpp:486
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:259
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:695
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void NodalToModal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
Array< OneD, LibUtilities::BasisSharedPtr > m_base
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
virtual void v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
Describes the specification for a Basis.
Definition: Basis.h:49
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:265
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)
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:295