Nektar++
StdNodalTetExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File StdNodalTetExp.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 tetrahedral routines built upon StdExpansion3D
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
38 #include <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
39 
40 namespace Nektar
41 {
42  namespace StdRegions
43  {
45  const LibUtilities::BasisKey &Ba,
46  const LibUtilities::BasisKey &Bb,
47  const LibUtilities::BasisKey &Bc,
49  StdExpansion (LibUtilities::StdTetData::getNumberOfCoefficients(
50  Ba.GetNumModes(),
51  Bb.GetNumModes(),
52  Bc.GetNumModes()),
53  3,Ba,Bb,Bc),
54  StdExpansion3D(LibUtilities::StdTetData::getNumberOfCoefficients(
55  Ba.GetNumModes(),
56  Bb.GetNumModes(),
57  Bc.GetNumModes()),
58  Ba,Bb,Bc),
59  StdTetExp (Ba,Bb,Bc),
60  m_nodalPointsKey(Ba.GetNumModes(),Ntype)
61  {
62  ASSERTL0(Ba.GetNumModes() <= Bb.GetNumModes(),
63  "order in 'a' direction is higher than order "
64  "in 'b' direction");
65  ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
66  "order in 'a' direction is higher than order "
67  "in 'c' direction");
68  ASSERTL0(Bb.GetNumModes() <= Bc.GetNumModes(),
69  "order in 'b' direction is higher than order "
70  "in 'c' direction");
71  }
72 
74  StdExpansion(T),
75  StdExpansion3D(T),
76  StdTetExp(T),
78  {
79  }
80 
82  {
83  }
84 
86  {
87  return true;
88  }
89 
90  //-------------------------------
91  // Nodal basis specific routines
92  //-------------------------------
93 
95  const Array<OneD, const NekDouble>& inarray,
96  Array<OneD, NekDouble>& outarray)
97  {
101  DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
102 
103  NekVector<NekDouble> nodal(m_ncoeffs,inarray,eWrapper);
104  NekVector<NekDouble> modal(m_ncoeffs,outarray,eWrapper);
105  modal = (*inv_vdm) * nodal;
106  }
107 
108  // Operate with transpose of NodalToModal transformation
110  const Array<OneD, const NekDouble>& inarray,
111  Array<OneD, NekDouble>& outarray)
112  {
116  DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
117 
118  NekVector<NekDouble> nodal(m_ncoeffs,inarray,eCopy);
119  NekVector<NekDouble> modal(m_ncoeffs,outarray,eWrapper);
120  modal = Transpose(*inv_vdm) * nodal;
121  }
122 
124  const Array<OneD, const NekDouble>& inarray,
125  Array<OneD, NekDouble>& outarray)
126  {
127  StdMatrixKey Nkey(eNBasisTrans, DetShapeType(), *this,
130  DNekMatSharedPtr vdm = GetStdMatrix(Nkey);
131 
132  // Multiply out matrix
133  NekVector<NekDouble> modal(m_ncoeffs,inarray,eWrapper);
134  NekVector<NekDouble> nodal(m_ncoeffs,outarray,eWrapper);
135  nodal = (*vdm)*modal;
136  }
137 
142  {
143  LibUtilities::PointsManager()[m_nodalPointsKey]->GetPoints(x,y,z);
144  }
145 
147  {
148  int i,j;
151  DNekMatSharedPtr Mat;
152 
155  GetNodalPoints(r,s,t);
156 
157  //Store the values of m_phys in a temporary array
158  int nqtot = GetTotPoints();
159  Array<OneD,NekDouble> tmp_phys(nqtot);
160 
161  for(i = 0; i < m_ncoeffs; ++i)
162  {
163  // fill physical space with mode i
164  StdTetExp::v_FillMode(i,tmp_phys);
165 
166  // interpolate mode i to the Nodal points 'j' and
167  // store in outarray
168  for(j = 0; j < m_ncoeffs; ++j)
169  {
170  c[0] = r[j];
171  c[1] = s[j];
172  c[2] = t[j];
173  (*Mat)(j,i) = StdTetExp::v_PhysEvaluate(c,tmp_phys);
174  }
175  }
176 
177  return Mat;
178  }
179 
180 
181  //---------------------------------------
182  // Transforms
183  //---------------------------------------
184 
186  const Array<OneD, const NekDouble>& inarray,
187  Array<OneD, NekDouble>& outarray)
188  {
189  v_BwdTrans_SumFac(inarray,outarray);
190  }
191 
193  const Array<OneD, const NekDouble>& inarray,
194  Array<OneD, NekDouble>& outarray)
195  {
197  NodalToModal(inarray,tmp);
198  StdTetExp::v_BwdTrans_SumFac(tmp,outarray);
199  }
200 
202  const Array<OneD, const NekDouble>& inarray,
203  Array<OneD, NekDouble>& outarray)
204  {
205  v_IProductWRTBase(inarray,outarray);
206 
207  // get Mass matrix inverse
208  StdMatrixKey masskey(eInvMass, DetShapeType(), *this,
211  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
212 
213  // copy inarray in case inarray == outarray
214  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
216 
217  out = (*matsys)*in;
218  }
219 
220 
221  //---------------------------------------
222  // Inner product functions
223  //---------------------------------------
224 
226  const Array<OneD, const NekDouble>& inarray,
227  Array<OneD, NekDouble>& outarray)
228  {
229  v_IProductWRTBase_SumFac(inarray,outarray);
230  }
231 
233  const Array<OneD, const NekDouble>& inarray,
234  Array<OneD, NekDouble>& outarray,
235  bool multiplybyweights)
236  {
237  StdTetExp::v_IProductWRTBase_SumFac(inarray,outarray,multiplybyweights);
238  NodalToModalTranspose(outarray,outarray);
239  }
240 
242  const int dir,
243  const Array<OneD, const NekDouble>& inarray,
244  Array<OneD, NekDouble>& outarray)
245  {
246  v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
247  }
248 
250  const int dir,
251  const Array<OneD, const NekDouble>& inarray,
252  Array<OneD, NekDouble>& outarray)
253  {
254  StdTetExp::v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
255  NodalToModalTranspose(outarray,outarray);
256  }
257 
258  //---------------------------------------
259  // Evaluation functions
260  //---------------------------------------
261 
263  const int mode,
264  Array<OneD, NekDouble> &outarray)
265  {
266  ASSERTL2(mode >= m_ncoeffs,
267  "calling argument mode is larger than total expansion order");
268 
269  Vmath::Zero(m_ncoeffs, outarray, 1);
270  outarray[mode] = 1.0;
271  v_BwdTrans(outarray,outarray);
272  }
273 
274 
275  //---------------------------------------
276  // Mapping functions
277  //---------------------------------------
278 
279  /*
280  void StdNodalTriExp::v_GetFaceToElementMap(
281  const int fid,
282  const FaceOrientation faceOrient,
283  Array<OneD, unsigned int> &maparray,
284  Array<OneD, int> &signarray,
285  int nummodesA,
286  int nummodesB)
287  {
288  int P, Q, i, j, k, idx = 0, nFaceCoeffs = 0;
289 
290  ASSERTL0(fid >= 0 && fid <= 3,
291  "Local face ID must be between 0 and 3");
292 
293  if (nummodesA == -1)
294  {
295  switch(fid)
296  {
297  case 0:
298  nummodesA = m_base[0]->GetNumModes();
299  nummodesB = m_base[1]->GetNumModes();
300  break;
301  case 1:
302  nummodesA = m_base[0]->GetNumModes();
303  nummodesB = m_base[2]->GetNumModes();
304  break;
305  case 2:
306  case 3:
307  nummodesA = m_base[1]->GetNumModes();
308  nummodesB = m_base[2]->GetNumModes();
309  break;
310  }
311  }
312 
313  P = nummodesA;
314  Q = nummodesB;
315  nFaceCoeffs = Q + ((P-1)*(1 + 2*(Q-1) - (P-1)))/2;
316 
317  if (maparray.num_elements() != nFaceCoeffs)
318  {
319  maparray = Array<OneD, unsigned int>(nFaceCoeffs);
320  }
321 
322  if (signarray.num_elements() != nFaceCoeffs)
323  {
324  signarray = Array<OneD, int>(nFaceCoeffs,1);
325  }
326  else
327  {
328  fill(signarray.get(), signarray.get()+nFaceCoeffs, 1);
329  }
330 
331  switch(fid)
332  {
333  case 0:
334  // Add vertices.
335  maparray[idx++] = 0;
336  maparray[idx++] = 1;
337  maparray[idx++] = 2;
338 
339  // Add edges.
340  for (i = 2; i < P; ++i)
341  {
342  maparray[idx++] = ;
343  }
344  }
345  }
346  */
347 
348  int StdNodalTetExp::v_GetVertexMap(const int localVertexId,
349  bool useCoeffPacking)
350  {
351  boost::ignore_unused(useCoeffPacking);
352  ASSERTL0(localVertexId >= 0 && localVertexId <= 3,
353  "Local Vertex ID must be between 0 and 3");
354  return localVertexId;
355  }
356 
358  Array<OneD, unsigned int>& outarray)
359  {
360  unsigned int i;
361  const unsigned int nBndryCoeff = NumBndryCoeffs();
362 
363  if (outarray.num_elements() != nBndryCoeff)
364  {
365  outarray = Array<OneD, unsigned int>(nBndryCoeff);
366  }
367 
368  for (i = 0; i < nBndryCoeff; i++)
369  {
370  outarray[i] = i;
371  }
372  }
373 
375  Array<OneD, unsigned int>& outarray)
376  {
377  unsigned int i;
378  const unsigned int nBndryCoeff = NumBndryCoeffs();
379 
380  if (outarray.num_elements() != m_ncoeffs-nBndryCoeff)
381  {
382  outarray = Array<OneD, unsigned int>(
383  m_ncoeffs-nBndryCoeff);
384  }
385 
386  for (i = nBndryCoeff; i < m_ncoeffs; i++)
387  {
388  outarray[i-nBndryCoeff] = i;
389  }
390  }
391 
392 
393  //---------------------------------------
394  // Wrapper functions
395  //---------------------------------------
396 
398  {
399  DNekMatSharedPtr Mat;
400 
401  switch(mkey.GetMatrixType())
402  {
403  case eNBasisTrans:
404  Mat = GenNBasisTransMatrix();
405  break;
406  default:
408  break;
409  }
410 
411  return Mat;
412  }
413 
415  const StdMatrixKey &mkey)
416  {
417  return StdNodalTetExp::v_GenMatrix(mkey);
418  }
419  } // end of namespace
420 } // end of namespace
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
PointsType GetPointsType() const
Definition: Points.h:112
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
void GetNodalPoints(Array< OneD, const NekDouble > &x, Array< OneD, const NekDouble > &y, Array< OneD, const NekDouble > &z)
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
LibUtilities::PointsKey m_nodalPointsKey
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdTetExp.cpp:549
LibUtilities::ShapeType DetShapeType() const
Definition: StdTetExp.h:69
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:81
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:714
virtual int v_GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
Definition: StdTetExp.cpp:905
The base class for all shapes.
Definition: StdExpansion.h:68
void NodalToModal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:86
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTetExp.cpp:324
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTetExp.cpp:720
virtual void v_IProductWRTBase(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)
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 DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274
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.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool mult=true)
void ModalToNodal(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
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
void NodalToModalTranspose(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Describes the specification for a Basis.
Definition: Basis.h:49
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:265
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:295