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 
37 #include <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
39 
40 namespace Nektar
41 {
42 namespace StdRegions
43 {
45 {
46 }
47 
49  const LibUtilities::BasisKey &Bb,
50  const LibUtilities::BasisKey &Bc,
52  : StdExpansion(LibUtilities::StdTetData::getNumberOfCoefficients(
53  Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
54  3, Ba, Bb, Bc),
55  StdExpansion3D(LibUtilities::StdTetData::getNumberOfCoefficients(
56  Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
57  Ba, Bb, Bc),
58  StdTetExp(Ba, Bb, Bc), m_nodalPointsKey(Ba.GetNumModes(), Ntype)
59 {
60  ASSERTL0(Ba.GetNumModes() <= Bb.GetNumModes(),
61  "order in 'a' direction is higher than order "
62  "in 'b' direction");
63  ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
64  "order in 'a' direction is higher than order "
65  "in 'c' direction");
66  ASSERTL0(Bb.GetNumModes() <= Bc.GetNumModes(),
67  "order in 'b' direction is higher than order "
68  "in 'c' direction");
69 }
70 
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  Array<OneD, NekDouble> &outarray)
92 {
96  DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
97 
98  NekVector<NekDouble> nodal(m_ncoeffs, inarray, eWrapper);
99  NekVector<NekDouble> modal(m_ncoeffs, outarray, eWrapper);
100  modal = (*inv_vdm) * nodal;
101 }
102 
103 // Operate with transpose of NodalToModal transformation
105  const Array<OneD, const NekDouble> &inarray,
106  Array<OneD, NekDouble> &outarray)
107 {
111  DNekMatSharedPtr inv_vdm = GetStdMatrix(Nkey);
112 
113  NekVector<NekDouble> nodal(m_ncoeffs, inarray, eCopy);
114  NekVector<NekDouble> modal(m_ncoeffs, outarray, eWrapper);
115  modal = Transpose(*inv_vdm) * nodal;
116 }
117 
119  Array<OneD, NekDouble> &outarray)
120 {
123  DNekMatSharedPtr vdm = GetStdMatrix(Nkey);
124 
125  // Multiply out matrix
126  NekVector<NekDouble> modal(m_ncoeffs, inarray, eWrapper);
127  NekVector<NekDouble> nodal(m_ncoeffs, outarray, eWrapper);
128  nodal = (*vdm) * modal;
129 }
130 
134 {
135  LibUtilities::PointsManager()[m_nodalPointsKey]->GetPoints(x, y, z);
136 }
137 
139 {
140  int i, j;
143  DNekMatSharedPtr Mat;
144 
146  GetNodalPoints(r, s, t);
147 
148  // Store the values of m_phys in a temporary array
149  int nqtot = GetTotPoints();
150  Array<OneD, NekDouble> tmp_phys(nqtot);
151 
152  for (i = 0; i < m_ncoeffs; ++i)
153  {
154  // fill physical space with mode i
155  StdTetExp::v_FillMode(i, tmp_phys);
156 
157  // interpolate mode i to the Nodal points 'j' and
158  // store in outarray
159  for (j = 0; j < m_ncoeffs; ++j)
160  {
161  c[0] = r[j];
162  c[1] = s[j];
163  c[2] = t[j];
164  (*Mat)(j, i) = StdExpansion3D::v_PhysEvaluate(c, tmp_phys);
165  }
166  }
167 
168  return Mat;
169 }
170 
171 //---------------------------------------
172 // Transforms
173 //---------------------------------------
174 
176  Array<OneD, NekDouble> &outarray)
177 {
178  v_BwdTrans_SumFac(inarray, outarray);
179 }
180 
182  const Array<OneD, const NekDouble> &inarray,
183  Array<OneD, NekDouble> &outarray)
184 {
186  NodalToModal(inarray, tmp);
187  StdTetExp::v_BwdTrans_SumFac(tmp, outarray);
188 }
189 
191  Array<OneD, NekDouble> &outarray)
192 {
193  v_IProductWRTBase(inarray, outarray);
194 
195  // get Mass matrix inverse
198  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
199 
200  // copy inarray in case inarray == outarray
201  NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
202  NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
203 
204  out = (*matsys) * in;
205 }
206 
207 //---------------------------------------
208 // Inner product functions
209 //---------------------------------------
210 
212  const Array<OneD, const NekDouble> &inarray,
213  Array<OneD, NekDouble> &outarray)
214 {
215  v_IProductWRTBase_SumFac(inarray, outarray);
216 }
217 
219  const Array<OneD, const NekDouble> &inarray,
220  Array<OneD, NekDouble> &outarray, bool multiplybyweights)
221 {
222  StdTetExp::v_IProductWRTBase_SumFac(inarray, outarray, multiplybyweights);
223  NodalToModalTranspose(outarray, outarray);
224 }
225 
227  const int dir, const Array<OneD, const NekDouble> &inarray,
228  Array<OneD, NekDouble> &outarray)
229 {
230  v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
231 }
232 
234  const int dir, const Array<OneD, const NekDouble> &inarray,
235  Array<OneD, NekDouble> &outarray)
236 {
237  StdTetExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
238  NodalToModalTranspose(outarray, outarray);
239 }
240 
241 //---------------------------------------
242 // Evaluation functions
243 //---------------------------------------
244 
245 void StdNodalTetExp::v_FillMode(const int mode,
246  Array<OneD, NekDouble> &outarray)
247 {
248  ASSERTL2(mode >= m_ncoeffs,
249  "calling argument mode is larger than total expansion order");
250 
251  Vmath::Zero(m_ncoeffs, outarray, 1);
252  outarray[mode] = 1.0;
253  v_BwdTrans(outarray, outarray);
254 }
255 
256 //---------------------------------------
257 // Mapping functions
258 //---------------------------------------
259 
260 /*
261 void StdNodalTriExp::v_GetFaceToElementMap(
262  const int fid,
263  const FaceOrientation faceOrient,
264  Array<OneD, unsigned int> &maparray,
265  Array<OneD, int> &signarray,
266  int nummodesA,
267  int nummodesB)
268 {
269  int P, Q, i, j, k, idx = 0, nFaceCoeffs = 0;
270 
271  ASSERTL0(fid >= 0 && fid <= 3,
272  "Local face ID must be between 0 and 3");
273 
274  if (nummodesA == -1)
275  {
276  switch(fid)
277  {
278  case 0:
279  nummodesA = m_base[0]->GetNumModes();
280  nummodesB = m_base[1]->GetNumModes();
281  break;
282  case 1:
283  nummodesA = m_base[0]->GetNumModes();
284  nummodesB = m_base[2]->GetNumModes();
285  break;
286  case 2:
287  case 3:
288  nummodesA = m_base[1]->GetNumModes();
289  nummodesB = m_base[2]->GetNumModes();
290  break;
291  }
292  }
293 
294  P = nummodesA;
295  Q = nummodesB;
296  nFaceCoeffs = Q + ((P-1)*(1 + 2*(Q-1) - (P-1)))/2;
297 
298  if (maparray.size() != nFaceCoeffs)
299  {
300  maparray = Array<OneD, unsigned int>(nFaceCoeffs);
301  }
302 
303  if (signarray.size() != nFaceCoeffs)
304  {
305  signarray = Array<OneD, int>(nFaceCoeffs,1);
306  }
307  else
308  {
309  fill(signarray.get(), signarray.get()+nFaceCoeffs, 1);
310  }
311 
312  switch(fid)
313  {
314  case 0:
315  // Add vertices.
316  maparray[idx++] = 0;
317  maparray[idx++] = 1;
318  maparray[idx++] = 2;
319 
320  // Add edges.
321  for (i = 2; i < P; ++i)
322  {
323  maparray[idx++] = ;
324  }
325  }
326 }
327 */
328 
329 int StdNodalTetExp::v_GetVertexMap(const int localVertexId,
330  bool useCoeffPacking)
331 {
332  boost::ignore_unused(useCoeffPacking);
333  ASSERTL0(localVertexId >= 0 && localVertexId <= 3,
334  "Local Vertex ID must be between 0 and 3");
335  return localVertexId;
336 }
337 
339 {
340  unsigned int i;
341  const unsigned int nBndryCoeff = NumBndryCoeffs();
342 
343  if (outarray.size() != nBndryCoeff)
344  {
345  outarray = Array<OneD, unsigned int>(nBndryCoeff);
346  }
347 
348  for (i = 0; i < nBndryCoeff; i++)
349  {
350  outarray[i] = i;
351  }
352 }
353 
355 {
356  unsigned int i;
357  const unsigned int nBndryCoeff = NumBndryCoeffs();
358 
359  if (outarray.size() != m_ncoeffs - nBndryCoeff)
360  {
361  outarray = Array<OneD, unsigned int>(m_ncoeffs - nBndryCoeff);
362  }
363 
364  for (i = nBndryCoeff; i < m_ncoeffs; i++)
365  {
366  outarray[i - nBndryCoeff] = i;
367  }
368 }
369 
370 //---------------------------------------
371 // Wrapper functions
372 //---------------------------------------
373 
375 {
376  DNekMatSharedPtr Mat;
377 
378  switch (mkey.GetMatrixType())
379  {
380  case eNBasisTrans:
381  Mat = GenNBasisTransMatrix();
382  break;
383  default:
385  break;
386  }
387 
388  return Mat;
389 }
390 
392 {
393  return StdNodalTetExp::v_GenMatrix(mkey);
394 }
395 } // namespace StdRegions
396 } // 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
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:83
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 GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:85
void NodalToModalTranspose(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
int v_GetVertexMap(const int localVertexId, bool useCoeffPacking=false) override
void NodalToModal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void ModalToNodal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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.
void GetNodalPoints(Array< OneD, const NekDouble > &x, Array< OneD, const NekDouble > &y, Array< OneD, const NekDouble > &z)
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool mult=true) override
LibUtilities::PointsKey m_nodalPointsKey
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:824
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:633
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdTetExp.cpp:509
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:309
LibUtilities::ShapeType DetShapeType() const
Definition: StdTetExp.h:64
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