Nektar++
PyrExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PyrExp.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: PyrExp routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
36#include <LocalRegions/PyrExp.h>
37
38using namespace std;
39
40namespace Nektar
41{
42namespace LocalRegions
43{
44
46 const LibUtilities::BasisKey &Bb,
47 const LibUtilities::BasisKey &Bc,
49 : StdExpansion(LibUtilities::StdPyrData::getNumberOfCoefficients(
50 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
51 3, Ba, Bb, Bc),
52 StdExpansion3D(LibUtilities::StdPyrData::getNumberOfCoefficients(
53 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
54 Ba, Bb, Bc),
55 StdPyrExp(Ba, Bb, Bc), Expansion(geom), Expansion3D(geom),
56 m_matrixManager(
57 std::bind(&Expansion3D::CreateMatrix, this, std::placeholders::_1),
58 std::string("PyrExpMatrix")),
59 m_staticCondMatrixManager(std::bind(&Expansion::CreateStaticCondMatrix,
60 this, std::placeholders::_1),
61 std::string("PyrExpStaticCondMatrix"))
62{
63}
64
66 : StdExpansion(T), StdExpansion3D(T), StdPyrExp(T), Expansion(T),
67 Expansion3D(T), m_matrixManager(T.m_matrixManager),
68 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
69{
70}
71
72//----------------------------
73// Integration Methods
74//----------------------------
75
76/**
77 * \brief Integrate the physical point list \a inarray over pyramidic
78 * region and return the value.
79 *
80 * Inputs:\n
81 *
82 * - \a inarray: definition of function to be returned at quadrature
83 * point of expansion.
84 *
85 * Outputs:\n
86 *
87 * - returns \f$\int^1_{-1}\int^1_{-1}\int^1_{-1} u(\bar \eta_1,
88 * \eta_2, \eta_3) J[i,j,k] d \bar \eta_1 d \eta_2 d \eta_3\f$ \n \f$=
89 * \sum_{i=0}^{Q_1 - 1} \sum_{j=0}^{Q_2 - 1} \sum_{k=0}^{Q_3 - 1}
90 * u(\bar \eta_{1i}^{0,0}, \eta_{2j}^{0,0},\eta_{3k}^{2,0})w_{i}^{0,0}
91 * w_{j}^{0,0} \hat w_{k}^{2,0} \f$ \n where \f$inarray[i,j, k] =
92 * u(\bar \eta_{1i},\eta_{2j}, \eta_{3k}) \f$, \n \f$\hat w_{k}^{2,0}
93 * = \frac {w^{2,0}} {2} \f$ \n and \f$ J[i,j,k] \f$ is the Jacobian
94 * evaluated at the quadrature point.
95 */
97{
98 int nquad0 = m_base[0]->GetNumPoints();
99 int nquad1 = m_base[1]->GetNumPoints();
100 int nquad2 = m_base[2]->GetNumPoints();
102 Array<OneD, NekDouble> tmp(nquad0 * nquad1 * nquad2);
103
104 // multiply inarray with Jacobian
105 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
106 {
107 Vmath::Vmul(nquad0 * nquad1 * nquad2, &jac[0], 1,
108 (NekDouble *)&inarray[0], 1, &tmp[0], 1);
109 }
110 else
111 {
112 Vmath::Smul(nquad0 * nquad1 * nquad2, (NekDouble)jac[0],
113 (NekDouble *)&inarray[0], 1, &tmp[0], 1);
114 }
115
116 // call StdPyrExp version;
117 return StdPyrExp::v_Integral(tmp);
118}
119
120//----------------------------
121// Differentiation Methods
122//----------------------------
123
128{
129 int nquad0 = m_base[0]->GetNumPoints();
130 int nquad1 = m_base[1]->GetNumPoints();
131 int nquad2 = m_base[2]->GetNumPoints();
133 m_metricinfo->GetDerivFactors(GetPointsKeys());
134 Array<OneD, NekDouble> diff0(nquad0 * nquad1 * nquad2);
135 Array<OneD, NekDouble> diff1(nquad0 * nquad1 * nquad2);
136 Array<OneD, NekDouble> diff2(nquad0 * nquad1 * nquad2);
137
138 StdPyrExp::v_PhysDeriv(inarray, diff0, diff1, diff2);
139
140 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
141 {
142 if (out_d0.size())
143 {
144 Vmath::Vmul(nquad0 * nquad1 * nquad2, &gmat[0][0], 1, &diff0[0], 1,
145 &out_d0[0], 1);
146 Vmath::Vvtvp(nquad0 * nquad1 * nquad2, &gmat[1][0], 1, &diff1[0], 1,
147 &out_d0[0], 1, &out_d0[0], 1);
148 Vmath::Vvtvp(nquad0 * nquad1 * nquad2, &gmat[2][0], 1, &diff2[0], 1,
149 &out_d0[0], 1, &out_d0[0], 1);
150 }
151
152 if (out_d1.size())
153 {
154 Vmath::Vmul(nquad0 * nquad1 * nquad2, &gmat[3][0], 1, &diff0[0], 1,
155 &out_d1[0], 1);
156 Vmath::Vvtvp(nquad0 * nquad1 * nquad2, &gmat[4][0], 1, &diff1[0], 1,
157 &out_d1[0], 1, &out_d1[0], 1);
158 Vmath::Vvtvp(nquad0 * nquad1 * nquad2, &gmat[5][0], 1, &diff2[0], 1,
159 &out_d1[0], 1, &out_d1[0], 1);
160 }
161
162 if (out_d2.size())
163 {
164 Vmath::Vmul(nquad0 * nquad1 * nquad2, &gmat[6][0], 1, &diff0[0], 1,
165 &out_d2[0], 1);
166 Vmath::Vvtvp(nquad0 * nquad1 * nquad2, &gmat[7][0], 1, &diff1[0], 1,
167 &out_d2[0], 1, &out_d2[0], 1);
168 Vmath::Vvtvp(nquad0 * nquad1 * nquad2, &gmat[8][0], 1, &diff2[0], 1,
169 &out_d2[0], 1, &out_d2[0], 1);
170 }
171 }
172 else // regular geometry
173 {
174 if (out_d0.size())
175 {
176 Vmath::Smul(nquad0 * nquad1 * nquad2, gmat[0][0], &diff0[0], 1,
177 &out_d0[0], 1);
178 Blas::Daxpy(nquad0 * nquad1 * nquad2, gmat[1][0], &diff1[0], 1,
179 &out_d0[0], 1);
180 Blas::Daxpy(nquad0 * nquad1 * nquad2, gmat[2][0], &diff2[0], 1,
181 &out_d0[0], 1);
182 }
183
184 if (out_d1.size())
185 {
186 Vmath::Smul(nquad0 * nquad1 * nquad2, gmat[3][0], &diff0[0], 1,
187 &out_d1[0], 1);
188 Blas::Daxpy(nquad0 * nquad1 * nquad2, gmat[4][0], &diff1[0], 1,
189 &out_d1[0], 1);
190 Blas::Daxpy(nquad0 * nquad1 * nquad2, gmat[5][0], &diff2[0], 1,
191 &out_d1[0], 1);
192 }
193
194 if (out_d2.size())
195 {
196 Vmath::Smul(nquad0 * nquad1 * nquad2, gmat[6][0], &diff0[0], 1,
197 &out_d2[0], 1);
198 Blas::Daxpy(nquad0 * nquad1 * nquad2, gmat[7][0], &diff1[0], 1,
199 &out_d2[0], 1);
200 Blas::Daxpy(nquad0 * nquad1 * nquad2, gmat[8][0], &diff2[0], 1,
201 &out_d2[0], 1);
202 }
203 }
204}
205
206//---------------------------------------
207// Transforms
208//---------------------------------------
209
210/**
211 * \brief Forward transform from physical quadrature space stored in
212 * \a inarray and evaluate the expansion coefficients and store in \a
213 * (this)->m_coeffs
214 *
215 * Inputs:\n
216 *
217 * - \a inarray: array of physical quadrature points to be transformed
218 *
219 * Outputs:\n
220 *
221 * - (this)->_coeffs: updated array of expansion coefficients.
222 */
224 Array<OneD, NekDouble> &outarray)
225{
226 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
227 m_base[2]->Collocation())
228 {
229 Vmath::Vcopy(GetNcoeffs(), &inarray[0], 1, &outarray[0], 1);
230 }
231 else
232 {
233 v_IProductWRTBase(inarray, outarray);
234
235 // get Mass matrix inverse
237 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
238
239 // copy inarray in case inarray == outarray
240 DNekVec in(m_ncoeffs, outarray);
241 DNekVec out(m_ncoeffs, outarray, eWrapper);
242
243 out = (*matsys) * in;
244 }
245}
246
247//---------------------------------------
248// Inner product functions
249//---------------------------------------
250
251/**
252 * \brief Calculate the inner product of inarray with respect to the
253 * basis B=base0*base1*base2 and put into outarray:
254 *
255 * \f$ \begin{array}{rcl} I_{pqr} = (\phi_{pqr}, u)_{\delta} & = &
256 * \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \sum_{k=0}^{nq_2} \psi_{p}^{a}
257 * (\bar \eta_{1i}) \psi_{q}^{a} (\eta_{2j}) \psi_{pqr}^{c}
258 * (\eta_{3k}) w_i w_j w_k u(\bar \eta_{1,i} \eta_{2,j} \eta_{3,k})
259 * J_{i,j,k}\\ & = & \sum_{i=0}^{nq_0} \psi_p^a(\bar \eta_{1,i})
260 * \sum_{j=0}^{nq_1} \psi_{q}^a(\eta_{2,j}) \sum_{k=0}^{nq_2}
261 * \psi_{pqr}^c u(\bar \eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k}
262 * \end{array} \f$ \n
263 *
264 * where
265 *
266 * \f$\phi_{pqr} (\xi_1 , \xi_2 , \xi_3) = \psi_p^a (\bar \eta_1)
267 * \psi_{q}^a (\eta_2) \psi_{pqr}^c (\eta_3) \f$ \n
268 *
269 * which can be implemented as \n \f$f_{pqr} (\xi_{3k}) =
270 * \sum_{k=0}^{nq_3} \psi_{pqr}^c u(\bar
271 * \eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} = {\bf B_3 U} \f$ \n \f$
272 * g_{pq} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{q}^a (\xi_{2j}) f_{pqr}
273 * (\xi_{3k}) = {\bf B_2 F} \f$ \n \f$ (\phi_{pqr}, u)_{\delta} =
274 * \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{pq} (\xi_{3k}) = {\bf
275 * B_1 G} \f$
276 */
277
279 Array<OneD, NekDouble> &outarray)
280{
281 v_IProductWRTBase_SumFac(inarray, outarray);
282}
283
285 const Array<OneD, const NekDouble> &inarray,
286 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
287{
288 const int nquad0 = m_base[0]->GetNumPoints();
289 const int nquad1 = m_base[1]->GetNumPoints();
290 const int nquad2 = m_base[2]->GetNumPoints();
291 const int order0 = m_base[0]->GetNumModes();
292 const int order1 = m_base[1]->GetNumModes();
293
294 Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
295
296 if (multiplybyweights)
297 {
298 Array<OneD, NekDouble> tmp(nquad0 * nquad1 * nquad2);
299
300 MultiplyByQuadratureMetric(inarray, tmp);
301
303 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
304 tmp, outarray, wsp, true, true, true);
305 }
306 else
307 {
309 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
310 inarray, outarray, wsp, true, true, true);
311 }
312}
313
314/**
315 * @brief Calculates the inner product \f$ I_{pqr} = (u,
316 * \partial_{x_i} \phi_{pqr}) \f$.
317 *
318 * The derivative of the basis functions is performed using the chain
319 * rule in order to incorporate the geometric factors. Assuming that
320 * the basis functions are a tensor product
321 * \f$\phi_{pqr}(\eta_1,\eta_2,\eta_3) =
322 * \phi_1(\eta_1)\phi_2(\eta_2)\phi_3(\eta_3)\f$, this yields the
323 * result
324 *
325 * \f[
326 * I_{pqr} = \sum_{j=1}^3 \left(u, \frac{\partial u}{\partial \eta_j}
327 * \frac{\partial \eta_j}{\partial x_i}\right)
328 * \f]
329 *
330 * In the pyramid element, we must also incorporate a second set
331 * of geometric factors which incorporate the collapsed co-ordinate
332 * system, so that
333 *
334 * \f[ \frac{\partial\eta_j}{\partial x_i} = \sum_{k=1}^3
335 * \frac{\partial\eta_j}{\partial\xi_k}\frac{\partial\xi_k}{\partial
336 * x_i} \f]
337 *
338 * These derivatives can be found on p152 of Sherwin & Karniadakis.
339 *
340 * @param dir Direction in which to take the derivative.
341 * @param inarray The function \f$ u \f$.
342 * @param outarray Value of the inner product.
343 */
345 const Array<OneD, const NekDouble> &inarray,
346 Array<OneD, NekDouble> &outarray)
347{
348 v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
349}
350
352 const int dir, const Array<OneD, const NekDouble> &inarray,
353 Array<OneD, NekDouble> &outarray)
354{
355 const int nquad0 = m_base[0]->GetNumPoints();
356 const int nquad1 = m_base[1]->GetNumPoints();
357 const int nquad2 = m_base[2]->GetNumPoints();
358 const int order0 = m_base[0]->GetNumModes();
359 const int order1 = m_base[1]->GetNumModes();
360 const int nqtot = nquad0 * nquad1 * nquad2;
361
362 Array<OneD, NekDouble> tmp1(nqtot);
363 Array<OneD, NekDouble> tmp2(nqtot);
364 Array<OneD, NekDouble> tmp3(nqtot);
365 Array<OneD, NekDouble> tmp4(nqtot);
368 std::max(nqtot, order0 * nquad2 * (nquad1 + order1)));
369
370 MultiplyByQuadratureMetric(inarray, tmp1);
371
373 tmp2D[0] = tmp2;
374 tmp2D[1] = tmp3;
375 tmp2D[2] = tmp4;
376
377 PyrExp::v_AlignVectorToCollapsedDir(dir, tmp1, tmp2D);
378
379 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
380 m_base[2]->GetBdata(), tmp2, outarray, wsp,
381 false, true, true);
382
383 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
384 m_base[2]->GetBdata(), tmp3, tmp6, wsp, true,
385 false, true);
386
387 Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
388
389 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
390 m_base[2]->GetDbdata(), tmp4, tmp6, wsp, true,
391 true, false);
392
393 Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
394}
395
397 const int dir, const Array<OneD, const NekDouble> &inarray,
399{
400 const int nquad0 = m_base[0]->GetNumPoints();
401 const int nquad1 = m_base[1]->GetNumPoints();
402 const int nquad2 = m_base[2]->GetNumPoints();
403 const int order0 = m_base[0]->GetNumModes();
404 const int order1 = m_base[1]->GetNumModes();
405 const int nqtot = nquad0 * nquad1 * nquad2;
406
407 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
408 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
409 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
410
411 Array<OneD, NekDouble> gfac0(nquad0);
412 Array<OneD, NekDouble> gfac1(nquad1);
413 Array<OneD, NekDouble> gfac2(nquad2);
414 Array<OneD, NekDouble> tmp5(nqtot);
416 std::max(nqtot, order0 * nquad2 * (nquad1 + order1)));
417
418 Array<OneD, NekDouble> tmp2 = outarray[0];
419 Array<OneD, NekDouble> tmp3 = outarray[1];
420 Array<OneD, NekDouble> tmp4 = outarray[2];
421
423 m_metricinfo->GetDerivFactors(GetPointsKeys());
424
426 tmp1 = inarray;
427
428 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
429 {
430 Vmath::Vmul(nqtot, &df[3 * dir][0], 1, tmp1.get(), 1, tmp2.get(), 1);
431 Vmath::Vmul(nqtot, &df[3 * dir + 1][0], 1, tmp1.get(), 1, tmp3.get(),
432 1);
433 Vmath::Vmul(nqtot, &df[3 * dir + 2][0], 1, tmp1.get(), 1, tmp4.get(),
434 1);
435 }
436 else
437 {
438 Vmath::Smul(nqtot, df[3 * dir][0], tmp1.get(), 1, tmp2.get(), 1);
439 Vmath::Smul(nqtot, df[3 * dir + 1][0], tmp1.get(), 1, tmp3.get(), 1);
440 Vmath::Smul(nqtot, df[3 * dir + 2][0], tmp1.get(), 1, tmp4.get(), 1);
441 }
442
443 // set up geometric factor: (1+z0)/2
444 for (int i = 0; i < nquad0; ++i)
445 {
446 gfac0[i] = 0.5 * (1 + z0[i]);
447 }
448
449 // set up geometric factor: (1+z1)/2
450 for (int i = 0; i < nquad1; ++i)
451 {
452 gfac1[i] = 0.5 * (1 + z1[i]);
453 }
454
455 // Set up geometric factor: 2/(1-z2)
456 for (int i = 0; i < nquad2; ++i)
457 {
458 gfac2[i] = 2.0 / (1 - z2[i]);
459 }
460
461 const int nq01 = nquad0 * nquad1;
462
463 for (int i = 0; i < nquad2; ++i)
464 {
465 Vmath::Smul(nq01, gfac2[i], &tmp2[0] + i * nq01, 1, &tmp2[0] + i * nq01,
466 1); // 2/(1-z2) for d/dxi_0
467 Vmath::Smul(nq01, gfac2[i], &tmp3[0] + i * nq01, 1, &tmp3[0] + i * nq01,
468 1); // 2/(1-z2) for d/dxi_1
469 Vmath::Smul(nq01, gfac2[i], &tmp4[0] + i * nq01, 1, &tmp5[0] + i * nq01,
470 1); // 2/(1-z2) for d/dxi_2
471 }
472
473 // (1+z0)/(1-z2) for d/d eta_0
474 for (int i = 0; i < nquad1 * nquad2; ++i)
475 {
476 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp5[0] + i * nquad0, 1,
477 &wsp[0] + i * nquad0, 1);
478 }
479
480 Vmath::Vadd(nqtot, &tmp2[0], 1, &wsp[0], 1, &tmp2[0], 1);
481
482 // (1+z1)/(1-z2) for d/d eta_1
483 for (int i = 0; i < nquad1 * nquad2; ++i)
484 {
485 Vmath::Smul(nquad0, gfac1[i % nquad1], &tmp5[0] + i * nquad0, 1,
486 &tmp5[0] + i * nquad0, 1);
487 }
488 Vmath::Vadd(nqtot, &tmp3[0], 1, &tmp5[0], 1, &tmp3[0], 1);
489}
490
491//---------------------------------------
492// Evaluation functions
493//---------------------------------------
494
496{
498 m_base[0]->GetBasisKey(), m_base[1]->GetBasisKey(),
499 m_base[2]->GetBasisKey());
500}
501
503{
505 m_base[0]->GetPointsKey());
507 m_base[1]->GetPointsKey());
509 m_base[2]->GetPointsKey());
510
512 bkey2);
513}
514
515/*
516 * @brief Get the coordinates #coords at the local coordinates
517 * #Lcoords
518 */
521{
522 int i;
523
524 ASSERTL1(Lcoords[0] <= -1.0 && Lcoords[0] >= 1.0 && Lcoords[1] <= -1.0 &&
525 Lcoords[1] >= 1.0 && Lcoords[2] <= -1.0 && Lcoords[2] >= 1.0,
526 "Local coordinates are not in region [-1,1]");
527
528 // m_geom->FillGeom(); // TODO: implement FillGeom()
529
530 for (i = 0; i < m_geom->GetCoordim(); ++i)
531 {
532 coords[i] = m_geom->GetCoord(i, Lcoords);
533 }
534}
535
537 Array<OneD, NekDouble> &coords_2,
538 Array<OneD, NekDouble> &coords_3)
539{
540 Expansion::v_GetCoords(coords_1, coords_2, coords_3);
541}
542
544 const NekDouble *data, const std::vector<unsigned int> &nummodes,
545 const int mode_offset, NekDouble *coeffs,
546 std::vector<LibUtilities::BasisType> &fromType)
547{
548 int data_order0 = nummodes[mode_offset];
549 int fillorder0 = min(m_base[0]->GetNumModes(), data_order0);
550 int data_order1 = nummodes[mode_offset + 1];
551 int order1 = m_base[1]->GetNumModes();
552 int fillorder1 = min(order1, data_order1);
553 int data_order2 = nummodes[mode_offset + 2];
554 int order2 = m_base[2]->GetNumModes();
555 int fillorder2 = min(order2, data_order2);
556
557 // Check if not same order or basis and if not make temp
558 // element to read in data
559 if (fromType[0] != m_base[0]->GetBasisType() ||
560 fromType[1] != m_base[1]->GetBasisType() ||
561 fromType[2] != m_base[2]->GetBasisType() || data_order0 != fillorder0 ||
562 data_order1 != fillorder1 || data_order2 != fillorder2)
563 {
564 // Construct a pyr with the appropriate basis type at our
565 // quadrature points, and one more to do a forwards
566 // transform. We can then copy the output to coeffs.
568 LibUtilities::BasisKey(fromType[0], data_order0,
569 m_base[0]->GetPointsKey()),
570 LibUtilities::BasisKey(fromType[1], data_order1,
571 m_base[1]->GetPointsKey()),
572 LibUtilities::BasisKey(fromType[2], data_order2,
573 m_base[2]->GetPointsKey()));
574
575 StdRegions::StdPyrExp tmpPyr2(m_base[0]->GetBasisKey(),
576 m_base[1]->GetBasisKey(),
577 m_base[2]->GetBasisKey());
578
579 Array<OneD, const NekDouble> tmpData(tmpPyr.GetNcoeffs(), data);
580 Array<OneD, NekDouble> tmpBwd(tmpPyr2.GetTotPoints());
581 Array<OneD, NekDouble> tmpOut(tmpPyr2.GetNcoeffs());
582
583 tmpPyr.BwdTrans(tmpData, tmpBwd);
584 tmpPyr2.FwdTrans(tmpBwd, tmpOut);
585 Vmath::Vcopy(tmpOut.size(), &tmpOut[0], 1, coeffs, 1);
586 }
587 else
588 {
589 Vmath::Vcopy(m_ncoeffs, &data[0], 1, coeffs, 1);
590 }
591}
592
593/**
594 * Given the local cartesian coordinate \a Lcoord evaluate the
595 * value of physvals at this point by calling through to the
596 * StdExpansion method
597 */
599 const Array<OneD, const NekDouble> &Lcoord,
600 const Array<OneD, const NekDouble> &physvals)
601{
602 // Evaluate point in local coordinates.
603 return StdExpansion3D::v_PhysEvaluate(Lcoord, physvals);
604}
605
607 const Array<OneD, const NekDouble> &physvals)
608{
609 Array<OneD, NekDouble> Lcoord(3);
610
611 ASSERTL0(m_geom, "m_geom not defined");
612
613 // TODO: check GetLocCoords()
614 m_geom->GetLocCoords(coord, Lcoord);
615
616 return StdExpansion3D::v_PhysEvaluate(Lcoord, physvals);
617}
618
620 const Array<OneD, const NekDouble> &inarray,
621 std::array<NekDouble, 3> &firstOrderDerivs)
622{
623 Array<OneD, NekDouble> Lcoord(3);
624 ASSERTL0(m_geom, "m_geom not defined");
625 m_geom->GetLocCoords(coord, Lcoord);
626 return StdPyrExp::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
627}
628
629//---------------------------------------
630// Helper functions
631//---------------------------------------
632
633void PyrExp::v_GetTracePhysMap(const int face, Array<OneD, int> &outarray)
634{
635 int nquad0 = m_base[0]->GetNumPoints();
636 int nquad1 = m_base[1]->GetNumPoints();
637 int nquad2 = m_base[2]->GetNumPoints();
638
639 int nq0 = 0;
640 int nq1 = 0;
641
642 switch (face)
643 {
644 case 0:
645 nq0 = nquad0;
646 nq1 = nquad1;
647 if (outarray.size() != nq0 * nq1)
648 {
649 outarray = Array<OneD, int>(nq0 * nq1);
650 }
651
652 // Directions A and B positive
653 for (int i = 0; i < nquad0 * nquad1; ++i)
654 {
655 outarray[i] = i;
656 }
657
658 break;
659 case 1:
660 nq0 = nquad0;
661 nq1 = nquad2;
662 if (outarray.size() != nq0 * nq1)
663 {
664 outarray = Array<OneD, int>(nq0 * nq1);
665 }
666
667 // Direction A and B positive
668 for (int k = 0; k < nquad2; k++)
669 {
670 for (int i = 0; i < nquad0; ++i)
671 {
672 outarray[k * nquad0 + i] = (nquad0 * nquad1 * k) + i;
673 }
674 }
675
676 break;
677 case 2:
678 nq0 = nquad1;
679 nq1 = nquad2;
680 if (outarray.size() != nq0 * nq1)
681 {
682 outarray = Array<OneD, int>(nq0 * nq1);
683 }
684
685 // Directions A and B positive
686 for (int j = 0; j < nquad1 * nquad2; ++j)
687 {
688 outarray[j] = nquad0 - 1 + j * nquad0;
689 }
690 break;
691 case 3:
692
693 nq0 = nquad0;
694 nq1 = nquad2;
695 if (outarray.size() != nq0 * nq1)
696 {
697 outarray = Array<OneD, int>(nq0 * nq1);
698 }
699
700 // Direction A and B positive
701 for (int k = 0; k < nquad2; k++)
702 {
703 for (int i = 0; i < nquad0; ++i)
704 {
705 outarray[k * nquad0 + i] =
706 nquad0 * (nquad1 - 1) + (nquad0 * nquad1 * k) + i;
707 }
708 }
709 break;
710 case 4:
711 nq0 = nquad1;
712 nq1 = nquad2;
713
714 if (outarray.size() != nq0 * nq1)
715 {
716 outarray = Array<OneD, int>(nq0 * nq1);
717 }
718
719 // Directions A and B positive
720 for (int j = 0; j < nquad1 * nquad2; ++j)
721 {
722 outarray[j] = j * nquad0;
723 }
724 break;
725 default:
726 ASSERTL0(false, "face value (> 4) is out of range");
727 break;
728 }
729}
730
732{
733 const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
734 GetGeom()->GetMetricInfo();
735
737 for (int i = 0; i < ptsKeys.size(); ++i)
738 {
739 // Need at least 2 points for computing normals
740 if (ptsKeys[i].GetNumPoints() == 1)
741 {
742 LibUtilities::PointsKey pKey(2, ptsKeys[i].GetPointsType());
743 ptsKeys[i] = pKey;
744 }
745 }
746
747 SpatialDomains::GeomType type = geomFactors->GetGtype();
749 geomFactors->GetDerivFactors(ptsKeys);
750 const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
751
752 LibUtilities::BasisKey tobasis0 = GetTraceBasisKey(face, 0);
753 LibUtilities::BasisKey tobasis1 = GetTraceBasisKey(face, 1);
754
755 // Number of quadrature points in face expansion.
756 int nq_face = tobasis0.GetNumPoints() * tobasis1.GetNumPoints();
757
758 int vCoordDim = GetCoordim();
759 int i;
760
763 for (i = 0; i < vCoordDim; ++i)
764 {
765 normal[i] = Array<OneD, NekDouble>(nq_face);
766 }
767
768 size_t nqb = nq_face;
769 size_t nbnd = face;
772
773 // Regular geometry case
774 if (type == SpatialDomains::eRegular ||
776 {
777 NekDouble fac;
778 // Set up normals
779 switch (face)
780 {
781 case 0:
782 {
783 for (i = 0; i < vCoordDim; ++i)
784 {
785 normal[i][0] = -df[3 * i + 2][0];
786 }
787 break;
788 }
789 case 1:
790 {
791 for (i = 0; i < vCoordDim; ++i)
792 {
793 normal[i][0] = -df[3 * i + 1][0];
794 }
795 break;
796 }
797 case 2:
798 {
799 for (i = 0; i < vCoordDim; ++i)
800 {
801 normal[i][0] = df[3 * i][0] + df[3 * i + 2][0];
802 }
803 break;
804 }
805 case 3:
806 {
807 for (i = 0; i < vCoordDim; ++i)
808 {
809 normal[i][0] = df[3 * i + 1][0] + df[3 * i + 2][0];
810 }
811 break;
812 }
813 case 4:
814 {
815 for (i = 0; i < vCoordDim; ++i)
816 {
817 normal[i][0] = -df[3 * i][0];
818 }
819 break;
820 }
821 default:
822 ASSERTL0(false, "face is out of range (face < 4)");
823 }
824
825 // Normalise resulting vector.
826 fac = 0.0;
827 for (i = 0; i < vCoordDim; ++i)
828 {
829 fac += normal[i][0] * normal[i][0];
830 }
831 fac = 1.0 / sqrt(fac);
832
833 Vmath::Fill(nqb, fac, length, 1);
834
835 for (i = 0; i < vCoordDim; ++i)
836 {
837 Vmath::Fill(nq_face, fac * normal[i][0], normal[i], 1);
838 }
839 }
840 else
841 {
842 // Set up deformed normals.
843 int j, k;
844
845 int nq0 = ptsKeys[0].GetNumPoints();
846 int nq1 = ptsKeys[1].GetNumPoints();
847 int nq2 = ptsKeys[2].GetNumPoints();
848 int nq01 = nq0 * nq1;
849 int nqtot;
850
851 // Determine number of quadrature points on the face.
852 if (face == 0)
853 {
854 nqtot = nq0 * nq1;
855 }
856 else if (face == 1 || face == 3)
857 {
858 nqtot = nq0 * nq2;
859 }
860 else
861 {
862 nqtot = nq1 * nq2;
863 }
864
867
868 Array<OneD, NekDouble> faceJac(nqtot);
869 Array<OneD, NekDouble> normals(vCoordDim * nqtot, 0.0);
870
871 // Extract Jacobian along face and recover local derivatives
872 // (dx/dr) for polynomial interpolation by multiplying m_gmat by
873 // jacobian
874 switch (face)
875 {
876 case 0:
877 {
878 for (j = 0; j < nq01; ++j)
879 {
880 normals[j] = -df[2][j] * jac[j];
881 normals[nqtot + j] = -df[5][j] * jac[j];
882 normals[2 * nqtot + j] = -df[8][j] * jac[j];
883 faceJac[j] = jac[j];
884 }
885
886 points0 = ptsKeys[0];
887 points1 = ptsKeys[1];
888 break;
889 }
890
891 case 1:
892 {
893 for (j = 0; j < nq0; ++j)
894 {
895 for (k = 0; k < nq2; ++k)
896 {
897 int tmp = j + nq01 * k;
898 normals[j + k * nq0] = -df[1][tmp] * jac[tmp];
899 normals[nqtot + j + k * nq0] = -df[4][tmp] * jac[tmp];
900 normals[2 * nqtot + j + k * nq0] =
901 -df[7][tmp] * jac[tmp];
902 faceJac[j + k * nq0] = jac[tmp];
903 }
904 }
905
906 points0 = ptsKeys[0];
907 points1 = ptsKeys[2];
908 break;
909 }
910
911 case 2:
912 {
913 for (j = 0; j < nq1; ++j)
914 {
915 for (k = 0; k < nq2; ++k)
916 {
917 int tmp = nq0 - 1 + nq0 * j + nq01 * k;
918 normals[j + k * nq1] =
919 (df[0][tmp] + df[2][tmp]) * jac[tmp];
920 normals[nqtot + j + k * nq1] =
921 (df[3][tmp] + df[5][tmp]) * jac[tmp];
922 normals[2 * nqtot + j + k * nq1] =
923 (df[6][tmp] + df[8][tmp]) * jac[tmp];
924 faceJac[j + k * nq1] = jac[tmp];
925 }
926 }
927
928 points0 = ptsKeys[1];
929 points1 = ptsKeys[2];
930 break;
931 }
932
933 case 3:
934 {
935 for (j = 0; j < nq0; ++j)
936 {
937 for (k = 0; k < nq2; ++k)
938 {
939 int tmp = nq0 * (nq1 - 1) + j + nq01 * k;
940 normals[j + k * nq0] =
941 (df[1][tmp] + df[2][tmp]) * jac[tmp];
942 normals[nqtot + j + k * nq0] =
943 (df[4][tmp] + df[5][tmp]) * jac[tmp];
944 normals[2 * nqtot + j + k * nq0] =
945 (df[7][tmp] + df[8][tmp]) * jac[tmp];
946 faceJac[j + k * nq0] = jac[tmp];
947 }
948 }
949
950 points0 = ptsKeys[0];
951 points1 = ptsKeys[2];
952 break;
953 }
954
955 case 4:
956 {
957 for (j = 0; j < nq1; ++j)
958 {
959 for (k = 0; k < nq2; ++k)
960 {
961 int tmp = j * nq0 + nq01 * k;
962 normals[j + k * nq1] = -df[0][tmp] * jac[tmp];
963 normals[nqtot + j + k * nq1] = -df[3][tmp] * jac[tmp];
964 normals[2 * nqtot + j + k * nq1] =
965 -df[6][tmp] * jac[tmp];
966 faceJac[j + k * nq1] = jac[tmp];
967 }
968 }
969
970 points0 = ptsKeys[1];
971 points1 = ptsKeys[2];
972 break;
973 }
974
975 default:
976 ASSERTL0(false, "face is out of range (face < 4)");
977 }
978
979 Array<OneD, NekDouble> work(nq_face, 0.0);
980 // Interpolate Jacobian and invert
981 LibUtilities::Interp2D(points0, points1, faceJac,
982 tobasis0.GetPointsKey(), tobasis1.GetPointsKey(),
983 work);
984 Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
985
986 // Interpolate normal and multiply by inverse Jacobian.
987 for (i = 0; i < vCoordDim; ++i)
988 {
989 LibUtilities::Interp2D(points0, points1, &normals[i * nqtot],
990 tobasis0.GetPointsKey(),
991 tobasis1.GetPointsKey(), &normal[i][0]);
992 Vmath::Vmul(nq_face, work, 1, normal[i], 1, normal[i], 1);
993 }
994
995 // Normalise to obtain unit normals.
996 Vmath::Zero(nq_face, work, 1);
997 for (i = 0; i < GetCoordim(); ++i)
998 {
999 Vmath::Vvtvp(nq_face, normal[i], 1, normal[i], 1, work, 1, work, 1);
1000 }
1001
1002 Vmath::Vsqrt(nq_face, work, 1, work, 1);
1003 Vmath::Sdiv(nq_face, 1.0, work, 1, work, 1);
1004
1005 Vmath::Vcopy(nqb, work, 1, length, 1);
1006
1007 for (i = 0; i < GetCoordim(); ++i)
1008 {
1009 Vmath::Vmul(nq_face, normal[i], 1, work, 1, normal[i], 1);
1010 }
1011 }
1012}
1013
1015 const StdRegions::StdMatrixKey &mkey)
1016{
1017 int nq = GetTotPoints();
1018
1019 // Calculate sqrt of the Jacobian
1021 Array<OneD, NekDouble> sqrt_jac(nq);
1022 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1023 {
1024 Vmath::Vsqrt(nq, jac, 1, sqrt_jac, 1);
1025 }
1026 else
1027 {
1028 Vmath::Fill(nq, sqrt(jac[0]), sqrt_jac, 1);
1029 }
1030
1031 // Multiply array by sqrt(Jac)
1032 Vmath::Vmul(nq, sqrt_jac, 1, array, 1, array, 1);
1033
1034 // Apply std region filter
1035 StdPyrExp::v_SVVLaplacianFilter(array, mkey);
1036
1037 // Divide by sqrt(Jac)
1038 Vmath::Vdiv(nq, array, 1, sqrt_jac, 1, array, 1);
1039}
1040
1041//---------------------------------------
1042// Matrix creation functions
1043//---------------------------------------
1044
1046{
1047 DNekMatSharedPtr returnval;
1048
1049 switch (mkey.GetMatrixType())
1050 {
1057 returnval = Expansion3D::v_GenMatrix(mkey);
1058 break;
1059 default:
1060 returnval = StdPyrExp::v_GenMatrix(mkey);
1061 }
1062
1063 return returnval;
1064}
1065
1067{
1068 LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
1069 LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
1070 LibUtilities::BasisKey bkey2 = m_base[2]->GetBasisKey();
1073
1074 return tmp->GetStdMatrix(mkey);
1075}
1076
1078{
1079 return m_matrixManager[mkey];
1080}
1081
1083{
1084 m_matrixManager.DeleteObject(mkey);
1085}
1086
1088{
1089 return m_staticCondMatrixManager[mkey];
1090}
1091
1093{
1094 m_staticCondMatrixManager.DeleteObject(mkey);
1095}
1096
1098{
1099 if (m_metrics.count(eMetricQuadrature) == 0)
1100 {
1102 }
1103
1104 int i, j;
1105 const unsigned int nqtot = GetTotPoints();
1106 const unsigned int dim = 3;
1107 const MetricType m[3][3] = {
1111
1112 for (unsigned int i = 0; i < dim; ++i)
1113 {
1114 for (unsigned int j = i; j < dim; ++j)
1115 {
1116 m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1117 }
1118 }
1119
1120 // Define shorthand synonyms for m_metrics storage
1121 Array<OneD, NekDouble> g0(m_metrics[m[0][0]]);
1122 Array<OneD, NekDouble> g1(m_metrics[m[1][1]]);
1123 Array<OneD, NekDouble> g2(m_metrics[m[2][2]]);
1124 Array<OneD, NekDouble> g3(m_metrics[m[0][1]]);
1125 Array<OneD, NekDouble> g4(m_metrics[m[0][2]]);
1126 Array<OneD, NekDouble> g5(m_metrics[m[1][2]]);
1127
1128 // Allocate temporary storage
1129 Array<OneD, NekDouble> alloc(9 * nqtot, 0.0);
1130 Array<OneD, NekDouble> h0(nqtot, alloc);
1131 Array<OneD, NekDouble> h1(nqtot, alloc + 1 * nqtot);
1132 Array<OneD, NekDouble> h2(nqtot, alloc + 2 * nqtot);
1133 Array<OneD, NekDouble> wsp1(nqtot, alloc + 3 * nqtot);
1134 Array<OneD, NekDouble> wsp2(nqtot, alloc + 4 * nqtot);
1135 Array<OneD, NekDouble> wsp3(nqtot, alloc + 5 * nqtot);
1136 Array<OneD, NekDouble> wsp4(nqtot, alloc + 6 * nqtot);
1137 Array<OneD, NekDouble> wsp5(nqtot, alloc + 7 * nqtot);
1138 Array<OneD, NekDouble> wsp6(nqtot, alloc + 8 * nqtot);
1139
1141 m_metricinfo->GetDerivFactors(GetPointsKeys());
1142 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
1143 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
1144 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
1145 const unsigned int nquad0 = m_base[0]->GetNumPoints();
1146 const unsigned int nquad1 = m_base[1]->GetNumPoints();
1147 const unsigned int nquad2 = m_base[2]->GetNumPoints();
1148
1149 // Populate collapsed coordinate arrays h0, h1 and h2.
1150 for (j = 0; j < nquad2; ++j)
1151 {
1152 for (i = 0; i < nquad1; ++i)
1153 {
1154 Vmath::Fill(nquad0, 2.0 / (1.0 - z2[j]),
1155 &h0[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1156 Vmath::Fill(nquad0, 1.0 / (1.0 - z2[j]),
1157 &h1[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1158 Vmath::Fill(nquad0, (1.0 + z1[i]) / (1.0 - z2[j]),
1159 &h2[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1160 }
1161 }
1162 for (i = 0; i < nquad0; i++)
1163 {
1164 Blas::Dscal(nquad1 * nquad2, 1 + z0[i], &h1[0] + i, nquad0);
1165 }
1166
1167 // Step 3. Construct combined metric terms for physical space to
1168 // collapsed coordinate system.
1169 // Order of construction optimised to minimise temporary storage
1170 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1171 {
1172 // f_{1k}
1173 Vmath::Vvtvvtp(nqtot, &df[0][0], 1, &h0[0], 1, &df[2][0], 1, &h1[0], 1,
1174 &wsp1[0], 1);
1175 Vmath::Vvtvvtp(nqtot, &df[3][0], 1, &h0[0], 1, &df[5][0], 1, &h1[0], 1,
1176 &wsp2[0], 1);
1177 Vmath::Vvtvvtp(nqtot, &df[6][0], 1, &h0[0], 1, &df[8][0], 1, &h1[0], 1,
1178 &wsp3[0], 1);
1179
1180 // g0
1181 Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp1[0], 1, &wsp2[0], 1, &wsp2[0],
1182 1, &g0[0], 1);
1183 Vmath::Vvtvp(nqtot, &wsp3[0], 1, &wsp3[0], 1, &g0[0], 1, &g0[0], 1);
1184
1185 // g4
1186 Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp1[0], 1, &df[5][0], 1, &wsp2[0],
1187 1, &g4[0], 1);
1188 Vmath::Vvtvp(nqtot, &df[8][0], 1, &wsp3[0], 1, &g4[0], 1, &g4[0], 1);
1189
1190 // f_{2k}
1191 Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &h0[0], 1, &df[2][0], 1, &h2[0], 1,
1192 &wsp4[0], 1);
1193 Vmath::Vvtvvtp(nqtot, &df[4][0], 1, &h0[0], 1, &df[5][0], 1, &h2[0], 1,
1194 &wsp5[0], 1);
1195 Vmath::Vvtvvtp(nqtot, &df[7][0], 1, &h0[0], 1, &df[8][0], 1, &h2[0], 1,
1196 &wsp6[0], 1);
1197
1198 // g1
1199 Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0],
1200 1, &g1[0], 1);
1201 Vmath::Vvtvp(nqtot, &wsp6[0], 1, &wsp6[0], 1, &g1[0], 1, &g1[0], 1);
1202
1203 // g3
1204 Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp4[0], 1, &wsp2[0], 1, &wsp5[0],
1205 1, &g3[0], 1);
1206 Vmath::Vvtvp(nqtot, &wsp3[0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1207
1208 // g5
1209 Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp4[0], 1, &df[5][0], 1, &wsp5[0],
1210 1, &g5[0], 1);
1211 Vmath::Vvtvp(nqtot, &df[8][0], 1, &wsp6[0], 1, &g5[0], 1, &g5[0], 1);
1212
1213 // g2
1214 Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &df[2][0], 1, &df[5][0], 1,
1215 &df[5][0], 1, &g2[0], 1);
1216 Vmath::Vvtvp(nqtot, &df[8][0], 1, &df[8][0], 1, &g2[0], 1, &g2[0], 1);
1217 }
1218 else
1219 {
1220 // f_{1k}
1221 Vmath::Svtsvtp(nqtot, df[0][0], &h0[0], 1, df[2][0], &h1[0], 1,
1222 &wsp1[0], 1);
1223 Vmath::Svtsvtp(nqtot, df[3][0], &h0[0], 1, df[5][0], &h1[0], 1,
1224 &wsp2[0], 1);
1225 Vmath::Svtsvtp(nqtot, df[6][0], &h0[0], 1, df[8][0], &h1[0], 1,
1226 &wsp3[0], 1);
1227
1228 // g0
1229 Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp1[0], 1, &wsp2[0], 1, &wsp2[0],
1230 1, &g0[0], 1);
1231 Vmath::Vvtvp(nqtot, &wsp3[0], 1, &wsp3[0], 1, &g0[0], 1, &g0[0], 1);
1232
1233 // g4
1234 Vmath::Svtsvtp(nqtot, df[2][0], &wsp1[0], 1, df[5][0], &wsp2[0], 1,
1235 &g4[0], 1);
1236 Vmath::Svtvp(nqtot, df[8][0], &wsp3[0], 1, &g4[0], 1, &g4[0], 1);
1237
1238 // f_{2k}
1239 Vmath::Svtsvtp(nqtot, df[1][0], &h0[0], 1, df[2][0], &h2[0], 1,
1240 &wsp4[0], 1);
1241 Vmath::Svtsvtp(nqtot, df[4][0], &h0[0], 1, df[5][0], &h2[0], 1,
1242 &wsp5[0], 1);
1243 Vmath::Svtsvtp(nqtot, df[7][0], &h0[0], 1, df[8][0], &h2[0], 1,
1244 &wsp6[0], 1);
1245
1246 // g1
1247 Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0],
1248 1, &g1[0], 1);
1249 Vmath::Vvtvp(nqtot, &wsp6[0], 1, &wsp6[0], 1, &g1[0], 1, &g1[0], 1);
1250
1251 // g3
1252 Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp4[0], 1, &wsp2[0], 1, &wsp5[0],
1253 1, &g3[0], 1);
1254 Vmath::Vvtvp(nqtot, &wsp3[0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1255
1256 // g5
1257 Vmath::Svtsvtp(nqtot, df[2][0], &wsp4[0], 1, df[5][0], &wsp5[0], 1,
1258 &g5[0], 1);
1259 Vmath::Svtvp(nqtot, df[8][0], &wsp6[0], 1, &g5[0], 1, &g5[0], 1);
1260
1261 // g2
1262 Vmath::Fill(nqtot,
1263 df[2][0] * df[2][0] + df[5][0] * df[5][0] +
1264 df[8][0] * df[8][0],
1265 &g2[0], 1);
1266 }
1267
1268 for (unsigned int i = 0; i < dim; ++i)
1269 {
1270 for (unsigned int j = i; j < dim; ++j)
1271 {
1273 }
1274 }
1275}
1276
1278 const Array<OneD, const NekDouble> &inarray,
1280{
1281 // This implementation is only valid when there are no coefficients
1282 // associated to the Laplacian operator
1283 if (m_metrics.count(eMetricLaplacian00) == 0)
1284 {
1286 }
1287
1288 int nquad0 = m_base[0]->GetNumPoints();
1289 int nquad1 = m_base[1]->GetNumPoints();
1290 int nq2 = m_base[2]->GetNumPoints();
1291 int nqtot = nquad0 * nquad1 * nq2;
1292
1293 ASSERTL1(wsp.size() >= 6 * nqtot, "Insufficient workspace size.");
1294 ASSERTL1(m_ncoeffs <= nqtot, "Workspace not set up for ncoeffs > nqtot");
1295
1296 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
1297 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
1298 const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
1299 const Array<OneD, const NekDouble> &dbase0 = m_base[0]->GetDbdata();
1300 const Array<OneD, const NekDouble> &dbase1 = m_base[1]->GetDbdata();
1301 const Array<OneD, const NekDouble> &dbase2 = m_base[2]->GetDbdata();
1302 const Array<OneD, const NekDouble> &metric00 =
1303 m_metrics[eMetricLaplacian00];
1304 const Array<OneD, const NekDouble> &metric01 =
1305 m_metrics[eMetricLaplacian01];
1306 const Array<OneD, const NekDouble> &metric02 =
1307 m_metrics[eMetricLaplacian02];
1308 const Array<OneD, const NekDouble> &metric11 =
1309 m_metrics[eMetricLaplacian11];
1310 const Array<OneD, const NekDouble> &metric12 =
1311 m_metrics[eMetricLaplacian12];
1312 const Array<OneD, const NekDouble> &metric22 =
1313 m_metrics[eMetricLaplacian22];
1314
1315 // Allocate temporary storage
1316 Array<OneD, NekDouble> wsp0(2 * nqtot, wsp);
1317 Array<OneD, NekDouble> wsp1(nqtot, wsp + 1 * nqtot);
1318 Array<OneD, NekDouble> wsp2(nqtot, wsp + 2 * nqtot);
1319 Array<OneD, NekDouble> wsp3(nqtot, wsp + 3 * nqtot);
1320 Array<OneD, NekDouble> wsp4(nqtot, wsp + 4 * nqtot);
1321 Array<OneD, NekDouble> wsp5(nqtot, wsp + 5 * nqtot);
1322
1323 // LAPLACIAN MATRIX OPERATION
1324 // wsp1 = du_dxi1 = D_xi1 * inarray = D_xi1 * u
1325 // wsp2 = du_dxi2 = D_xi2 * inarray = D_xi2 * u
1326 // wsp2 = du_dxi3 = D_xi3 * inarray = D_xi3 * u
1327 StdExpansion3D::PhysTensorDeriv(inarray, wsp0, wsp1, wsp2);
1328
1329 // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1330 // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1331 // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1332 // especially for this purpose
1333 Vmath::Vvtvvtp(nqtot, &metric00[0], 1, &wsp0[0], 1, &metric01[0], 1,
1334 &wsp1[0], 1, &wsp3[0], 1);
1335 Vmath::Vvtvp(nqtot, &metric02[0], 1, &wsp2[0], 1, &wsp3[0], 1, &wsp3[0], 1);
1336 Vmath::Vvtvvtp(nqtot, &metric01[0], 1, &wsp0[0], 1, &metric11[0], 1,
1337 &wsp1[0], 1, &wsp4[0], 1);
1338 Vmath::Vvtvp(nqtot, &metric12[0], 1, &wsp2[0], 1, &wsp4[0], 1, &wsp4[0], 1);
1339 Vmath::Vvtvvtp(nqtot, &metric02[0], 1, &wsp0[0], 1, &metric12[0], 1,
1340 &wsp1[0], 1, &wsp5[0], 1);
1341 Vmath::Vvtvp(nqtot, &metric22[0], 1, &wsp2[0], 1, &wsp5[0], 1, &wsp5[0], 1);
1342
1343 // outarray = m = (D_xi1 * B)^T * k
1344 // wsp1 = n = (D_xi2 * B)^T * l
1345 IProductWRTBase_SumFacKernel(dbase0, base1, base2, wsp3, outarray, wsp0,
1346 false, true, true);
1347 IProductWRTBase_SumFacKernel(base0, dbase1, base2, wsp4, wsp2, wsp0, true,
1348 false, true);
1349 Vmath::Vadd(m_ncoeffs, wsp2.get(), 1, outarray.get(), 1, outarray.get(), 1);
1350 IProductWRTBase_SumFacKernel(base0, base1, dbase2, wsp5, wsp2, wsp0, true,
1351 true, false);
1352 Vmath::Vadd(m_ncoeffs, wsp2.get(), 1, outarray.get(), 1, outarray.get(), 1);
1353}
1354
1355/** @brief: This method gets all of the factors which are
1356 required as part of the Gradient Jump Penalty
1357 stabilisation and involves the product of the normal and
1358 geometric factors along the element trace.
1359*/
1360void PyrExp::v_NormalTraceDerivFactors(
1361 Array<OneD, Array<OneD, NekDouble>> &d0factors,
1362 Array<OneD, Array<OneD, NekDouble>> &d1factors,
1363 Array<OneD, Array<OneD, NekDouble>> &d2factors)
1364{
1365 int nquad0 = GetNumPoints(0);
1366 int nquad1 = GetNumPoints(1);
1367 int nquad2 = GetNumPoints(2);
1368
1369 const Array<TwoD, const NekDouble> &df =
1370 m_metricinfo->GetDerivFactors(GetPointsKeys());
1371
1372 if (d0factors.size() != 5)
1373 {
1374 d0factors = Array<OneD, Array<OneD, NekDouble>>(5);
1375 d1factors = Array<OneD, Array<OneD, NekDouble>>(5);
1376 d2factors = Array<OneD, Array<OneD, NekDouble>>(5);
1377 }
1378
1379 if (d0factors[0].size() != nquad0 * nquad1)
1380 {
1381 d0factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1382 d1factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1383 d2factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1384 }
1385
1386 if (d0factors[1].size() != nquad0 * nquad2)
1387 {
1388 d0factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1389 d0factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1390 d1factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1391 d1factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1392 d2factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1393 d2factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1394 }
1395
1396 if (d0factors[2].size() != nquad1 * nquad2)
1397 {
1398 d0factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1399 d0factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1400 d1factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1401 d1factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1402 d2factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1403 d2factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1404 }
1405
1406 // Outwards normals
1407 const Array<OneD, const Array<OneD, NekDouble>> &normal_0 =
1408 GetTraceNormal(0);
1409 const Array<OneD, const Array<OneD, NekDouble>> &normal_1 =
1410 GetTraceNormal(1);
1411 const Array<OneD, const Array<OneD, NekDouble>> &normal_2 =
1412 GetTraceNormal(2);
1413 const Array<OneD, const Array<OneD, NekDouble>> &normal_3 =
1414 GetTraceNormal(3);
1415 const Array<OneD, const Array<OneD, NekDouble>> &normal_4 =
1416 GetTraceNormal(4);
1417
1418 int ncoords = normal_0.size();
1419
1420 // first gather together standard cartesian inner products
1421 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1422 {
1423 // face 0
1424 for (int i = 0; i < nquad0 * nquad1; ++i)
1425 {
1426 d0factors[0][i] = df[0][i] * normal_0[0][i];
1427 d1factors[0][i] = df[1][i] * normal_0[0][i];
1428 d2factors[0][i] = df[2][i] * normal_0[0][i];
1429 }
1430
1431 for (int n = 1; n < ncoords; ++n)
1432 {
1433 for (int i = 0; i < nquad0 * nquad1; ++i)
1434 {
1435 d0factors[0][i] += df[3 * n][i] * normal_0[n][i];
1436 d1factors[0][i] += df[3 * n + 1][i] * normal_0[n][i];
1437 d2factors[0][i] += df[3 * n + 2][i] * normal_0[n][i];
1438 }
1439 }
1440
1441 // faces 1 and 3
1442 for (int j = 0; j < nquad2; ++j)
1443 {
1444 for (int i = 0; i < nquad0; ++i)
1445 {
1446 d0factors[1][j * nquad0 + i] = df[0][j * nquad0 * nquad1 + i] *
1447 normal_1[0][j * nquad0 + i];
1448 d1factors[1][j * nquad0 + i] = df[1][j * nquad0 * nquad1 + i] *
1449 normal_1[0][j * nquad0 + i];
1450 d2factors[1][j * nquad0 + i] = df[2][j * nquad0 * nquad1 + i] *
1451 normal_1[0][j * nquad0 + i];
1452
1453 d0factors[3][j * nquad0 + i] =
1454 df[0][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1455 normal_3[0][j * nquad0 + i];
1456 d1factors[3][j * nquad0 + i] =
1457 df[1][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1458 normal_3[0][j * nquad0 + i];
1459 d2factors[3][j * nquad0 + i] =
1460 df[2][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1461 normal_3[0][j * nquad0 + i];
1462 }
1463 }
1464
1465 for (int n = 1; n < ncoords; ++n)
1466 {
1467 for (int j = 0; j < nquad2; ++j)
1468 {
1469 for (int i = 0; i < nquad0; ++i)
1470 {
1471 d0factors[1][j * nquad0 + i] +=
1472 df[3 * n][j * nquad0 * nquad1 + i] *
1473 normal_1[0][j * nquad0 + i];
1474 d1factors[1][j * nquad0 + i] +=
1475 df[3 * n + 1][j * nquad0 * nquad1 + i] *
1476 normal_1[0][j * nquad0 + i];
1477 d2factors[1][j * nquad0 + i] +=
1478 df[3 * n + 2][j * nquad0 * nquad1 + i] *
1479 normal_1[0][j * nquad0 + i];
1480
1481 d0factors[3][j * nquad0 + i] +=
1482 df[3 * n][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1483 normal_3[0][j * nquad0 + i];
1484 d1factors[3][j * nquad0 + i] +=
1485 df[3 * n + 1][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1486 normal_3[0][j * nquad0 + i];
1487 d2factors[3][j * nquad0 + i] +=
1488 df[3 * n + 2][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1489 normal_3[0][j * nquad0 + i];
1490 }
1491 }
1492 }
1493
1494 // faces 2 and 4
1495 for (int j = 0; j < nquad2; ++j)
1496 {
1497 for (int i = 0; i < nquad1; ++i)
1498 {
1499 d0factors[2][j * nquad1 + i] =
1500 df[0][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1501 normal_2[0][j * nquad1 + i];
1502 d1factors[2][j * nquad1 + i] =
1503 df[1][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1504 normal_2[0][j * nquad1 + i];
1505 d2factors[2][j * nquad1 + i] =
1506 df[2][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1507 normal_2[0][j * nquad1 + i];
1508
1509 d0factors[4][j * nquad1 + i] =
1510 df[0][j * nquad0 * nquad1 + i * nquad0] *
1511 normal_4[0][j * nquad1 + i];
1512 d1factors[4][j * nquad1 + i] =
1513 df[1][j * nquad0 * nquad1 + i * nquad0] *
1514 normal_4[0][j * nquad1 + i];
1515 d2factors[4][j * nquad1 + i] =
1516 df[2][j * nquad0 * nquad1 + i * nquad0] *
1517 normal_4[0][j * nquad1 + i];
1518 }
1519 }
1520
1521 for (int n = 1; n < ncoords; ++n)
1522 {
1523 for (int j = 0; j < nquad2; ++j)
1524 {
1525 for (int i = 0; i < nquad1; ++i)
1526 {
1527 d0factors[2][j * nquad1 + i] +=
1528 df[3 * n][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1529 normal_2[n][j * nquad1 + i];
1530 d1factors[2][j * nquad1 + i] +=
1531 df[3 * n + 1]
1532 [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1533 normal_2[n][j * nquad1 + i];
1534 d2factors[2][j * nquad1 + i] +=
1535 df[3 * n + 2]
1536 [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1537 normal_2[n][j * nquad1 + i];
1538
1539 d0factors[4][j * nquad1 + i] +=
1540 df[3 * n][i * nquad0 + j * nquad0 * nquad1] *
1541 normal_4[n][j * nquad1 + i];
1542 d1factors[4][j * nquad1 + i] +=
1543 df[3 * n + 1][i * nquad0 + j * nquad0 * nquad1] *
1544 normal_4[n][j * nquad1 + i];
1545 d2factors[4][j * nquad1 + i] +=
1546 df[3 * n + 2][i * nquad0 + j * nquad0 * nquad1] *
1547 normal_4[n][j * nquad1 + i];
1548 }
1549 }
1550 }
1551 }
1552 else
1553 {
1554 // Face 0
1555 for (int i = 0; i < nquad0 * nquad1; ++i)
1556 {
1557 d0factors[0][i] = df[0][0] * normal_0[0][i];
1558 d1factors[0][i] = df[1][0] * normal_0[0][i];
1559 d2factors[0][i] = df[2][0] * normal_0[0][i];
1560 }
1561
1562 for (int n = 1; n < ncoords; ++n)
1563 {
1564 for (int i = 0; i < nquad0 * nquad1; ++i)
1565 {
1566 d0factors[0][i] += df[3 * n][0] * normal_0[n][i];
1567 d1factors[0][i] += df[3 * n + 1][0] * normal_0[n][i];
1568 d2factors[0][i] += df[3 * n + 2][0] * normal_0[n][i];
1569 }
1570 }
1571
1572 // faces 1 and 3
1573 for (int i = 0; i < nquad0 * nquad2; ++i)
1574 {
1575 d0factors[1][i] = df[0][0] * normal_1[0][i];
1576 d0factors[3][i] = df[0][0] * normal_3[0][i];
1577
1578 d1factors[1][i] = df[1][0] * normal_1[0][i];
1579 d1factors[3][i] = df[1][0] * normal_3[0][i];
1580
1581 d2factors[1][i] = df[2][0] * normal_1[0][i];
1582 d2factors[3][i] = df[2][0] * normal_3[0][i];
1583 }
1584
1585 for (int n = 1; n < ncoords; ++n)
1586 {
1587 for (int i = 0; i < nquad0 * nquad2; ++i)
1588 {
1589 d0factors[1][i] += df[3 * n][0] * normal_1[n][i];
1590 d0factors[3][i] += df[3 * n][0] * normal_3[n][i];
1591
1592 d1factors[1][i] += df[3 * n + 1][0] * normal_1[n][i];
1593 d1factors[3][i] += df[3 * n + 1][0] * normal_3[n][i];
1594
1595 d2factors[1][i] += df[3 * n + 2][0] * normal_1[n][i];
1596 d2factors[3][i] += df[3 * n + 2][0] * normal_3[n][i];
1597 }
1598 }
1599
1600 // faces 2 and 4
1601 for (int i = 0; i < nquad1 * nquad2; ++i)
1602 {
1603 d0factors[2][i] = df[0][0] * normal_2[0][i];
1604 d0factors[4][i] = df[0][0] * normal_4[0][i];
1605
1606 d1factors[2][i] = df[1][0] * normal_2[0][i];
1607 d1factors[4][i] = df[1][0] * normal_4[0][i];
1608
1609 d2factors[2][i] = df[2][0] * normal_2[0][i];
1610 d2factors[4][i] = df[2][0] * normal_4[0][i];
1611 }
1612
1613 for (int n = 1; n < ncoords; ++n)
1614 {
1615 for (int i = 0; i < nquad1 * nquad2; ++i)
1616 {
1617 d0factors[2][i] += df[3 * n][0] * normal_2[n][i];
1618 d0factors[4][i] += df[3 * n][0] * normal_4[n][i];
1619
1620 d1factors[2][i] += df[3 * n + 1][0] * normal_2[n][i];
1621 d1factors[4][i] += df[3 * n + 1][0] * normal_4[n][i];
1622
1623 d2factors[2][i] += df[3 * n + 2][0] * normal_2[n][i];
1624 d2factors[4][i] += df[3 * n + 2][0] * normal_4[n][i];
1625 }
1626 }
1627 }
1628}
1629
1630} // namespace LocalRegions
1631} // 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
Describes the specification for a Basis.
Definition: Basis.h:47
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:122
PointsKey GetPointsKey() const
Return distribution of points.
Definition: Basis.h:139
Defines a specification for a set of points.
Definition: Points.h:55
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
std::map< int, NormalVector > m_traceNormals
Definition: Expansion.h:278
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
Definition: Expansion.h:288
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:171
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:275
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: Expansion.cpp:535
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:276
virtual void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType) override
Definition: PyrExp.cpp:543
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrate the physical point list inarray over pyramidic region and return the value.
Definition: PyrExp.cpp:96
virtual void v_ComputeLaplacianMetric() override
Definition: PyrExp.cpp:1097
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: PyrExp.h:176
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: PyrExp.cpp:536
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdRegions::StdMatrixKey &mkey) override
Definition: PyrExp.cpp:1066
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
Definition: PyrExp.cpp:223
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey) override
Definition: PyrExp.cpp:1014
void v_ComputeTraceNormal(const int face) override
Definition: PyrExp.cpp:731
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
Definition: PyrExp.cpp:519
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: PyrExp.h:178
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculates the inner product .
Definition: PyrExp.cpp:344
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const override
Definition: PyrExp.cpp:495
void v_DropLocMatrix(const MatrixKey &mkey) override
Definition: PyrExp.cpp:1082
virtual void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp) override
Definition: PyrExp.cpp:1277
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: PyrExp.cpp:351
PyrExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, const SpatialDomains::PyrGeomSharedPtr &geom)
Constructor using BasisKey class for quadrature points and order definition.
Definition: PyrExp.cpp:45
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
Definition: PyrExp.cpp:1045
virtual StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const override
Definition: PyrExp.cpp:502
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const MatrixKey &mkey) override
Definition: PyrExp.cpp:1087
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
Definition: PyrExp.cpp:396
virtual void v_GetTracePhysMap(const int face, Array< OneD, int > &outarray) override
Definition: PyrExp.cpp:633
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
Calculate the derivative of the physical points.
Definition: PyrExp.cpp:124
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals) override
This function evaluates the expansion at a single (arbitrary) point of the domain.
Definition: PyrExp.cpp:606
void v_DropLocStaticCondMatrix(const MatrixKey &mkey) override
Definition: PyrExp.cpp:1092
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculate the inner product of inarray with respect to the basis B=base0*base1*base2 and put into out...
Definition: PyrExp.cpp:278
virtual DNekScalMatSharedPtr v_GetLocMatrix(const MatrixKey &mkey) override
Definition: PyrExp.cpp:1077
NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
Definition: PyrExp.cpp:598
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: PyrExp.cpp:284
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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)
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
const LibUtilities::PointsKeyVector GetPointsKeys() const
const LibUtilities::BasisKey GetTraceBasisKey(const int i, int k=-1) const
This function returns the basis key belonging to the i-th trace.
Definition: StdExpansion.h:305
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
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
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:87
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:151
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:137
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
Definition: Interp.cpp:103
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:236
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< PyrGeom > PyrGeomSharedPtr
Definition: PyrGeom.h:77
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::shared_ptr< StdPyrExp > StdPyrExpSharedPtr
Definition: StdPyrExp.h:258
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:529
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
svtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:746
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:207
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:617
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:569
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:354
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:245
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
Definition: Vmath.cpp:319
void Vdiv(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:280
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:43
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.cpp:687
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294