Nektar++
StdTetExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: StdTetExp.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: Header field for tetrahedral routines built upon
32 // StdExpansion3D
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <boost/core/ignore_unused.hpp>
37 
38 #include <StdRegions/StdTetExp.h>
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44 namespace StdRegions
45 {
46 StdTetExp::StdTetExp()
47 {
48 }
49 
50 StdTetExp::StdTetExp(const LibUtilities::BasisKey &Ba,
51  const LibUtilities::BasisKey &Bb,
52  const LibUtilities::BasisKey &Bc)
53  : StdExpansion(LibUtilities::StdTetData::getNumberOfCoefficients(
54  Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
55  3, Ba, Bb, Bc),
56  StdExpansion3D(LibUtilities::StdTetData::getNumberOfCoefficients(
57  Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
58  Ba, Bb, Bc)
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 
72 {
73 }
74 
76 {
77 }
78 
79 //----------------------------
80 // Differentiation Methods
81 //----------------------------
82 
83 /**
84  * \brief Calculate the derivative of the physical points
85  *
86  * The derivative is evaluated at the nodal physical points.
87  * Derivatives with respect to the local Cartesian coordinates
88  *
89  * \f$\begin{Bmatrix} \frac {\partial} {\partial \xi_1} \\ \frac
90  * {\partial} {\partial \xi_2} \\ \frac {\partial} {\partial \xi_3}
91  * \end{Bmatrix} = \begin{Bmatrix} \frac 4 {(1-\eta_2)(1-\eta_3)}
92  * \frac \partial {\partial \eta_1} \ \ \frac {2(1+\eta_1)}
93  * {(1-\eta_2)(1-\eta_3)} \frac \partial {\partial \eta_1} + \frac 2
94  * {1-\eta_3} \frac \partial {\partial \eta_3} \\ \frac {2(1 +
95  * \eta_1)} {2(1 - \eta_2)(1-\eta_3)} \frac \partial {\partial \eta_1}
96  * + \frac {1 + \eta_2} {1 - \eta_3} \frac \partial {\partial \eta_2}
97  * + \frac \partial {\partial \eta_3} \end{Bmatrix}\f$
98  **/
100  Array<OneD, NekDouble> &out_dxi0,
101  Array<OneD, NekDouble> &out_dxi1,
102  Array<OneD, NekDouble> &out_dxi2)
103 {
104  int Q0 = m_base[0]->GetNumPoints();
105  int Q1 = m_base[1]->GetNumPoints();
106  int Q2 = m_base[2]->GetNumPoints();
107  int Qtot = Q0 * Q1 * Q2;
108 
109  // Compute the physical derivative
110  Array<OneD, NekDouble> out_dEta0(3 * Qtot, 0.0);
111  Array<OneD, NekDouble> out_dEta1 = out_dEta0 + Qtot;
112  Array<OneD, NekDouble> out_dEta2 = out_dEta1 + Qtot;
113 
114  bool Do_2 = (out_dxi2.size() > 0) ? true : false;
115  bool Do_1 = (out_dxi1.size() > 0) ? true : false;
116 
117  if (Do_2) // Need all local derivatives
118  {
119  PhysTensorDeriv(inarray, out_dEta0, out_dEta1, out_dEta2);
120  }
121  else if (Do_1) // Need 0 and 1 derivatives
122  {
123  PhysTensorDeriv(inarray, out_dEta0, out_dEta1, NullNekDouble1DArray);
124  }
125  else // Only need Eta0 derivaitve
126  {
127  PhysTensorDeriv(inarray, out_dEta0, NullNekDouble1DArray,
129  }
130 
131  Array<OneD, const NekDouble> eta_0, eta_1, eta_2;
132  eta_0 = m_base[0]->GetZ();
133  eta_1 = m_base[1]->GetZ();
134  eta_2 = m_base[2]->GetZ();
135 
136  // calculate 2.0/((1-eta_1)(1-eta_2)) Out_dEta0
137 
138  NekDouble *dEta0 = &out_dEta0[0];
139  NekDouble fac;
140  for (int k = 0; k < Q2; ++k)
141  {
142  for (int j = 0; j < Q1; ++j, dEta0 += Q0)
143  {
144  Vmath::Smul(Q0, 2.0 / (1.0 - eta_1[j]), dEta0, 1, dEta0, 1);
145  }
146  fac = 1.0 / (1.0 - eta_2[k]);
147  Vmath::Smul(Q0 * Q1, fac, &out_dEta0[0] + k * Q0 * Q1, 1,
148  &out_dEta0[0] + k * Q0 * Q1, 1);
149  }
150 
151  if (out_dxi0.size() > 0)
152  {
153  // out_dxi0 = 4.0/((1-eta_1)(1-eta_2)) Out_dEta0
154  Vmath::Smul(Qtot, 2.0, out_dEta0, 1, out_dxi0, 1);
155  }
156 
157  if (Do_1 || Do_2)
158  {
159  Array<OneD, NekDouble> Fac0(Q0);
160  Vmath::Sadd(Q0, 1.0, eta_0, 1, Fac0, 1);
161 
162  // calculate 2.0*(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0
163  for (int k = 0; k < Q1 * Q2; ++k)
164  {
165  Vmath::Vmul(Q0, &Fac0[0], 1, &out_dEta0[0] + k * Q0, 1,
166  &out_dEta0[0] + k * Q0, 1);
167  }
168  // calculate 2/(1.0-eta_2) out_dEta1
169  for (int k = 0; k < Q2; ++k)
170  {
171  Vmath::Smul(Q0 * Q1, 2.0 / (1.0 - eta_2[k]),
172  &out_dEta1[0] + k * Q0 * Q1, 1,
173  &out_dEta1[0] + k * Q0 * Q1, 1);
174  }
175 
176  if (Do_1)
177  {
178  // calculate out_dxi1 = 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0
179  // + 2/(1.0-eta_2) out_dEta1
180  Vmath::Vadd(Qtot, out_dEta0, 1, out_dEta1, 1, out_dxi1, 1);
181  }
182 
183  if (Do_2)
184  {
185  // calculate (1 + eta_1)/(1 -eta_2)*out_dEta1
186  NekDouble *dEta1 = &out_dEta1[0];
187  for (int k = 0; k < Q2; ++k)
188  {
189  for (int j = 0; j < Q1; ++j, dEta1 += Q0)
190  {
191  Vmath::Smul(Q0, (1.0 + eta_1[j]) / 2.0, dEta1, 1, dEta1, 1);
192  }
193  }
194 
195  // calculate out_dxi2 =
196  // 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0 +
197  // (1 + eta_1)/(1 -eta_2)*out_dEta1 + out_dEta2
198  Vmath::Vadd(Qtot, out_dEta0, 1, out_dEta1, 1, out_dxi2, 1);
199  Vmath::Vadd(Qtot, out_dEta2, 1, out_dxi2, 1, out_dxi2, 1);
200  }
201  }
202 }
203 
204 /**
205  * @param dir Direction in which to compute derivative.
206  * Valid values are 0, 1, 2.
207  * @param inarray Input array.
208  * @param outarray Output array.
209  */
210 void StdTetExp::v_PhysDeriv(const int dir,
211  const Array<OneD, const NekDouble> &inarray,
212  Array<OneD, NekDouble> &outarray)
213 {
214  switch (dir)
215  {
216  case 0:
217  {
218  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
220  break;
221  }
222  case 1:
223  {
224  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
226  break;
227  }
228  case 2:
229  {
231  outarray);
232  break;
233  }
234  default:
235  {
236  ASSERTL1(false, "input dir is out of range");
237  }
238  break;
239  }
240 }
241 
243  Array<OneD, NekDouble> &out_d0,
244  Array<OneD, NekDouble> &out_d1,
245  Array<OneD, NekDouble> &out_d2)
246 {
247  StdTetExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
248 }
249 
250 void StdTetExp::v_StdPhysDeriv(const int dir,
251  const Array<OneD, const NekDouble> &inarray,
252  Array<OneD, NekDouble> &outarray)
253 {
254  StdTetExp::v_PhysDeriv(dir, inarray, outarray);
255 }
256 
257 //---------------------------------------
258 // Transforms
259 //---------------------------------------
260 
261 /**
262  * @note 'r' (base[2]) runs fastest in this element
263  *
264  * \f$ u^{\delta} (\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{m(pqr)} \hat
265  * u_{pqr} \phi_{pqr} (\xi_{1i}, \xi_{2j}, \xi_{3k})\f$
266  *
267  * Backward transformation is three dimensional tensorial expansion
268  * \f$ u (\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_p^a
269  * (\xi_{1i}) \lbrace { \sum_{q=0}^{Q_y} \psi_{pq}^b (\xi_{2j})
270  * \lbrace { \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c (\xi_{3k})
271  * \rbrace} \rbrace}. \f$ And sumfactorizing step of the form is as:\\
272  *
273  * \f$ f_{pq} (\xi_{3k}) = \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c
274  * (\xi_{3k}),\\
275  *
276  * g_{p} (\xi_{2j}, \xi_{3k}) = \sum_{r=0}^{Q_y} \psi_{pq}^b
277  * (\xi_{2j}) f_{pq} (\xi_{3k}),\\
278  *
279  * u(\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_{p}^a
280  * (\xi_{1i}) g_{p} (\xi_{2j}, \xi_{3k}). \f$
281  */
283  Array<OneD, NekDouble> &outarray)
284 {
287  "Basis[1] is not a general tensor type");
288 
291  "Basis[2] is not a general tensor type");
292 
293  if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
294  m_base[2]->Collocation())
295  {
297  m_base[2]->GetNumPoints(),
298  inarray, 1, outarray, 1);
299  }
300  else
301  {
302  StdTetExp::v_BwdTrans_SumFac(inarray, outarray);
303  }
304 }
305 
306 /**
307  * Sum-factorisation implementation of the BwdTrans operation.
308  */
310  Array<OneD, NekDouble> &outarray)
311 {
312  int nquad1 = m_base[1]->GetNumPoints();
313  int nquad2 = m_base[2]->GetNumPoints();
314  int order0 = m_base[0]->GetNumModes();
315  int order1 = m_base[1]->GetNumModes();
316 
317  Array<OneD, NekDouble> wsp(nquad2 * order0 * (2 * order1 - order0 + 1) / 2 +
318  nquad2 * nquad1 * order0);
319 
320  BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
321  m_base[2]->GetBdata(), inarray, outarray, wsp, true,
322  true, true);
323 }
324 
325 /**
326  * @param base0 x-dirn basis matrix
327  * @param base1 y-dirn basis matrix
328  * @param base2 z-dirn basis matrix
329  * @param inarray Input vector of modes.
330  * @param outarray Output vector of physical space data.
331  * @param wsp Workspace of size Q_x*P_z*(P_y+Q_y)
332  * @param doCheckCollDir0 Check for collocation of basis.
333  * @param doCheckCollDir1 Check for collocation of basis.
334  * @param doCheckCollDir2 Check for collocation of basis.
335  * @todo Account for some directions being collocated. See
336  * StdQuadExp as an example.
337  */
339  const Array<OneD, const NekDouble> &base0,
340  const Array<OneD, const NekDouble> &base1,
341  const Array<OneD, const NekDouble> &base2,
342  const Array<OneD, const NekDouble> &inarray,
344  bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
345 {
346  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
347 
348  int nquad0 = m_base[0]->GetNumPoints();
349  int nquad1 = m_base[1]->GetNumPoints();
350  int nquad2 = m_base[2]->GetNumPoints();
351 
352  int order0 = m_base[0]->GetNumModes();
353  int order1 = m_base[1]->GetNumModes();
354  int order2 = m_base[2]->GetNumModes();
355 
356  Array<OneD, NekDouble> tmp = wsp;
358  tmp + nquad2 * order0 * (2 * order1 - order0 + 1) / 2;
359 
360  int i, j, mode, mode1, cnt;
361 
362  // Perform summation over '2' direction
363  mode = mode1 = cnt = 0;
364  for (i = 0; i < order0; ++i)
365  {
366  for (j = 0; j < order1 - i; ++j, ++cnt)
367  {
368  Blas::Dgemv('N', nquad2, order2 - i - j, 1.0,
369  base2.get() + mode * nquad2, nquad2,
370  inarray.get() + mode1, 1, 0.0, tmp.get() + cnt * nquad2,
371  1);
372  mode += order2 - i - j;
373  mode1 += order2 - i - j;
374  }
375  // increment mode in case order1!=order2
376  for (j = order1 - i; j < order2 - i; ++j)
377  {
378  mode += order2 - i - j;
379  }
380  }
381 
382  // fix for modified basis by adding split of top singular
383  // vertex mode - currently (1+c)/2 x (1-b)/2 x (1-a)/2
384  // component is evaluated
386  {
387  // top singular vertex - (1+c)/2 x (1+b)/2 x (1-a)/2 component
388  Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
389  &tmp[0] + nquad2, 1);
390 
391  // top singular vertex - (1+c)/2 x (1-b)/2 x (1+a)/2 component
392  Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
393  &tmp[0] + order1 * nquad2, 1);
394  }
395 
396  // Perform summation over '1' direction
397  mode = 0;
398  for (i = 0; i < order0; ++i)
399  {
400  Blas::Dgemm('N', 'T', nquad1, nquad2, order1 - i, 1.0,
401  base1.get() + mode * nquad1, nquad1,
402  tmp.get() + mode * nquad2, nquad2, 0.0,
403  tmp1.get() + i * nquad1 * nquad2, nquad1);
404  mode += order1 - i;
405  }
406 
407  // fix for modified basis by adding additional split of
408  // top and base singular vertex modes as well as singular
409  // edge
411  {
412  // use tmp to sort out singular vertices and
413  // singular edge components with (1+b)/2 (1+a)/2 form
414  for (i = 0; i < nquad2; ++i)
415  {
416  Blas::Daxpy(nquad1, tmp[nquad2 + i], base1.get() + nquad1, 1,
417  &tmp1[nquad1 * nquad2] + i * nquad1, 1);
418  }
419  }
420 
421  // Perform summation over '0' direction
422  Blas::Dgemm('N', 'T', nquad0, nquad1 * nquad2, order0, 1.0, base0.get(),
423  nquad0, tmp1.get(), nquad1 * nquad2, 0.0, outarray.get(),
424  nquad0);
425 }
426 
427 /**
428  * @param inarray array of physical quadrature points to be
429  * transformed.
430  * @param outarray updated array of expansion coefficients.
431  */
433  Array<OneD, NekDouble> &outarray)
434 { // int numMax = nmodes0;
435  v_IProductWRTBase(inarray, outarray);
436 
437  // get Mass matrix inverse
438  StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
439  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
440 
441  // copy inarray in case inarray == outarray
442  DNekVec in(m_ncoeffs, outarray);
443  DNekVec out(m_ncoeffs, outarray, eWrapper);
444 
445  out = (*matsys) * in;
446 }
447 
448 //---------------------------------------
449 // Inner product functions
450 //---------------------------------------
451 
452 /**
453  * \f$ \begin{array}{rcl} I_{pqr} = (\phi_{pqr}, u)_{\delta} & = &
454  * \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \sum_{k=0}^{nq_2} \psi_{p}^{a}
455  * (\eta_{1i}) \psi_{pq}^{b} (\eta_{2j}) \psi_{pqr}^{c} (\eta_{3k})
456  * w_i w_j w_k u(\eta_{1,i} \eta_{2,j} \eta_{3,k}) J_{i,j,k}\\ & = &
457  * \sum_{i=0}^{nq_0} \psi_p^a(\eta_{1,i}) \sum_{j=0}^{nq_1}
458  * \psi_{pq}^b(\eta_{2,j}) \sum_{k=0}^{nq_2} \psi_{pqr}^c
459  * u(\eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} \end{array} \f$ \n
460  *
461  * where
462  *
463  * \f$ \phi_{pqr} (\xi_1 , \xi_2 , \xi_3) = \psi_p^a (\eta_1)
464  * \psi_{pq}^b (\eta_2) \psi_{pqr}^c (\eta_3) \f$
465  *
466  * which can be implemented as \n \f$f_{pqr} (\xi_{3k}) =
467  * \sum_{k=0}^{nq_3} \psi_{pqr}^c u(\eta_{1i},\eta_{2j},\eta_{3k})
468  *
469  * J_{i,j,k} = {\bf B_3 U} \f$ \n
470  *
471  * \f$ g_{pq} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{pq}^b (\xi_{2j})
472  * f_{pqr} (\xi_{3k}) = {\bf B_2 F} \f$ \n
473  *
474  * \f$ (\phi_{pqr}, u)_{\delta} = \sum_{k=0}^{nq_0} \psi_{p}^a
475  * (\xi_{3k}) g_{pq} (\xi_{3k}) = {\bf B_1 G} \f$
476  *
477  * @param inarray Function evaluated at physical collocation
478  * points.
479  * @param outarray Inner product with respect to each basis
480  * function over the element.
481  */
483  Array<OneD, NekDouble> &outarray)
484 {
487  "Basis[1] is not a general tensor type");
488 
491  "Basis[2] is not a general tensor type");
492 
493  if (m_base[0]->Collocation() && m_base[1]->Collocation())
494  {
495  MultiplyByQuadratureMetric(inarray, outarray);
496  }
497  else
498  {
499  StdTetExp::v_IProductWRTBase_SumFac(inarray, outarray);
500  }
501 }
502 
503 /**
504  * @param inarray Function evaluated at physical collocation
505  * points.
506  * @param outarray Inner product with respect to each basis
507  * function over the element.
508  */
510  const Array<OneD, const NekDouble> &inarray,
511  Array<OneD, NekDouble> &outarray, bool multiplybyweights)
512 {
513  int nquad0 = m_base[0]->GetNumPoints();
514  int nquad1 = m_base[1]->GetNumPoints();
515  int nquad2 = m_base[2]->GetNumPoints();
516  int order0 = m_base[0]->GetNumModes();
517  int order1 = m_base[1]->GetNumModes();
518 
519  Array<OneD, NekDouble> wsp(nquad1 * nquad2 * order0 +
520  nquad2 * order0 * (2 * order1 - order0 + 1) / 2);
521 
522  if (multiplybyweights)
523  {
524  Array<OneD, NekDouble> tmp(nquad0 * nquad1 * nquad2);
525  MultiplyByQuadratureMetric(inarray, tmp);
526 
528  m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
529  tmp, outarray, wsp, true, true, true);
530  }
531  else
532  {
534  m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
535  inarray, outarray, wsp, true, true, true);
536  }
537 }
538 
540  const Array<OneD, const NekDouble> &base0,
541  const Array<OneD, const NekDouble> &base1,
542  const Array<OneD, const NekDouble> &base2,
543  const Array<OneD, const NekDouble> &inarray,
545  bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
546 {
547  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
548 
549  int nquad0 = m_base[0]->GetNumPoints();
550  int nquad1 = m_base[1]->GetNumPoints();
551  int nquad2 = m_base[2]->GetNumPoints();
552 
553  int order0 = m_base[0]->GetNumModes();
554  int order1 = m_base[1]->GetNumModes();
555  int order2 = m_base[2]->GetNumModes();
556 
557  Array<OneD, NekDouble> tmp1 = wsp;
558  Array<OneD, NekDouble> tmp2 = wsp + nquad1 * nquad2 * order0;
559 
560  int i, j, mode, mode1, cnt;
561 
562  // Inner product with respect to the '0' direction
563  Blas::Dgemm('T', 'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.get(),
564  nquad0, base0.get(), nquad0, 0.0, tmp1.get(), nquad1 * nquad2);
565 
566  // Inner product with respect to the '1' direction
567  for (mode = i = 0; i < order0; ++i)
568  {
569  Blas::Dgemm('T', 'N', nquad2, order1 - i, nquad1, 1.0,
570  tmp1.get() + i * nquad1 * nquad2, nquad1,
571  base1.get() + mode * nquad1, nquad1, 0.0,
572  tmp2.get() + mode * nquad2, nquad2);
573  mode += order1 - i;
574  }
575 
576  // fix for modified basis for base singular vertex
578  {
579  // base singular vertex and singular edge (1+b)/2
580  //(1+a)/2 components (makes tmp[nquad2] entry into (1+b)/2)
581  Blas::Dgemv('T', nquad1, nquad2, 1.0, tmp1.get() + nquad1 * nquad2,
582  nquad1, base1.get() + nquad1, 1, 1.0, tmp2.get() + nquad2,
583  1);
584  }
585 
586  // Inner product with respect to the '2' direction
587  mode = mode1 = cnt = 0;
588  for (i = 0; i < order0; ++i)
589  {
590  for (j = 0; j < order1 - i; ++j, ++cnt)
591  {
592  Blas::Dgemv('T', nquad2, order2 - i - j, 1.0,
593  base2.get() + mode * nquad2, nquad2,
594  tmp2.get() + cnt * nquad2, 1, 0.0,
595  outarray.get() + mode1, 1);
596  mode += order2 - i - j;
597  mode1 += order2 - i - j;
598  }
599  // increment mode in case order1!=order2
600  for (j = order1 - i; j < order2 - i; ++j)
601  {
602  mode += order2 - i - j;
603  }
604  }
605 
606  // fix for modified basis for top singular vertex component
607  // Already have evaluated (1+c)/2 (1-b)/2 (1-a)/2
609  {
610  // add in (1+c)/2 (1+b)/2 component
611  outarray[1] +=
612  Blas::Ddot(nquad2, base2.get() + nquad2, 1, &tmp2[nquad2], 1);
613 
614  // add in (1+c)/2 (1-b)/2 (1+a)/2 component
615  outarray[1] += Blas::Ddot(nquad2, base2.get() + nquad2, 1,
616  &tmp2[nquad2 * order1], 1);
617  }
618 }
619 
621  const int dir, const Array<OneD, const NekDouble> &inarray,
622  Array<OneD, NekDouble> &outarray)
623 {
624  StdTetExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
625 }
626 
627 /**
628  * @param inarray Function evaluated at physical collocation
629  * points.
630  * @param outarray Inner product with respect to each basis
631  * function over the element.
632  */
634  const int dir, const Array<OneD, const NekDouble> &inarray,
635  Array<OneD, NekDouble> &outarray)
636 {
637  int i;
638  int nquad0 = m_base[0]->GetNumPoints();
639  int nquad1 = m_base[1]->GetNumPoints();
640  int nquad2 = m_base[2]->GetNumPoints();
641  int nqtot = nquad0 * nquad1 * nquad2;
642  int nmodes0 = m_base[0]->GetNumModes();
643  int nmodes1 = m_base[1]->GetNumModes();
644  int wspsize = nquad0 + nquad1 + nquad2 + max(nqtot, m_ncoeffs) +
645  nquad1 * nquad2 * nmodes0 +
646  nquad2 * nmodes0 * (2 * nmodes1 - nmodes0 + 1) / 2;
647 
648  Array<OneD, NekDouble> gfac0(wspsize);
649  Array<OneD, NekDouble> gfac1(gfac0 + nquad0);
650  Array<OneD, NekDouble> gfac2(gfac1 + nquad1);
651  Array<OneD, NekDouble> tmp0(gfac2 + nquad2);
652  Array<OneD, NekDouble> wsp(tmp0 + max(nqtot, m_ncoeffs));
653 
654  const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
655  const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
656  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
657 
658  // set up geometric factor: (1+z0)/2
659  for (i = 0; i < nquad0; ++i)
660  {
661  gfac0[i] = 0.5 * (1 + z0[i]);
662  }
663 
664  // set up geometric factor: 2/(1-z1)
665  for (i = 0; i < nquad1; ++i)
666  {
667  gfac1[i] = 2.0 / (1 - z1[i]);
668  }
669 
670  // Set up geometric factor: 2/(1-z2)
671  for (i = 0; i < nquad2; ++i)
672  {
673  gfac2[i] = 2.0 / (1 - z2[i]);
674  }
675 
676  // Derivative in first direction is always scaled as follows
677  for (i = 0; i < nquad1 * nquad2; ++i)
678  {
679  Vmath::Smul(nquad0, gfac1[i % nquad1], &inarray[0] + i * nquad0, 1,
680  &tmp0[0] + i * nquad0, 1);
681  }
682  for (i = 0; i < nquad2; ++i)
683  {
684  Vmath::Smul(nquad0 * nquad1, gfac2[i], &tmp0[0] + i * nquad0 * nquad1,
685  1, &tmp0[0] + i * nquad0 * nquad1, 1);
686  }
687 
688  MultiplyByQuadratureMetric(tmp0, tmp0);
689 
690  switch (dir)
691  {
692  case 0:
693  {
695  m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
696  m_base[2]->GetBdata(), tmp0, outarray, wsp, false, true, true);
697  }
698  break;
699  case 1:
700  {
702 
703  for (i = 0; i < nquad1 * nquad2; ++i)
704  {
705  Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
706  &tmp0[0] + i * nquad0, 1);
707  }
708 
710  m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
711  m_base[2]->GetBdata(), tmp0, tmp3, wsp, false, true, true);
712 
713  for (i = 0; i < nquad2; ++i)
714  {
715  Vmath::Smul(nquad0 * nquad1, gfac2[i],
716  &inarray[0] + i * nquad0 * nquad1, 1,
717  &tmp0[0] + i * nquad0 * nquad1, 1);
718  }
719  MultiplyByQuadratureMetric(tmp0, tmp0);
721  m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
722  m_base[2]->GetBdata(), tmp0, outarray, wsp, true, false, true);
723  Vmath::Vadd(m_ncoeffs, &tmp3[0], 1, &outarray[0], 1, &outarray[0],
724  1);
725  }
726  break;
727  case 2:
728  {
731  for (i = 0; i < nquad1; ++i)
732  {
733  gfac1[i] = (1 + z1[i]) / 2;
734  }
735 
736  for (i = 0; i < nquad1 * nquad2; ++i)
737  {
738  Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
739  &tmp0[0] + i * nquad0, 1);
740  }
742  m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
743  m_base[2]->GetBdata(), tmp0, tmp3, wsp, false, true, true);
744 
745  for (i = 0; i < nquad2; ++i)
746  {
747  Vmath::Smul(nquad0 * nquad1, gfac2[i],
748  &inarray[0] + i * nquad0 * nquad1, 1,
749  &tmp0[0] + i * nquad0 * nquad1, 1);
750  }
751  for (i = 0; i < nquad1 * nquad2; ++i)
752  {
753  Vmath::Smul(nquad0, gfac1[i % nquad1], &tmp0[0] + i * nquad0, 1,
754  &tmp0[0] + i * nquad0, 1);
755  }
756  MultiplyByQuadratureMetric(tmp0, tmp0);
758  m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
759  m_base[2]->GetBdata(), tmp0, tmp4, wsp, true, false, true);
760 
761  MultiplyByQuadratureMetric(inarray, tmp0);
763  m_base[0]->GetBdata(), m_base[1]->GetBdata(),
764  m_base[2]->GetDbdata(), tmp0, outarray, wsp, true, true, false);
765 
766  Vmath::Vadd(m_ncoeffs, &tmp3[0], 1, &outarray[0], 1, &outarray[0],
767  1);
768  Vmath::Vadd(m_ncoeffs, &tmp4[0], 1, &outarray[0], 1, &outarray[0],
769  1);
770  }
771  break;
772  default:
773  {
774  ASSERTL1(false, "input dir is out of range");
775  }
776  break;
777  }
778 }
779 
780 //---------------------------------------
781 // Evaluation functions
782 //---------------------------------------
783 
786 {
787  NekDouble d2 = 1.0 - xi[2];
788  NekDouble d12 = -xi[1] - xi[2];
789  if (fabs(d2) < NekConstants::kNekZeroTol)
790  {
791  if (d2 >= 0.)
792  {
794  }
795  else
796  {
798  }
799  }
800  if (fabs(d12) < NekConstants::kNekZeroTol)
801  {
802  if (d12 >= 0.)
803  {
805  }
806  else
807  {
809  }
810  }
811  eta[0] = 2.0 * (1.0 + xi[0]) / d12 - 1.0;
812  eta[1] = 2.0 * (1.0 + xi[1]) / d2 - 1.0;
813  eta[2] = xi[2];
814 }
815 
818 {
819  xi[1] = (1.0 + eta[0]) * (1.0 - eta[2]) * 0.5 - 1.0;
820  xi[0] = (1.0 + eta[0]) * (-xi[1] - eta[2]) * 0.5 - 1.0;
821  xi[2] = eta[2];
822 }
823 
824 void StdTetExp::v_FillMode(const int mode, Array<OneD, NekDouble> &outarray)
825 {
827  tmp[mode] = 1.0;
828  StdTetExp::v_BwdTrans(tmp, outarray);
829 }
830 
832  const Array<OneD, const NekDouble> &coords, int mode)
833 {
834  Array<OneD, NekDouble> coll(3);
835  LocCoordToLocCollapsed(coords, coll);
836 
837  const int nm1 = m_base[1]->GetNumModes();
838  const int nm2 = m_base[2]->GetNumModes();
839 
840  const int b = 2 * nm2 + 1;
841  const int mode0 = floor(0.5 * (b - sqrt(b * b - 8.0 * mode / nm1)));
842  const int tmp =
843  mode - nm1 * (mode0 * (nm2 - 1) + 1 - (mode0 - 2) * (mode0 - 1) / 2);
844  const int mode1 = tmp / (nm2 - mode0);
845  const int mode2 = tmp % (nm2 - mode0);
846 
848  {
849  // Handle the collapsed vertices and edges in the modified
850  // basis.
851  if (mode == 1)
852  {
853  // Collapsed top vertex
854  return StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
855  }
856  else if (mode0 == 0 && mode2 == 1)
857  {
858  return StdExpansion::BaryEvaluateBasis<1>(coll[1], 0) *
859  StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
860  }
861  else if (mode0 == 1 && mode1 == 1 && mode2 == 0)
862  {
863  return StdExpansion::BaryEvaluateBasis<0>(coll[0], 0) *
864  StdExpansion::BaryEvaluateBasis<1>(coll[1], 1);
865  }
866  }
867 
868  return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
869  StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
870  StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
871 }
872 
874  const Array<OneD, const NekDouble> &inarray,
875  std::array<NekDouble, 3> &firstOrderDerivs)
876 {
877  // Collapse coordinates
878  Array<OneD, NekDouble> coll(3, 0.0);
879  LocCoordToLocCollapsed(coord, coll);
880 
881  // If near singularity do the old interpolation matrix method
882  if ((1 - coll[1]) < 1e-5 || (1 - coll[2]) < 1e-5)
883  {
884  int totPoints = GetTotPoints();
885  Array<OneD, NekDouble> EphysDeriv0(totPoints), EphysDeriv1(totPoints),
886  EphysDeriv2(totPoints);
887  PhysDeriv(inarray, EphysDeriv0, EphysDeriv1, EphysDeriv2);
888 
890  I[0] = GetBase()[0]->GetI(coll);
891  I[1] = GetBase()[1]->GetI(coll + 1);
892  I[2] = GetBase()[2]->GetI(coll + 2);
893 
894  firstOrderDerivs[0] = PhysEvaluate(I, EphysDeriv0);
895  firstOrderDerivs[1] = PhysEvaluate(I, EphysDeriv1);
896  firstOrderDerivs[2] = PhysEvaluate(I, EphysDeriv2);
897  return PhysEvaluate(I, inarray);
898  }
899 
900  std::array<NekDouble, 3> interDeriv;
901  NekDouble val = BaryTensorDeriv(coll, inarray, interDeriv);
902 
903  // calculate 2.0/((1-eta_1)(1-eta_2)) * Out_dEta0
904  NekDouble temp = 2.0 / ((1 - coll[1]) * (1 - coll[2]));
905  interDeriv[0] *= temp;
906 
907  // out_dxi0 = 4.0/((1-eta_1)(1-eta_2)) * Out_dEta0
908  firstOrderDerivs[0] = 2 * interDeriv[0];
909 
910  // fac0 = 1 + eta_0
911  NekDouble fac0;
912  fac0 = 1 + coll[0];
913 
914  // calculate 2.0*(1+eta_0)/((1-eta_1)(1-eta_2)) * Out_dEta0
915  interDeriv[0] *= fac0;
916 
917  // calculate 2/(1.0-eta_2) * out_dEta1
918  fac0 = 2 / (1 - coll[2]);
919  interDeriv[1] *= fac0;
920 
921  // calculate out_dxi1 = 2.0(1+eta_0)/((1-eta_1)(1-eta_2))
922  // * Out_dEta0 + 2/(1.0-eta_2) out_dEta1
923  firstOrderDerivs[1] = interDeriv[0] + interDeriv[1];
924 
925  // calculate (1 + eta_1)/(1 -eta_2)*out_dEta1
926  fac0 = (1 + coll[1]) / 2;
927  interDeriv[1] *= fac0;
928 
929  // calculate out_dxi2 =
930  // 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0 +
931  // (1 + eta_1)/(1 -eta_2)*out_dEta1 + out_dEta2
932  firstOrderDerivs[2] = interDeriv[0] + interDeriv[1] + interDeriv[2];
933 
934  return val;
935 }
936 
937 void StdTetExp::v_GetTraceNumModes(const int fid,
938 
939  int &numModes0, int &numModes1,
940  Orientation faceOrient)
941 {
942  boost::ignore_unused(faceOrient);
943 
944  int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
945  m_base[2]->GetNumModes()};
946  switch (fid)
947  {
948  case 0:
949  {
950  numModes0 = nummodes[0];
951  numModes1 = nummodes[1];
952  }
953  break;
954  case 1:
955  {
956  numModes0 = nummodes[0];
957  numModes1 = nummodes[2];
958  }
959  break;
960  case 2:
961  case 3:
962  {
963  numModes0 = nummodes[1];
964  numModes1 = nummodes[2];
965  }
966  break;
967  }
968 }
969 
970 //---------------------------
971 // Helper functions
972 //---------------------------
973 
975 {
976  return 4;
977 }
978 
980 {
981  return 6;
982 }
983 
985 {
986  return 4;
987 }
988 
990 {
991  return DetShapeType();
992 }
993 
995 {
998  "BasisType is not a boundary interior form");
1001  "BasisType is not a boundary interior form");
1004  "BasisType is not a boundary interior form");
1005 
1006  int P = m_base[0]->GetNumModes();
1007  int Q = m_base[1]->GetNumModes();
1008  int R = m_base[2]->GetNumModes();
1009 
1011 }
1012 
1014 {
1017  "BasisType is not a boundary interior form");
1020  "BasisType is not a boundary interior form");
1023  "BasisType is not a boundary interior form");
1024 
1025  int P = m_base[0]->GetNumModes() - 1;
1026  int Q = m_base[1]->GetNumModes() - 1;
1027  int R = m_base[2]->GetNumModes() - 1;
1028 
1029  return (Q + 1) + P * (1 + 2 * Q - P) / 2 // base face
1030  + (R + 1) + P * (1 + 2 * R - P) / 2 // front face
1031  + 2 * (R + 1) + Q * (1 + 2 * R - Q); // back two faces
1032 }
1033 
1034 int StdTetExp::v_GetTraceNcoeffs(const int i) const
1035 {
1036  ASSERTL2((i >= 0) && (i <= 3), "face id is out of range");
1037  int nFaceCoeffs = 0;
1038  int nummodesA, nummodesB, P, Q;
1039  if (i == 0)
1040  {
1041  nummodesA = GetBasisNumModes(0);
1042  nummodesB = GetBasisNumModes(1);
1043  }
1044  else if ((i == 1) || (i == 2))
1045  {
1046  nummodesA = GetBasisNumModes(0);
1047  nummodesB = GetBasisNumModes(2);
1048  }
1049  else
1050  {
1051  nummodesA = GetBasisNumModes(1);
1052  nummodesB = GetBasisNumModes(2);
1053  }
1054  P = nummodesA - 1;
1055  Q = nummodesB - 1;
1056  nFaceCoeffs = Q + 1 + (P * (1 + 2 * Q - P)) / 2;
1057  return nFaceCoeffs;
1058 }
1059 
1060 int StdTetExp::v_GetTraceIntNcoeffs(const int i) const
1061 {
1062  ASSERTL2((i >= 0) && (i <= 3), "face id is out of range");
1063  int Pi = m_base[0]->GetNumModes() - 2;
1064  int Qi = m_base[1]->GetNumModes() - 2;
1065  int Ri = m_base[2]->GetNumModes() - 2;
1066 
1067  if ((i == 0))
1068  {
1069  return Pi * (2 * Qi - Pi - 1) / 2;
1070  }
1071  else if ((i == 1))
1072  {
1073  return Pi * (2 * Ri - Pi - 1) / 2;
1074  }
1075  else
1076  {
1077  return Qi * (2 * Ri - Qi - 1) / 2;
1078  }
1079 }
1080 
1081 int StdTetExp::v_GetTraceNumPoints(const int i) const
1082 {
1083  ASSERTL2(i >= 0 && i <= 3, "face id is out of range");
1084 
1085  if (i == 0)
1086  {
1087  return m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints();
1088  }
1089  else if (i == 1)
1090  {
1091  return m_base[0]->GetNumPoints() * m_base[2]->GetNumPoints();
1092  }
1093  else
1094  {
1095  return m_base[1]->GetNumPoints() * m_base[2]->GetNumPoints();
1096  }
1097 }
1098 
1099 int StdTetExp::v_GetEdgeNcoeffs(const int i) const
1100 {
1101  ASSERTL2((i >= 0) && (i <= 5), "edge id is out of range");
1102  int P = m_base[0]->GetNumModes();
1103  int Q = m_base[1]->GetNumModes();
1104  int R = m_base[2]->GetNumModes();
1105 
1106  if (i == 0)
1107  {
1108  return P;
1109  }
1110  else if (i == 1 || i == 2)
1111  {
1112  return Q;
1113  }
1114  else
1115  {
1116  return R;
1117  }
1118 }
1119 
1121  const int j) const
1122 {
1123  ASSERTL2(i >= 0 && i <= 3, "face id is out of range");
1124  ASSERTL2(j == 0 || j == 1, "face direction is out of range");
1125 
1126  if (i == 0)
1127  {
1128  return m_base[j]->GetPointsKey();
1129  }
1130  else if (i == 1)
1131  {
1132  return m_base[2 * j]->GetPointsKey();
1133  }
1134  else
1135  {
1136  return m_base[j + 1]->GetPointsKey();
1137  }
1138 }
1139 
1141  const std::vector<unsigned int> &nummodes, int &modes_offset)
1142 {
1144  nummodes[modes_offset], nummodes[modes_offset + 1],
1145  nummodes[modes_offset + 2]);
1146  modes_offset += 3;
1147 
1148  return nmodes;
1149 }
1150 
1152  const int k) const
1153 {
1154  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
1155  ASSERTL2(k == 0 || k == 1, "face direction out of range");
1156 
1157  int dir = k;
1158  switch (i)
1159  {
1160  case 0:
1161  dir = k;
1162  break;
1163  case 1:
1164  dir = 2 * k;
1165  break;
1166  case 2:
1167  case 3:
1168  dir = k + 1;
1169  break;
1170  }
1171 
1172  return EvaluateTriFaceBasisKey(k, m_base[dir]->GetBasisType(),
1173  m_base[dir]->GetNumPoints(),
1174  m_base[dir]->GetNumModes());
1175 
1176  // Should not get here.
1178 }
1179 
1181  Array<OneD, NekDouble> &xi_y,
1182  Array<OneD, NekDouble> &xi_z)
1183 {
1184  Array<OneD, const NekDouble> eta_x = m_base[0]->GetZ();
1185  Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
1186  Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
1187  int Qx = GetNumPoints(0);
1188  int Qy = GetNumPoints(1);
1189  int Qz = GetNumPoints(2);
1190 
1191  // Convert collapsed coordinates into cartesian coordinates: eta
1192  // --> xi
1193  for (int k = 0; k < Qz; ++k)
1194  {
1195  for (int j = 0; j < Qy; ++j)
1196  {
1197  for (int i = 0; i < Qx; ++i)
1198  {
1199  int s = i + Qx * (j + Qy * k);
1200  xi_x[s] =
1201  (eta_x[i] + 1.0) * (1.0 - eta_y[j]) * (1.0 - eta_z[k]) / 4 -
1202  1.0;
1203  xi_y[s] = (eta_y[j] + 1.0) * (1.0 - eta_z[k]) / 2 - 1.0;
1204  xi_z[s] = eta_z[k];
1205  }
1206  }
1207  }
1208 }
1209 
1211 {
1212  return (m_base[0]->GetBasisType() == LibUtilities::eModified_A) &&
1213  (m_base[1]->GetBasisType() == LibUtilities::eModified_B) &&
1215 }
1216 
1217 //--------------------------
1218 // Mappings
1219 //--------------------------
1220 int StdTetExp::v_GetVertexMap(const int localVertexId, bool useCoeffPacking)
1221 {
1225  "Mapping not defined for this type of basis");
1226 
1227  int localDOF = 0;
1228  if (useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
1229  {
1230  switch (localVertexId)
1231  {
1232  case 0:
1233  {
1234  localDOF = GetMode(0, 0, 0);
1235  break;
1236  }
1237  case 1:
1238  {
1239  localDOF = GetMode(0, 0, 1);
1240  break;
1241  }
1242  case 2:
1243  {
1244  localDOF = GetMode(0, 1, 0);
1245  break;
1246  }
1247  case 3:
1248  {
1249  localDOF = GetMode(1, 0, 0);
1250  break;
1251  }
1252  default:
1253  {
1254  ASSERTL0(false, "Vertex ID must be between 0 and 3");
1255  break;
1256  }
1257  }
1258  }
1259  else
1260  {
1261  switch (localVertexId)
1262  {
1263  case 0:
1264  {
1265  localDOF = GetMode(0, 0, 0);
1266  break;
1267  }
1268  case 1:
1269  {
1270  localDOF = GetMode(1, 0, 0);
1271  break;
1272  }
1273  case 2:
1274  {
1275  localDOF = GetMode(0, 1, 0);
1276  break;
1277  }
1278  case 3:
1279  {
1280  localDOF = GetMode(0, 0, 1);
1281  break;
1282  }
1283  default:
1284  {
1285  ASSERTL0(false, "Vertex ID must be between 0 and 3");
1286  break;
1287  }
1288  }
1289  }
1290 
1291  return localDOF;
1292 }
1293 
1294 /**
1295  * Maps interior modes of an edge to the elemental modes.
1296  */
1297 
1298 /**
1299  * List of all interior modes in the expansion.
1300  */
1302 {
1305  "BasisType is not a boundary interior form");
1308  "BasisType is not a boundary interior form");
1311  "BasisType is not a boundary interior form");
1312 
1313  int P = m_base[0]->GetNumModes();
1314  int Q = m_base[1]->GetNumModes();
1315  int R = m_base[2]->GetNumModes();
1316 
1317  int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
1318 
1319  if (outarray.size() != nIntCoeffs)
1320  {
1321  outarray = Array<OneD, unsigned int>(nIntCoeffs);
1322  }
1323 
1324  int idx = 0;
1325  for (int i = 2; i < P - 2; ++i)
1326  {
1327  for (int j = 1; j < Q - i - 1; ++j)
1328  {
1329  for (int k = 1; k < R - i - j; ++k)
1330  {
1331  outarray[idx++] = GetMode(i, j, k);
1332  }
1333  }
1334  }
1335 }
1336 
1337 /**
1338  * List of all boundary modes in the the expansion.
1339  */
1341 {
1344  "BasisType is not a boundary interior form");
1347  "BasisType is not a boundary interior form");
1350  "BasisType is not a boundary interior form");
1351 
1352  int P = m_base[0]->GetNumModes();
1353  int Q = m_base[1]->GetNumModes();
1354  int R = m_base[2]->GetNumModes();
1355 
1356  int i, j, k;
1357  int idx = 0;
1358 
1359  int nBnd = NumBndryCoeffs();
1360 
1361  if (outarray.size() != nBnd)
1362  {
1363  outarray = Array<OneD, unsigned int>(nBnd);
1364  }
1365 
1366  for (i = 0; i < P; ++i)
1367  {
1368  // First two Q-R planes are entirely boundary modes
1369  if (i < 2)
1370  {
1371  for (j = 0; j < Q - i; j++)
1372  {
1373  for (k = 0; k < R - i - j; ++k)
1374  {
1375  outarray[idx++] = GetMode(i, j, k);
1376  }
1377  }
1378  }
1379  // Remaining Q-R planes contain boundary modes on bottom and
1380  // left edge.
1381  else
1382  {
1383  for (k = 0; k < R - i; ++k)
1384  {
1385  outarray[idx++] = GetMode(i, 0, k);
1386  }
1387  for (j = 1; j < Q - i; ++j)
1388  {
1389  outarray[idx++] = GetMode(i, j, 0);
1390  }
1391  }
1392  }
1393 }
1394 
1395 void StdTetExp::v_GetTraceCoeffMap(const unsigned int fid,
1396  Array<OneD, unsigned int> &maparray)
1397 {
1398  int i, j, k;
1399  int P = 0, Q = 0, idx = 0;
1400  int nFaceCoeffs = 0;
1401 
1402  switch (fid)
1403  {
1404  case 0:
1405  P = m_base[0]->GetNumModes();
1406  Q = m_base[1]->GetNumModes();
1407  break;
1408  case 1:
1409  P = m_base[0]->GetNumModes();
1410  Q = m_base[2]->GetNumModes();
1411  break;
1412  case 2:
1413  case 3:
1414  P = m_base[1]->GetNumModes();
1415  Q = m_base[2]->GetNumModes();
1416  break;
1417  default:
1418  ASSERTL0(false, "fid must be between 0 and 3");
1419  }
1420 
1421  nFaceCoeffs = P * (2 * Q - P + 1) / 2;
1422 
1423  if (maparray.size() != nFaceCoeffs)
1424  {
1425  maparray = Array<OneD, unsigned int>(nFaceCoeffs);
1426  }
1427 
1428  switch (fid)
1429  {
1430  case 0:
1431  idx = 0;
1432  for (i = 0; i < P; ++i)
1433  {
1434  for (j = 0; j < Q - i; ++j)
1435  {
1436  maparray[idx++] = GetMode(i, j, 0);
1437  }
1438  }
1439  break;
1440  case 1:
1441  idx = 0;
1442  for (i = 0; i < P; ++i)
1443  {
1444  for (k = 0; k < Q - i; ++k)
1445  {
1446  maparray[idx++] = GetMode(i, 0, k);
1447  }
1448  }
1449  break;
1450  case 2:
1451  idx = 0;
1452  for (j = 0; j < P - 1; ++j)
1453  {
1454  for (k = 0; k < Q - 1 - j; ++k)
1455  {
1456  maparray[idx++] = GetMode(1, j, k);
1457  // Incorporate modes from zeroth plane where needed.
1458  if (j == 0 && k == 0)
1459  {
1460  maparray[idx++] = GetMode(0, 0, 1);
1461  }
1462  if (j == 0 && k == Q - 2)
1463  {
1464  for (int r = 0; r < Q - 1; ++r)
1465  {
1466  maparray[idx++] = GetMode(0, 1, r);
1467  }
1468  }
1469  }
1470  }
1471  break;
1472  case 3:
1473  idx = 0;
1474  for (j = 0; j < P; ++j)
1475  {
1476  for (k = 0; k < Q - j; ++k)
1477  {
1478  maparray[idx++] = GetMode(0, j, k);
1479  }
1480  }
1481  break;
1482  default:
1483  ASSERTL0(false, "Element map not available.");
1484  }
1485 }
1486 
1487 void StdTetExp::v_GetElmtTraceToTraceMap(const unsigned int fid,
1488  Array<OneD, unsigned int> &maparray,
1489  Array<OneD, int> &signarray,
1490  Orientation faceOrient, int P, int Q)
1491 {
1492  int nummodesA = 0, nummodesB = 0, i, j, k, idx;
1493 
1495  "Method only implemented for Modified_A BasisType (x "
1496  "direction), Modified_B BasisType (y direction), and "
1497  "Modified_C BasisType(z direction)");
1498 
1499  int nFaceCoeffs = 0;
1500 
1501  switch (fid)
1502  {
1503  case 0:
1504  nummodesA = m_base[0]->GetNumModes();
1505  nummodesB = m_base[1]->GetNumModes();
1506  break;
1507  case 1:
1508  nummodesA = m_base[0]->GetNumModes();
1509  nummodesB = m_base[2]->GetNumModes();
1510  break;
1511  case 2:
1512  case 3:
1513  nummodesA = m_base[1]->GetNumModes();
1514  nummodesB = m_base[2]->GetNumModes();
1515  break;
1516  default:
1517  ASSERTL0(false, "fid must be between 0 and 3");
1518  }
1519 
1520  if (P == -1)
1521  {
1522  P = nummodesA;
1523  Q = nummodesB;
1524  }
1525 
1526  nFaceCoeffs = P * (2 * Q - P + 1) / 2;
1527 
1528  // Allocate the map array and sign array; set sign array to ones (+)
1529  if (maparray.size() != nFaceCoeffs)
1530  {
1531  maparray = Array<OneD, unsigned int>(nFaceCoeffs, 1);
1532  }
1533 
1534  if (signarray.size() != nFaceCoeffs)
1535  {
1536  signarray = Array<OneD, int>(nFaceCoeffs, 1);
1537  }
1538  else
1539  {
1540  fill(signarray.get(), signarray.get() + nFaceCoeffs, 1);
1541  }
1542 
1543  // zero signmap and set maparray to zero if elemental
1544  // modes are not as large as face modesl
1545  idx = 0;
1546  int cnt = 0;
1547  int minPA = min(nummodesA, P);
1548  int minQB = min(nummodesB, Q);
1549 
1550  for (j = 0; j < minPA; ++j)
1551  {
1552  // set maparray
1553  for (k = 0; k < minQB - j; ++k, ++cnt)
1554  {
1555  maparray[idx++] = cnt;
1556  }
1557 
1558  cnt += nummodesB - minQB;
1559 
1560  for (k = nummodesB - j; k < Q - j; ++k)
1561  {
1562  signarray[idx] = 0.0;
1563  maparray[idx++] = maparray[0];
1564  }
1565  }
1566 
1567  for (j = nummodesA; j < P; ++j)
1568  {
1569  for (k = 0; k < Q - j; ++k)
1570  {
1571  signarray[idx] = 0.0;
1572  maparray[idx++] = maparray[0];
1573  }
1574  }
1575 
1576  if (faceOrient == eDir1BwdDir1_Dir2FwdDir2)
1577  {
1578  idx = 0;
1579  for (i = 0; i < P; ++i)
1580  {
1581  for (j = 0; j < Q - i; ++j, idx++)
1582  {
1583  if (i > 1)
1584  {
1585  signarray[idx] = (i % 2 ? -1 : 1);
1586  }
1587  }
1588  }
1589 
1590  swap(maparray[0], maparray[Q]);
1591 
1592  for (i = 1; i < Q - 1; ++i)
1593  {
1594  swap(maparray[i + 1], maparray[Q + i]);
1595  }
1596  }
1597 }
1598 
1599 /**
1600  * Maps interior modes of an edge to the elemental modes.
1601  */
1603  const int eid, Array<OneD, unsigned int> &maparray,
1604  Array<OneD, int> &signarray, const Orientation edgeOrient)
1605 {
1606  int i;
1607  const int P = m_base[0]->GetNumModes();
1608  const int Q = m_base[1]->GetNumModes();
1609  const int R = m_base[2]->GetNumModes();
1610 
1611  const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid) - 2;
1612 
1613  if (maparray.size() != nEdgeIntCoeffs)
1614  {
1615  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1616  }
1617  else
1618  {
1619  fill(maparray.get(), maparray.get() + nEdgeIntCoeffs, 0);
1620  }
1621 
1622  if (signarray.size() != nEdgeIntCoeffs)
1623  {
1624  signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1625  }
1626  else
1627  {
1628  fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1629  }
1630 
1631  switch (eid)
1632  {
1633  case 0:
1634  for (i = 0; i < P - 2; ++i)
1635  {
1636  maparray[i] = GetMode(i + 2, 0, 0);
1637  }
1638  if (edgeOrient == eBackwards)
1639  {
1640  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1641  {
1642  signarray[i] = -1;
1643  }
1644  }
1645  break;
1646  case 1:
1647  for (i = 0; i < Q - 2; ++i)
1648  {
1649  maparray[i] = GetMode(1, i + 1, 0);
1650  }
1651  if (edgeOrient == eBackwards)
1652  {
1653  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1654  {
1655  signarray[i] = -1;
1656  }
1657  }
1658  break;
1659  case 2:
1660  for (i = 0; i < Q - 2; ++i)
1661  {
1662  maparray[i] = GetMode(0, i + 2, 0);
1663  }
1664  if (edgeOrient == eBackwards)
1665  {
1666  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1667  {
1668  signarray[i] = -1;
1669  }
1670  }
1671  break;
1672  case 3:
1673  for (i = 0; i < R - 2; ++i)
1674  {
1675  maparray[i] = GetMode(0, 0, i + 2);
1676  }
1677  if (edgeOrient == eBackwards)
1678  {
1679  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1680  {
1681  signarray[i] = -1;
1682  }
1683  }
1684  break;
1685  case 4:
1686  for (i = 0; i < R - 2; ++i)
1687  {
1688  maparray[i] = GetMode(1, 0, i + 1);
1689  }
1690  if (edgeOrient == eBackwards)
1691  {
1692  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1693  {
1694  signarray[i] = -1;
1695  }
1696  }
1697  break;
1698  case 5:
1699  for (i = 0; i < R - 2; ++i)
1700  {
1701  maparray[i] = GetMode(0, 1, i + 1);
1702  }
1703  if (edgeOrient == eBackwards)
1704  {
1705  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1706  {
1707  signarray[i] = -1;
1708  }
1709  }
1710  break;
1711  default:
1712  ASSERTL0(false, "Edge not defined.");
1713  break;
1714  }
1715 }
1716 
1718  const int fid, Array<OneD, unsigned int> &maparray,
1719  Array<OneD, int> &signarray, const Orientation faceOrient)
1720 {
1721  int i, j, idx, k;
1722  const int P = m_base[0]->GetNumModes();
1723  const int Q = m_base[1]->GetNumModes();
1724  const int R = m_base[2]->GetNumModes();
1725 
1726  const int nFaceIntCoeffs = v_GetTraceIntNcoeffs(fid);
1727 
1728  if (maparray.size() != nFaceIntCoeffs)
1729  {
1730  maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
1731  }
1732 
1733  if (signarray.size() != nFaceIntCoeffs)
1734  {
1735  signarray = Array<OneD, int>(nFaceIntCoeffs, 1);
1736  }
1737  else
1738  {
1739  fill(signarray.get(), signarray.get() + nFaceIntCoeffs, 1);
1740  }
1741 
1742  switch (fid)
1743  {
1744  case 0:
1745  idx = 0;
1746  for (i = 2; i < P - 1; ++i)
1747  {
1748  for (j = 1; j < Q - i; ++j)
1749  {
1750  if ((int)faceOrient == 7)
1751  {
1752  signarray[idx] = (i % 2 ? -1 : 1);
1753  }
1754  maparray[idx++] = GetMode(i, j, 0);
1755  }
1756  }
1757  break;
1758  case 1:
1759  idx = 0;
1760  for (i = 2; i < P; ++i)
1761  {
1762  for (k = 1; k < R - i; ++k)
1763  {
1764  if ((int)faceOrient == 7)
1765  {
1766  signarray[idx] = (i % 2 ? -1 : 1);
1767  }
1768  maparray[idx++] = GetMode(i, 0, k);
1769  }
1770  }
1771  break;
1772  case 2:
1773  idx = 0;
1774  for (j = 1; j < Q - 2; ++j)
1775  {
1776  for (k = 1; k < R - 1 - j; ++k)
1777  {
1778  if ((int)faceOrient == 7)
1779  {
1780  signarray[idx] = ((j + 1) % 2 ? -1 : 1);
1781  }
1782  maparray[idx++] = GetMode(1, j, k);
1783  }
1784  }
1785  break;
1786  case 3:
1787  idx = 0;
1788  for (j = 2; j < Q - 1; ++j)
1789  {
1790  for (k = 1; k < R - j; ++k)
1791  {
1792  if ((int)faceOrient == 7)
1793  {
1794  signarray[idx] = (j % 2 ? -1 : 1);
1795  }
1796  maparray[idx++] = GetMode(0, j, k);
1797  }
1798  }
1799  break;
1800  default:
1801  ASSERTL0(false, "Face interior map not available.");
1802  break;
1803  }
1804 }
1805 //---------------------------------------
1806 // Wrapper functions
1807 //---------------------------------------
1809 {
1810 
1811  MatrixType mtype = mkey.GetMatrixType();
1812 
1813  DNekMatSharedPtr Mat;
1814 
1815  switch (mtype)
1816  {
1818  {
1819  int nq0 = m_base[0]->GetNumPoints();
1820  int nq1 = m_base[1]->GetNumPoints();
1821  int nq2 = m_base[2]->GetNumPoints();
1822  int nq;
1823 
1824  // take definition from key
1825  if (mkey.ConstFactorExists(eFactorConst))
1826  {
1827  nq = (int)mkey.GetConstFactor(eFactorConst);
1828  }
1829  else
1830  {
1831  nq = max(nq0, max(nq1, nq2));
1832  }
1833 
1834  int neq =
1837  Array<OneD, NekDouble> coll(3);
1839  Array<OneD, NekDouble> tmp(nq0);
1840 
1841  Mat =
1842  MemoryManager<DNekMat>::AllocateSharedPtr(neq, nq0 * nq1 * nq2);
1843  int cnt = 0;
1844 
1845  for (int i = 0; i < nq; ++i)
1846  {
1847  for (int j = 0; j < nq - i; ++j)
1848  {
1849  for (int k = 0; k < nq - i - j; ++k, ++cnt)
1850  {
1851  coords[cnt] = Array<OneD, NekDouble>(3);
1852  coords[cnt][0] = -1.0 + 2 * k / (NekDouble)(nq - 1);
1853  coords[cnt][1] = -1.0 + 2 * j / (NekDouble)(nq - 1);
1854  coords[cnt][2] = -1.0 + 2 * i / (NekDouble)(nq - 1);
1855  }
1856  }
1857  }
1858 
1859  for (int i = 0; i < neq; ++i)
1860  {
1861  LocCoordToLocCollapsed(coords[i], coll);
1862 
1863  I[0] = m_base[0]->GetI(coll);
1864  I[1] = m_base[1]->GetI(coll + 1);
1865  I[2] = m_base[2]->GetI(coll + 2);
1866 
1867  // interpolate first coordinate direction
1868  NekDouble fac;
1869  for (int k = 0; k < nq2; ++k)
1870  {
1871  for (int j = 0; j < nq1; ++j)
1872  {
1873 
1874  fac = (I[1]->GetPtr())[j] * (I[2]->GetPtr())[k];
1875  Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1876 
1877  Vmath::Vcopy(nq0, &tmp[0], 1,
1878  Mat->GetRawPtr() + k * nq0 * nq1 * neq +
1879  j * nq0 * neq + i,
1880  neq);
1881  }
1882  }
1883  }
1884  }
1885  break;
1886  default:
1887  {
1889  }
1890  break;
1891  }
1892 
1893  return Mat;
1894 }
1895 
1897 {
1898  return v_GenMatrix(mkey);
1899 }
1900 
1901 //---------------------------------------
1902 // Private helper functions
1903 //---------------------------------------
1904 
1905 /**
1906  * @brief Compute the mode number in the expansion for a particular
1907  * tensorial combination.
1908  *
1909  * Modes are numbered with the r index travelling fastest, followed by
1910  * q and then p, and each q-r plane is of size
1911  * (Q+1)*(Q+2)/2+max(0,R-Q-p)*Q. For example, when P=2, Q=3 and R=4
1912  * the indexing inside each q-r plane (with r increasing upwards and q
1913  * to the right) is:
1914  *
1915  * p = 0: p = 1: p = 2:
1916  * ----------------------------------
1917  * 4
1918  * 3 8 17
1919  * 2 7 11 16 20 26
1920  * 1 6 10 13 15 19 22 25 28
1921  * 0 5 9 12 14 18 21 23 24 27 29
1922  *
1923  * Note that in this element, we must have that \f$ P \leq Q \leq
1924  * R\f$.
1925  */
1926 int StdTetExp::GetMode(const int I, const int J, const int K)
1927 {
1928  const int Q = m_base[1]->GetNumModes();
1929  const int R = m_base[2]->GetNumModes();
1930 
1931  int i, j, q_hat, k_hat;
1932  int cnt = 0;
1933 
1934  // Traverse to q-r plane number I
1935  for (i = 0; i < I; ++i)
1936  {
1937  // Size of triangle part
1938  q_hat = min(Q, R - i);
1939  // Size of rectangle part
1940  k_hat = max(R - Q - i, 0);
1941  cnt += q_hat * (q_hat + 1) / 2 + k_hat * Q;
1942  }
1943 
1944  // Traverse to q column J
1945  q_hat = R - I;
1946  for (j = 0; j < J; ++j)
1947  {
1948  cnt += q_hat;
1949  q_hat--;
1950  }
1951 
1952  // Traverse up stacks to K
1953  cnt += K;
1954 
1955  return cnt;
1956 }
1957 
1959  const Array<OneD, const NekDouble> &inarray,
1960  Array<OneD, NekDouble> &outarray)
1961 {
1962  int i, j;
1963 
1964  int nquad0 = m_base[0]->GetNumPoints();
1965  int nquad1 = m_base[1]->GetNumPoints();
1966  int nquad2 = m_base[2]->GetNumPoints();
1967 
1968  const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
1969  const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
1970  const Array<OneD, const NekDouble> &w2 = m_base[2]->GetW();
1971 
1972  const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
1973  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
1974 
1975  // multiply by integration constants
1976  for (i = 0; i < nquad1 * nquad2; ++i)
1977  {
1978  Vmath::Vmul(nquad0, (NekDouble *)&inarray[0] + i * nquad0, 1, w0.get(),
1979  1, &outarray[0] + i * nquad0, 1);
1980  }
1981 
1982  switch (m_base[1]->GetPointsType())
1983  {
1984  // (1,0) Jacobi Inner product.
1985  case LibUtilities::eGaussRadauMAlpha1Beta0:
1986  for (j = 0; j < nquad2; ++j)
1987  {
1988  for (i = 0; i < nquad1; ++i)
1989  {
1990  Blas::Dscal(nquad0, 0.5 * w1[i],
1991  &outarray[0] + i * nquad0 + j * nquad0 * nquad1,
1992  1);
1993  }
1994  }
1995  break;
1996 
1997  default:
1998  for (j = 0; j < nquad2; ++j)
1999  {
2000  for (i = 0; i < nquad1; ++i)
2001  {
2002  Blas::Dscal(nquad0, 0.5 * (1 - z1[i]) * w1[i],
2003  &outarray[0] + i * nquad0 + j * nquad0 * nquad1,
2004  1);
2005  }
2006  }
2007  break;
2008  }
2009 
2010  switch (m_base[2]->GetPointsType())
2011  {
2012  // (2,0) Jacobi inner product.
2013  case LibUtilities::eGaussRadauMAlpha2Beta0:
2014  for (i = 0; i < nquad2; ++i)
2015  {
2016  Blas::Dscal(nquad0 * nquad1, 0.25 * w2[i],
2017  &outarray[0] + i * nquad0 * nquad1, 1);
2018  }
2019  break;
2020  // (1,0) Jacobi inner product.
2021  case LibUtilities::eGaussRadauMAlpha1Beta0:
2022  for (i = 0; i < nquad2; ++i)
2023  {
2024  Blas::Dscal(nquad0 * nquad1, 0.25 * (1 - z2[i]) * w2[i],
2025  &outarray[0] + i * nquad0 * nquad1, 1);
2026  }
2027  break;
2028  default:
2029  for (i = 0; i < nquad2; ++i)
2030  {
2031  Blas::Dscal(nquad0 * nquad1,
2032  0.25 * (1 - z2[i]) * (1 - z2[i]) * w2[i],
2033  &outarray[0] + i * nquad0 * nquad1, 1);
2034  }
2035  break;
2036  }
2037 }
2038 
2040  const StdMatrixKey &mkey)
2041 {
2042  // To do : 1) add a test to ensure 0 \leq SvvCutoff \leq 1.
2043  // 2) check if the transfer function needs an analytical
2044  // Fourier transform.
2045  // 3) if it doesn't : find a transfer function that renders
2046  // the if( cutoff_a ...) useless to reduce computational
2047  // cost.
2048  // 4) add SVVDiffCoef to both models!!
2049 
2050  int qa = m_base[0]->GetNumPoints();
2051  int qb = m_base[1]->GetNumPoints();
2052  int qc = m_base[2]->GetNumPoints();
2053  int nmodes_a = m_base[0]->GetNumModes();
2054  int nmodes_b = m_base[1]->GetNumModes();
2055  int nmodes_c = m_base[2]->GetNumModes();
2056 
2057  // Declare orthogonal basis.
2061 
2065 
2066  StdTetExp OrthoExp(Ba, Bb, Bc);
2067 
2068  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2069  int i, j, k, cnt = 0;
2070 
2071  // project onto physical space.
2072  OrthoExp.FwdTrans(array, orthocoeffs);
2073 
2075  {
2076  // Rodrigo's power kernel
2078  NekDouble SvvDiffCoeff =
2081 
2082  for (int i = 0; i < nmodes_a; ++i)
2083  {
2084  for (int j = 0; j < nmodes_b - j; ++j)
2085  {
2086  NekDouble fac1 = std::max(
2087  pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2088  pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2089 
2090  for (int k = 0; k < nmodes_c - i - j; ++k)
2091  {
2092  NekDouble fac =
2093  std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2094  cutoff * nmodes_c));
2095 
2096  orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2097  cnt++;
2098  }
2099  }
2100  }
2101  }
2102  else if (mkey.ConstFactorExists(
2103  eFactorSVVDGKerDiffCoeff)) // Rodrigo/Mansoor's DG Kernel
2104  {
2105  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff) *
2107 
2108  int max_abc = max(nmodes_a - kSVVDGFiltermodesmin,
2109  nmodes_b - kSVVDGFiltermodesmin);
2110  max_abc = max(max_abc, nmodes_c - kSVVDGFiltermodesmin);
2111  // clamp max_abc
2112  max_abc = max(max_abc, 0);
2113  max_abc = min(max_abc, kSVVDGFiltermodesmax - kSVVDGFiltermodesmin);
2114 
2115  for (int i = 0; i < nmodes_a; ++i)
2116  {
2117  for (int j = 0; j < nmodes_b - j; ++j)
2118  {
2119  int maxij = max(i, j);
2120 
2121  for (int k = 0; k < nmodes_c - i - j; ++k)
2122  {
2123  int maxijk = max(maxij, k);
2124  maxijk = min(maxijk, kSVVDGFiltermodesmax - 1);
2125 
2126  orthocoeffs[cnt] *=
2127  SvvDiffCoeff * kSVVDGFilter[max_abc][maxijk];
2128  cnt++;
2129  }
2130  }
2131  }
2132  }
2133  else
2134  {
2135 
2136  // SVV filter paramaters (how much added diffusion
2137  // relative to physical one and fraction of modes from
2138  // which you start applying this added diffusion)
2139 
2140  NekDouble SvvDiffCoeff =
2142  NekDouble SVVCutOff =
2144 
2145  // Defining the cut of mode
2146  int cutoff_a = (int)(SVVCutOff * nmodes_a);
2147  int cutoff_b = (int)(SVVCutOff * nmodes_b);
2148  int cutoff_c = (int)(SVVCutOff * nmodes_c);
2149  int nmodes = min(min(nmodes_a, nmodes_b), nmodes_c);
2150  NekDouble cutoff = min(min(cutoff_a, cutoff_b), cutoff_c);
2151  NekDouble epsilon = 1;
2152 
2153  //------"New" Version August 22nd '13--------------------
2154  for (i = 0; i < nmodes_a; ++i)
2155  {
2156  for (j = 0; j < nmodes_b - i; ++j)
2157  {
2158  for (k = 0; k < nmodes_c - i - j; ++k)
2159  {
2160  if (i + j + k >= cutoff)
2161  {
2162  orthocoeffs[cnt] *= ((SvvDiffCoeff)*exp(
2163  -(i + j + k - nmodes) * (i + j + k - nmodes) /
2164  ((NekDouble)((i + j + k - cutoff + epsilon) *
2165  (i + j + k - cutoff + epsilon)))));
2166  }
2167  else
2168  {
2169  orthocoeffs[cnt] *= 0.0;
2170  }
2171  cnt++;
2172  }
2173  }
2174  }
2175  }
2176 
2177  // backward transform to physical space
2178  OrthoExp.BwdTrans(orthocoeffs, array);
2179 }
2180 
2182  const Array<OneD, const NekDouble> &inarray,
2183  Array<OneD, NekDouble> &outarray)
2184 {
2185  int nquad0 = m_base[0]->GetNumPoints();
2186  int nquad1 = m_base[1]->GetNumPoints();
2187  int nquad2 = m_base[2]->GetNumPoints();
2188  int nqtot = nquad0 * nquad1 * nquad2;
2189  int nmodes0 = m_base[0]->GetNumModes();
2190  int nmodes1 = m_base[1]->GetNumModes();
2191  int nmodes2 = m_base[2]->GetNumModes();
2192  int numMax = nmodes0;
2193 
2195  Array<OneD, NekDouble> coeff_tmp1(m_ncoeffs, 0.0);
2196  Array<OneD, NekDouble> coeff_tmp2(m_ncoeffs, 0.0);
2197  Array<OneD, NekDouble> phys_tmp(nqtot, 0.0);
2198  Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
2199 
2200  Vmath::Vcopy(m_ncoeffs, inarray, 1, coeff_tmp2, 1);
2201 
2202  const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
2203  const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
2204  const LibUtilities::PointsKey Pkey2 = m_base[2]->GetPointsKey();
2205 
2206  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
2207  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_B, nmodes1, Pkey1);
2208  LibUtilities::BasisKey bortho2(LibUtilities::eOrtho_C, nmodes2, Pkey2);
2209 
2210  Vmath::Zero(m_ncoeffs, coeff_tmp2, 1);
2211 
2212  StdRegions::StdTetExpSharedPtr OrthoTetExp;
2214  bortho0, bortho1, bortho2);
2215 
2216  BwdTrans(inarray, phys_tmp);
2217  OrthoTetExp->FwdTrans(phys_tmp, coeff);
2218 
2219  Vmath::Zero(m_ncoeffs, outarray, 1);
2220 
2221  // filtering
2222  int cnt = 0;
2223  for (int u = 0; u < numMin; ++u)
2224  {
2225  for (int i = 0; i < numMin - u; ++i)
2226  {
2227  Vmath::Vcopy(numMin - u - i, tmp = coeff + cnt, 1,
2228  tmp2 = coeff_tmp1 + cnt, 1);
2229  cnt += numMax - u - i;
2230  }
2231  for (int i = numMin; i < numMax - u; ++i)
2232  {
2233  cnt += numMax - u - i;
2234  }
2235  }
2236 
2237  OrthoTetExp->BwdTrans(coeff_tmp1, phys_tmp);
2238  FwdTrans(phys_tmp, outarray);
2239 }
2240 
2242  bool standard)
2243 {
2244  boost::ignore_unused(standard);
2245 
2246  int np0 = m_base[0]->GetNumPoints();
2247  int np1 = m_base[1]->GetNumPoints();
2248  int np2 = m_base[2]->GetNumPoints();
2249  int np = max(np0, max(np1, np2));
2250 
2251  conn = Array<OneD, int>(4 * (np - 1) * (np - 1) * (np - 1));
2252 
2253  int row = 0;
2254  int rowp1 = 0;
2255  int plane = 0;
2256  int row1 = 0;
2257  int row1p1 = 0;
2258  int planep1 = 0;
2259  int cnt = 0;
2260  for (int i = 0; i < np - 1; ++i)
2261  {
2262  planep1 += (np - i) * (np - i + 1) / 2;
2263  row = 0; // current plane row offset
2264  rowp1 = 0; // current plane row plus one offset
2265  row1 = 0; // next plane row offset
2266  row1p1 = 0; // nex plane row plus one offset
2267  for (int j = 0; j < np - i - 1; ++j)
2268  {
2269  rowp1 += np - i - j;
2270  row1p1 += np - i - j - 1;
2271  for (int k = 0; k < np - i - j - 2; ++k)
2272  {
2273  conn[cnt++] = plane + row + k + 1;
2274  conn[cnt++] = plane + row + k;
2275  conn[cnt++] = plane + rowp1 + k;
2276  conn[cnt++] = planep1 + row1 + k;
2277 
2278  conn[cnt++] = plane + row + k + 1;
2279  conn[cnt++] = plane + rowp1 + k + 1;
2280  conn[cnt++] = planep1 + row1 + k + 1;
2281  conn[cnt++] = planep1 + row1 + k;
2282 
2283  conn[cnt++] = plane + rowp1 + k + 1;
2284  conn[cnt++] = plane + row + k + 1;
2285  conn[cnt++] = plane + rowp1 + k;
2286  conn[cnt++] = planep1 + row1 + k;
2287 
2288  conn[cnt++] = planep1 + row1 + k;
2289  conn[cnt++] = planep1 + row1p1 + k;
2290  conn[cnt++] = plane + rowp1 + k;
2291  conn[cnt++] = plane + rowp1 + k + 1;
2292 
2293  conn[cnt++] = planep1 + row1 + k;
2294  conn[cnt++] = planep1 + row1p1 + k;
2295  conn[cnt++] = planep1 + row1 + k + 1;
2296  conn[cnt++] = plane + rowp1 + k + 1;
2297 
2298  if (k < np - i - j - 3)
2299  {
2300  conn[cnt++] = plane + rowp1 + k + 1;
2301  conn[cnt++] = planep1 + row1p1 + k + 1;
2302  conn[cnt++] = planep1 + row1 + k + 1;
2303  conn[cnt++] = planep1 + row1p1 + k;
2304  }
2305  }
2306 
2307  conn[cnt++] = plane + row + np - i - j - 1;
2308  conn[cnt++] = plane + row + np - i - j - 2;
2309  conn[cnt++] = plane + rowp1 + np - i - j - 2;
2310  conn[cnt++] = planep1 + row1 + np - i - j - 2;
2311 
2312  row += np - i - j;
2313  row1 += np - i - j - 1;
2314  }
2315  plane += (np - i) * (np - i + 1) / 2;
2316  }
2317 }
2318 
2319 } // namespace StdRegions
2320 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#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
Defines a specification for a set of points.
Definition: Points.h:59
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points.
The base class for all shapes.
Definition: StdExpansion.h:71
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:162
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
This function gets the shared point to basis.
Definition: StdExpansion.h:106
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
NekDouble 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.
Definition: StdExpansion.h:918
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
Definition: StdExpansion.h:430
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:211
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:224
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:729
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:175
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:848
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:85
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:126
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:135
virtual int v_GetNtraces() const override
Definition: StdTetExp.cpp:984
virtual int v_GetTraceNcoeffs(const int i) const override
Definition: StdTetExp.cpp:1034
void v_GetTraceNumModes(const int fid, int &numModes0, int &numModes1, Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
Definition: StdTetExp.cpp:937
NekDouble v_PhysEvaluate(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
Definition: StdTetExp.cpp:873
virtual int v_GetTraceIntNcoeffs(const int i) const override
Definition: StdTetExp.cpp:1060
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:824
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:432
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
Definition: StdTetExp.cpp:1301
virtual int v_GetNedges() const override
Definition: StdTetExp.cpp:979
int GetMode(const int i, const int j, const int k)
Compute the mode number in the expansion for a particular tensorial combination.
Definition: StdTetExp.cpp:1926
virtual int v_GetEdgeNcoeffs(const int i) const override
Definition: StdTetExp.cpp:1099
virtual bool v_IsBoundaryInteriorExpansion() const override
Definition: StdTetExp.cpp:1210
virtual LibUtilities::PointsKey v_GetTracePointsKey(const int i, const int j) const override
Definition: StdTetExp.cpp:1120
virtual const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int k) const override
Definition: StdTetExp.cpp:1151
void v_GetEdgeInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
Definition: StdTetExp.cpp:1602
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:1958
virtual int v_GetTraceNumPoints(const int i) const override
Definition: StdTetExp.cpp:1081
NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) final override
Definition: StdTetExp.cpp:831
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_BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2) override
Definition: StdTetExp.cpp:338
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
Definition: StdTetExp.cpp:2039
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
Definition: StdTetExp.cpp:1896
virtual void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
Definition: StdTetExp.cpp:784
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
Definition: StdTetExp.cpp:1340
virtual void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true) override
Definition: StdTetExp.cpp:2241
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:482
virtual int v_NumBndryCoeffs() const override
Definition: StdTetExp.cpp:994
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdTetExp.cpp:509
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:282
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:620
virtual ~StdTetExp() override
Definition: StdTetExp.cpp:75
void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:2181
virtual void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
Definition: StdTetExp.cpp:242
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z) override
Definition: StdTetExp.cpp:1180
virtual void v_GetTraceCoeffMap(const unsigned int fid, Array< OneD, unsigned int > &maparray) override
Definition: StdTetExp.cpp:1395
void v_GetTraceInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
Definition: StdTetExp.cpp:1717
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:309
void v_GetElmtTraceToTraceMap(const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1) override
Definition: StdTetExp.cpp:1487
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
Definition: StdTetExp.cpp:1220
virtual int v_GetNverts() const override
Definition: StdTetExp.cpp:974
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dx, Array< OneD, NekDouble > &out_dy, Array< OneD, NekDouble > &out_dz) override
Calculate the derivative of the physical points.
Definition: StdTetExp.cpp:99
virtual LibUtilities::ShapeType v_DetShapeType() const override
Definition: StdTetExp.cpp:989
virtual void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
Definition: StdTetExp.cpp:816
LibUtilities::ShapeType DetShapeType() const
Definition: StdTetExp.h:64
virtual int v_NumDGBndryCoeffs() const override
Definition: StdTetExp.cpp:1013
virtual void v_IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2) override
Definition: StdTetExp.cpp:539
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) override
Definition: StdTetExp.cpp:1140
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
Definition: StdTetExp.cpp:1808
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:246
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:168
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:182
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:368
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:154
int getNumberOfBndCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:215
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:192
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:52
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:48
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
static const NekDouble kNekZeroTol
LibUtilities::BasisKey EvaluateTriFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
std::shared_ptr< StdTetExp > StdTetExpSharedPtr
Definition: StdTetExp.h:255
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:469
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:470
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:472
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:359
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:384
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294