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