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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Nodal triangle routines built upon StdExpansion2D
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 #include <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
38 
39 namespace Nektar
40 {
41  namespace StdRegions
42  {
44  StdTriExp(),
45  m_nodalPointsKey()
46  {
47  }
48 
50  const LibUtilities::BasisKey &Ba,
51  const LibUtilities::BasisKey &Bb,
53  StdExpansion (LibUtilities::StdTriData::getNumberOfCoefficients(
54  Ba.GetNumModes(),
55  Bb.GetNumModes()),
56  2,Ba,Bb),
57  StdExpansion2D(LibUtilities::StdTriData::getNumberOfCoefficients(
58  Ba.GetNumModes(),
59  Bb.GetNumModes()),
60  Ba,Bb),
61  StdTriExp (Ba,Bb),
62  m_nodalPointsKey(Ba.GetNumModes(),Ntype)
63  {
64  ASSERTL0(m_base[0]->GetNumModes() == m_base[1]->GetNumModes(),
65  "Nodal basis initiated with different orders in the a "
66  "and b directions");
67  }
68 
70  StdExpansion(T),
71  StdExpansion2D(T),
72  StdTriExp(T),
73  m_nodalPointsKey(T.m_nodalPointsKey)
74  {
75  }
76 
78  {
79  }
80 
82  {
83  return true;
84  }
85 
86  //-------------------------------
87  // Nodal basis specific routines
88  //-------------------------------
89 
91  const Array<OneD, const NekDouble>& inarray,
92  Array<OneD, NekDouble>& outarray)
93  {
97  DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
98 
100  NekVector<NekDouble> modal(m_ncoeffs,outarray,eWrapper);
101  modal = (*inv_vdm) * nodal;
102  }
103 
104  // Operate with transpose of NodalToModal transformation
106  const Array<OneD, const NekDouble>& inarray,
107  Array<OneD, NekDouble>& outarray)
108  {
112  DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
113 
114  NekVector<NekDouble> nodal(m_ncoeffs,inarray,eCopy);
115  NekVector<NekDouble> modal(m_ncoeffs,outarray,eWrapper);
116  modal = Transpose(*inv_vdm) * nodal;
117  }
118 
120  const Array<OneD, const NekDouble>& inarray,
121  Array<OneD, NekDouble>& outarray)
122  {
123  StdMatrixKey Nkey(eNBasisTrans, DetShapeType(), *this,
126  DNekMatSharedPtr vdm = GetStdMatrix(Nkey);
127 
128  // Multiply out matrix
129  NekVector<NekDouble> modal(m_ncoeffs,inarray,eWrapper);
130  NekVector<NekDouble> nodal(m_ncoeffs,outarray,eWrapper);
131  nodal = (*vdm)*modal;
132  }
133 
137  {
139  }
140 
142  {
143  int i,j;
146  DNekMatSharedPtr Mat;
147 
150  GetNodalPoints(r,s);
151 
152  //Store the values of m_phys in a temporary array
153  int nqtot = GetTotPoints();
154  Array<OneD,NekDouble> phys(nqtot);
155 
156  for(i = 0; i < m_ncoeffs; ++i)
157  {
158  // fill physical space with mode i
159  StdTriExp::v_FillMode(i,phys);
160 
161  // interpolate mode i to the Nodal points 'j' and
162  // store in outarray
163  for(j = 0; j < m_ncoeffs; ++j)
164  {
165  c[0] = r[j];
166  c[1] = s[j];
167  (*Mat)(j,i) = StdTriExp::v_PhysEvaluate(c,phys);
168  }
169  }
170  return Mat;
171  }
172 
173 
174  //---------------------------------------
175  // Transforms
176  //---------------------------------------
177 
179  const Array<OneD, const NekDouble>& inarray,
180  Array<OneD, NekDouble>& outarray)
181  {
182  v_BwdTrans_SumFac(inarray,outarray);
183  }
184 
186  const Array<OneD, const NekDouble>& inarray,
187  Array<OneD, NekDouble>& outarray)
188  {
190  NodalToModal(inarray,tmp);
191  StdTriExp::v_BwdTrans_SumFac(tmp,outarray);
192  }
193 
195  const Array<OneD, const NekDouble>& inarray,
196  Array<OneD, NekDouble>& outarray)
197  {
198  v_IProductWRTBase(inarray,outarray);
199 
200  // get Mass matrix inverse
201  StdMatrixKey masskey(eInvMass, DetShapeType(), *this,
204  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
205 
206  // copy inarray in case inarray == outarray
207  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
209 
210  out = (*matsys)*in;
211  }
212 
213 
214  //---------------------------------------
215  // Inner product functions
216  //---------------------------------------
217 
219  const Array<OneD, const NekDouble>& inarray,
220  Array<OneD, NekDouble>& outarray)
221  {
222  v_IProductWRTBase_SumFac(inarray,outarray);
223  }
224 
226  const Array<OneD, const NekDouble>& inarray,
227  Array<OneD, NekDouble>& outarray,
228  bool multiplybyweights)
229  {
230  StdTriExp::v_IProductWRTBase_SumFac(inarray,outarray,multiplybyweights);
231  NodalToModalTranspose(outarray,outarray);
232  }
233 
235  const int dir,
236  const Array<OneD, const NekDouble>& inarray,
237  Array<OneD, NekDouble>& outarray)
238  {
239  v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
240  }
241 
243  const int dir,
244  const Array<OneD, const NekDouble>& inarray,
245  Array<OneD, NekDouble>& outarray)
246  {
247  StdTriExp::v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
248  NodalToModalTranspose(outarray,outarray);
249  }
250 
251  //---------------------------------------
252  // Evaluation functions
253  //---------------------------------------
254 
256  const int mode,
257  Array<OneD, NekDouble> &outarray)
258  {
259  ASSERTL2(mode >= m_ncoeffs,
260  "calling argument mode is larger than total expansion order");
261 
262  Vmath::Zero(m_ncoeffs, outarray, 1);
263  outarray[mode] = 1.0;
264  v_BwdTrans(outarray,outarray);
265  }
266 
267 
268  //---------------------------
269  // Helper functions
270  //---------------------------
271 
273  {
274  return 3 + (GetBasisNumModes(0)-2) + 2*(GetBasisNumModes(1)-2);
275  }
276 
277  //--------------------------
278  // Mappings
279  //--------------------------
280 
282  const int eid,
283  const Orientation edgeOrient,
284  Array<OneD, unsigned int> &maparray,
285  Array<OneD, int> &signarray)
286  {
287  ASSERTL0(eid >= 0 && eid <= 2,
288  "Local Edge ID must be between 0 and 2");
289 
290  const int nEdgeCoeffs = GetEdgeNcoeffs(eid);
291 
292  if (maparray.num_elements() != nEdgeCoeffs)
293  {
294  maparray = Array<OneD, unsigned int>(nEdgeCoeffs);
295  }
296 
297  if (signarray.num_elements() != nEdgeCoeffs)
298  {
299  signarray = Array<OneD, int>(nEdgeCoeffs,1);
300  }
301  else
302  {
303  fill(signarray.get(), signarray.get()+nEdgeCoeffs, 1);
304  }
305 
306  maparray[0] = eid;
307  maparray[1] = eid == 2 ? 0 : eid+1;
308  for (int i = 2; i < nEdgeCoeffs; i++)
309  {
310  maparray[i] = eid*(nEdgeCoeffs-2)+1+i;
311  }
312 
313  if (edgeOrient == eBackwards)
314  {
315  reverse(maparray.get(), maparray.get()+nEdgeCoeffs);
316  }
317  }
318 
319  int StdNodalTriExp::v_GetVertexMap(const int localVertexId,
320  bool useCoeffPacking)
321  {
322  ASSERTL0(localVertexId >= 0 && localVertexId <= 2,
323  "Local Vertex ID must be between 0 and 2");
324  return localVertexId;
325  }
326 
328  const int eid,
329  const Orientation edgeOrient,
330  Array<OneD, unsigned int> &maparray,
331  Array<OneD, int> &signarray)
332  {
333  ASSERTL0(eid >= 0 && eid <= 2,
334  "Local Edge ID must be between 0 and 2");
335 
336  const int nEdgeIntCoeffs = GetEdgeNcoeffs(eid)-2;
337 
338  if (maparray.num_elements() != nEdgeIntCoeffs)
339  {
340  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
341  }
342 
343  if (signarray.num_elements() != nEdgeIntCoeffs)
344  {
345  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
346  }
347  else
348  {
349  fill(signarray.get(), signarray.get()+nEdgeIntCoeffs, 1);
350  }
351 
352  for (int i = 0; i < nEdgeIntCoeffs; i++)
353  {
354  maparray[i] = eid*nEdgeIntCoeffs+3+i;
355  }
356 
357  if (edgeOrient == eBackwards)
358  {
359  reverse(maparray.get(), maparray.get()+nEdgeIntCoeffs);
360  }
361  }
362 
364  Array<OneD, unsigned int>& outarray)
365  {
366  unsigned int i;
367  if (outarray.num_elements() != GetNcoeffs()-NumBndryCoeffs())
368  {
369  outarray = Array<OneD, unsigned int>(
371  }
372 
373  for (i = NumBndryCoeffs(); i < GetNcoeffs(); i++)
374  {
375  outarray[i-NumBndryCoeffs()] = i;
376  }
377  }
378 
380  Array<OneD, unsigned int>& outarray)
381  {
382  unsigned int i;
383  if (outarray.num_elements()!=NumBndryCoeffs())
384  {
386  }
387 
388  for (i = 0; i < NumBndryCoeffs(); i++)
389  {
390  outarray[i] = i;
391  }
392  }
393 
394 
395  //---------------------------------------
396  // Wrapper functions
397  //---------------------------------------
398 
400  {
401  DNekMatSharedPtr Mat;
402 
403  switch(mkey.GetMatrixType())
404  {
405  case eNBasisTrans:
406  Mat = GenNBasisTransMatrix();
407  break;
408  default:
410  break;
411  }
412 
413  return Mat;
414  }
415 
417  const StdMatrixKey &mkey)
418  {
419  return StdNodalTriExp::v_GenMatrix(mkey);
420  }
421 
422 
423  //---------------------------------------
424  // Operator evaluation functions
425  //---------------------------------------
426 
428  const Array<OneD, const NekDouble> &inarray,
429  Array<OneD, NekDouble> &outarray,
430  const StdMatrixKey &mkey)
431  {
432  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
433  }
434 
436  const Array<OneD, const NekDouble> &inarray,
437  Array<OneD, NekDouble> &outarray,
438  const StdMatrixKey &mkey)
439  {
441  inarray,outarray,mkey);
442  }
443 
445  const int k1,
446  const int k2,
447  const Array<OneD, const NekDouble> &inarray,
448  Array<OneD, NekDouble> &outarray,
449  const StdMatrixKey &mkey)
450 
451  {
453  k1,k2,inarray,outarray,mkey);
454  }
455 
457  const int i,
458  const Array<OneD, const NekDouble> &inarray,
459  Array<OneD, NekDouble> &outarray,
460  const StdMatrixKey &mkey)
461  {
462  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
463  }
464 
466  const Array<OneD, const NekDouble> &inarray,
467  Array<OneD, NekDouble> &outarray,
468  const StdMatrixKey &mkey)
469  {
471  inarray,outarray,mkey);
472  }
473 
474 
475  //---------------------------------------
476  // Private helper functions
477  //---------------------------------------
478 
479  } // end of namespace
480 } // end of namespace
481 
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:82
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...
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:593
int GetEdgeNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge.
Definition: StdExpansion.h:287
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:684
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.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
LibUtilities::PointsKey m_nodalPointsKey
The base class for all shapes.
Definition: StdExpansion.h:69
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
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)
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:483
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:258
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:213
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:690
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:131
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:359
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:50
PointsType GetPointsType() const
Definition: Points.h:111
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:226
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:249