Nektar++
Loading...
Searching...
No Matches
PrismExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PrismExp.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: PrismExp routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
39
40using namespace std;
41
43{
44
46 const LibUtilities::BasisKey &Bb,
47 const LibUtilities::BasisKey &Bc,
49 : StdExpansion(LibUtilities::StdPrismData::getNumberOfCoefficients(
50 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
51 3, Ba, Bb, Bc),
52 StdExpansion3D(LibUtilities::StdPrismData::getNumberOfCoefficients(
53 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
54 Ba, Bb, Bc),
55 StdPrismExp(Ba, Bb, Bc), Expansion(geom), Expansion3D(geom),
56 m_matrixManager(
57 std::bind(&Expansion3D::CreateMatrix, this, std::placeholders::_1),
58 std::string("PrismExpMatrix")),
59 m_staticCondMatrixManager(std::bind(&Expansion::CreateStaticCondMatrix,
60 this, std::placeholders::_1),
61 std::string("PrismExpStaticCondMatrix"))
62{
63}
64
66 : StdExpansion(T), StdExpansion3D(T), StdPrismExp(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 prismatic
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 * \xi_2, \xi_3) J[i,j,k] d \bar \eta_1 d \xi_2 d \xi_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}, \xi_{2j}^{0,0},\xi_{3k}^{1,0})w_{i}^{0,0}
91 * w_{j}^{0,0} \hat w_{k}^{1,0} \f$ \n where \f$ inarray[i,j, k] =
92 * u(\bar \eta_{1i}^{0,0}, \xi_{2j}^{0,0},\xi_{3k}^{1,0}) \f$, \n
93 * \f$\hat w_{i}^{1,0} = \frac {w_{j}^{1,0}} {2} \f$ \n and \f$
94 * J[i,j,k] \f$ is the Jacobian 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 StdPrismExp version.
117 return StdPrismExp::v_Integral(tmp);
118}
119
120//----------------------------
121// Differentiation Methods
122//----------------------------
127{
128 int nqtot = GetTotPoints();
129
131 m_metricinfo->GetDerivFactors(GetPointsKeys());
132 Array<OneD, NekDouble> diff0(nqtot);
133 Array<OneD, NekDouble> diff1(nqtot);
134 Array<OneD, NekDouble> diff2(nqtot);
135
136 StdPrismExp::v_PhysDeriv(inarray, diff0, diff1, diff2);
137
138 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
139 {
140 if (out_d0.size())
141 {
142 Vmath::Vmul(nqtot, &df[0][0], 1, &diff0[0], 1, &out_d0[0], 1);
143 Vmath::Vvtvp(nqtot, &df[1][0], 1, &diff1[0], 1, &out_d0[0], 1,
144 &out_d0[0], 1);
145 Vmath::Vvtvp(nqtot, &df[2][0], 1, &diff2[0], 1, &out_d0[0], 1,
146 &out_d0[0], 1);
147 }
148
149 if (out_d1.size())
150 {
151 Vmath::Vmul(nqtot, &df[3][0], 1, &diff0[0], 1, &out_d1[0], 1);
152 Vmath::Vvtvp(nqtot, &df[4][0], 1, &diff1[0], 1, &out_d1[0], 1,
153 &out_d1[0], 1);
154 Vmath::Vvtvp(nqtot, &df[5][0], 1, &diff2[0], 1, &out_d1[0], 1,
155 &out_d1[0], 1);
156 }
157
158 if (out_d2.size())
159 {
160 Vmath::Vmul(nqtot, &df[6][0], 1, &diff0[0], 1, &out_d2[0], 1);
161 Vmath::Vvtvp(nqtot, &df[7][0], 1, &diff1[0], 1, &out_d2[0], 1,
162 &out_d2[0], 1);
163 Vmath::Vvtvp(nqtot, &df[8][0], 1, &diff2[0], 1, &out_d2[0], 1,
164 &out_d2[0], 1);
165 }
166 }
167 else // regular geometry
168 {
169 if (out_d0.size())
170 {
171 Vmath::Smul(nqtot, df[0][0], &diff0[0], 1, &out_d0[0], 1);
172 Blas::Daxpy(nqtot, df[1][0], &diff1[0], 1, &out_d0[0], 1);
173 Blas::Daxpy(nqtot, df[2][0], &diff2[0], 1, &out_d0[0], 1);
174 }
175
176 if (out_d1.size())
177 {
178 Vmath::Smul(nqtot, df[3][0], &diff0[0], 1, &out_d1[0], 1);
179 Blas::Daxpy(nqtot, df[4][0], &diff1[0], 1, &out_d1[0], 1);
180 Blas::Daxpy(nqtot, df[5][0], &diff2[0], 1, &out_d1[0], 1);
181 }
182
183 if (out_d2.size())
184 {
185 Vmath::Smul(nqtot, df[6][0], &diff0[0], 1, &out_d2[0], 1);
186 Blas::Daxpy(nqtot, df[7][0], &diff1[0], 1, &out_d2[0], 1);
187 Blas::Daxpy(nqtot, df[8][0], &diff2[0], 1, &out_d2[0], 1);
188 }
189 }
190}
191
192//---------------------------------------
193// Transforms
194//---------------------------------------
195
196/**
197 * \brief Forward transform from physical quadrature space stored in
198 * \a inarray and evaluate the expansion coefficients and store in \a
199 * (this)->m_coeffs
200 *
201 * Inputs:\n
202 *
203 * - \a inarray: array of physical quadrature points to be transformed
204 *
205 * Outputs:\n
206 *
207 * - (this)->_coeffs: updated array of expansion coefficients.
208 */
210 Array<OneD, NekDouble> &outarray)
211{
212 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
213 m_base[2]->Collocation())
214 {
215 Vmath::Vcopy(GetNcoeffs(), &inarray[0], 1, &outarray[0], 1);
216 }
217 else
218 {
219 v_IProductWRTBase(inarray, outarray);
220
221 // get Mass matrix inverse
223 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
224
225 // copy inarray in case inarray == outarray
226 DNekVec in(m_ncoeffs, outarray);
227 DNekVec out(m_ncoeffs, outarray, eWrapper);
228
229 out = (*matsys) * in;
230 }
231}
232
233//---------------------------------------
234// Inner product functions
235//---------------------------------------
236
237/**
238 * \brief Calculate the inner product of inarray with respect to the
239 * basis B=base0*base1*base2 and put into outarray:
240 *
241 * \f$ \begin{array}{rcl} I_{pqr} = (\phi_{pqr}, u)_{\delta} & = &
242 * \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \sum_{k=0}^{nq_2} \psi_{p}^{a}
243 * (\bar \eta_{1i}) \psi_{q}^{a} (\xi_{2j}) \psi_{pr}^{b} (\xi_{3k})
244 * w_i w_j w_k u(\bar \eta_{1,i} \xi_{2,j} \xi_{3,k}) J_{i,j,k}\\ & =
245 * & \sum_{i=0}^{nq_0} \psi_p^a(\bar \eta_{1,i}) \sum_{j=0}^{nq_1}
246 * \psi_{q}^a(\xi_{2,j}) \sum_{k=0}^{nq_2} \psi_{pr}^b u(\bar
247 * \eta_{1i},\xi_{2j},\xi_{3k}) J_{i,j,k} \end{array} \f$ \n
248 *
249 * where
250 *
251 * \f$ \phi_{pqr} (\xi_1 , \xi_2 , \xi_3) = \psi_p^a (\bar \eta_1)
252 * \psi_{q}^a (\xi_2) \psi_{pr}^b (\xi_3) \f$ \n
253 *
254 * which can be implemented as \n \f$f_{pr} (\xi_{3k}) =
255 * \sum_{k=0}^{nq_3} \psi_{pr}^b u(\bar \eta_{1i},\xi_{2j},\xi_{3k})
256 * J_{i,j,k} = {\bf B_3 U} \f$ \n \f$ g_{q} (\xi_{3k}) =
257 * \sum_{j=0}^{nq_1} \psi_{q}^a (\xi_{2j}) f_{pr} (\xi_{3k}) = {\bf
258 * B_2 F} \f$ \n \f$ (\phi_{pqr}, u)_{\delta} = \sum_{k=0}^{nq_0}
259 * \psi_{p}^a (\xi_{3k}) g_{q} (\xi_{3k}) = {\bf B_1 G} \f$
260 */
262 Array<OneD, NekDouble> &outarray)
263{
264 v_IProductWRTBase_SumFac(inarray, outarray);
265}
266
268 const Array<OneD, const NekDouble> &inarray,
269 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
270{
271 const int nquad0 = m_base[0]->GetNumPoints();
272 const int nquad1 = m_base[1]->GetNumPoints();
273 const int nquad2 = m_base[2]->GetNumPoints();
274 const int order0 = m_base[0]->GetNumModes();
275 const int order1 = m_base[1]->GetNumModes();
276
277 Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
278
279 if (multiplybyweights)
280 {
281 Array<OneD, NekDouble> tmp(nquad0 * nquad1 * nquad2);
282
283 MultiplyByQuadratureMetric(inarray, tmp);
284
286 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
287 tmp, outarray, wsp, true, true, true);
288 }
289 else
290 {
292 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
293 inarray, outarray, wsp, true, true, true);
294 }
295}
296
297/**
298 * @brief Calculates the inner product \f$ I_{pqr} = (u,
299 * \partial_{x_i} \phi_{pqr}) \f$.
300 *
301 * The derivative of the basis functions is performed using the chain
302 * rule in order to incorporate the geometric factors. Assuming that
303 * the basis functions are a tensor product
304 * \f$\phi_{pqr}(\eta_1,\eta_2,\eta_3) =
305 * \phi_1(\eta_1)\phi_2(\eta_2)\phi_3(\eta_3)\f$, this yields the
306 * result
307 *
308 * \f[
309 * I_{pqr} = \sum_{j=1}^3 \left(u, \frac{\partial u}{\partial \eta_j}
310 * \frac{\partial \eta_j}{\partial x_i}\right)
311 * \f]
312 *
313 * In the tetrahedral element, we must also incorporate a second set
314 * of geometric factors which incorporate the collapsed co-ordinate
315 * system, so that
316 *
317 * \f[ \frac{\partial\eta_j}{\partial x_i} = \sum_{k=1}^3
318 * \frac{\partial\eta_j}{\partial\xi_k}\frac{\partial\xi_k}{\partial
319 * x_i} \f]
320 *
321 * These derivatives can be found on p152 of Sherwin & Karniadakis.
322 *
323 * @param dir Direction in which to take the derivative.
324 * @param inarray The function \f$ u \f$.
325 * @param outarray Value of the inner product.
326 */
328 const int dir, const Array<OneD, const NekDouble> &inarray,
329 Array<OneD, NekDouble> &outarray)
330{
331 v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
332}
333
335 const int dir, const Array<OneD, const NekDouble> &inarray,
336 Array<OneD, NekDouble> &outarray)
337{
338 const int nquad0 = m_base[0]->GetNumPoints();
339 const int nquad1 = m_base[1]->GetNumPoints();
340 const int nquad2 = m_base[2]->GetNumPoints();
341 const int order0 = m_base[0]->GetNumModes();
342 const int order1 = m_base[1]->GetNumModes();
343 const int nqtot = nquad0 * nquad1 * nquad2;
344
345 Array<OneD, NekDouble> tmp1(nqtot);
346 Array<OneD, NekDouble> tmp2(nqtot);
347 Array<OneD, NekDouble> tmp3(nqtot);
348 Array<OneD, NekDouble> tmp4(nqtot);
350 Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
351
352 MultiplyByQuadratureMetric(inarray, tmp1);
353
355 tmp2D[0] = tmp2;
356 tmp2D[1] = tmp3;
357 tmp2D[2] = tmp4;
358
360
361 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
362 m_base[2]->GetBdata(), tmp2, outarray, wsp,
363 true, true, true);
364
365 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
366 m_base[2]->GetBdata(), tmp3, tmp6, wsp, true,
367 true, true);
368
369 Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
370
371 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
372 m_base[2]->GetDbdata(), tmp4, tmp6, wsp, true,
373 true, true);
374
375 Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
376}
377
379 const int dir, const Array<OneD, const NekDouble> &inarray,
381{
382 const int nquad0 = m_base[0]->GetNumPoints();
383 const int nquad1 = m_base[1]->GetNumPoints();
384 const int nquad2 = m_base[2]->GetNumPoints();
385 const int nqtot = nquad0 * nquad1 * nquad2;
386
387 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
388 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
389
390 Array<OneD, NekDouble> tmp1(nqtot);
391
392 Array<OneD, NekDouble> tmp2 = outarray[0];
393 Array<OneD, NekDouble> tmp3 = outarray[1];
394 Array<OneD, NekDouble> tmp4 = outarray[2];
395
397 m_metricinfo->GetDerivFactors(GetPointsKeys());
398
399 Vmath::Vcopy(nqtot, inarray, 1, tmp1, 1); // Dir3 metric
400
401 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
402 {
403 Vmath::Vmul(nqtot, &df[3 * dir][0], 1, tmp1.data(), 1, tmp2.data(), 1);
404 Vmath::Vmul(nqtot, &df[3 * dir + 1][0], 1, tmp1.data(), 1, tmp3.data(),
405 1);
406 Vmath::Vmul(nqtot, &df[3 * dir + 2][0], 1, tmp1.data(), 1, tmp4.data(),
407 1);
408 }
409 else
410 {
411 Vmath::Smul(nqtot, df[3 * dir][0], tmp1.data(), 1, tmp2.data(), 1);
412 Vmath::Smul(nqtot, df[3 * dir + 1][0], tmp1.data(), 1, tmp3.data(), 1);
413 Vmath::Smul(nqtot, df[3 * dir + 2][0], tmp1.data(), 1, tmp4.data(), 1);
414 }
415
416 int cnt = 0;
417 int i, j;
418
419 NekDouble g0, g2, g02;
420 for (int k = 0; k < nquad2; ++k)
421 {
422 g2 = 2.0 / (1.0 - z2[k]);
423
424 for (j = 0; j < nquad1; ++j)
425 {
426 for (i = 0; i < nquad0; ++i, ++cnt)
427 {
428 g0 = 0.5 * (1.0 + z0[i]);
429 g02 = g0 * g2;
430 tmp2[cnt] = g2 * tmp2[cnt] + g02 * tmp4[cnt];
431 }
432 }
433 }
434}
435
436//---------------------------------------
437// Evaluation functions
438//---------------------------------------
439
441{
443 m_base[0]->GetBasisKey(), m_base[1]->GetBasisKey(),
444 m_base[2]->GetBasisKey());
445}
446
448{
450 m_base[0]->GetPointsKey());
452 m_base[1]->GetPointsKey());
454 m_base[2]->GetPointsKey());
455
457 bkey0, bkey1, bkey2);
458}
459
460/**
461 * @brief Get the coordinates #coords at the local coordinates
462 * #Lcoords.
463 */
466{
467 int i;
468
469 ASSERTL1(Lcoords[0] <= -1.0 && Lcoords[0] >= 1.0 && Lcoords[1] <= -1.0 &&
470 Lcoords[1] >= 1.0 && Lcoords[2] <= -1.0 && Lcoords[2] >= 1.0,
471 "Local coordinates are not in region [-1,1]");
472
473 m_geom->FillGeom();
474
475 for (i = 0; i < m_geom->GetCoordim(); ++i)
476 {
477 coords[i] = m_geom->GetCoord(i, Lcoords);
478 }
479}
480
482 Array<OneD, NekDouble> &coords_1,
483 Array<OneD, NekDouble> &coords_2)
484{
485 Expansion::v_GetCoords(coords_0, coords_1, coords_2);
486}
487
488/**
489 * Given the local cartesian coordinate \a Lcoord evaluate the
490 * value of physvals at this point by calling through to the
491 * StdExpansion method
492 */
494 const Array<OneD, const NekDouble> &Lcoord,
495 const Array<OneD, const NekDouble> &physvals)
496{
497 // Evaluate point in local (eta) coordinates.
498 return StdExpansion3D::v_PhysEvaluate(Lcoord, physvals);
499}
500
502 const Array<OneD, const NekDouble> &physvals)
503{
504 Array<OneD, NekDouble> Lcoord(3);
505
506 ASSERTL0(m_geom, "m_geom not defined");
507
508 m_geom->GetLocCoords(coord, Lcoord);
509
510 return StdExpansion3D::v_PhysEvaluate(Lcoord, physvals);
511}
512
514 const Array<OneD, NekDouble> &coord,
515 const Array<OneD, const NekDouble> &inarray,
516 std::array<NekDouble, 3> &firstOrderDerivs)
517{
518 Array<OneD, NekDouble> Lcoord(3);
519 ASSERTL0(m_geom, "m_geom not defined");
520 m_geom->GetLocCoords(coord, Lcoord);
521 return StdPrismExp::v_PhysEvalFirstDeriv(Lcoord, inarray, firstOrderDerivs);
522}
523
524//---------------------------------------
525// Helper functions
526//---------------------------------------
527
529 const NekDouble *data, const std::vector<unsigned int> &nummodes,
530 const int mode_offset, NekDouble *coeffs,
531 [[maybe_unused]] std::vector<LibUtilities::BasisType> &fromType)
532{
533 int data_order0 = nummodes[mode_offset];
534 int fillorder0 = min(m_base[0]->GetNumModes(), data_order0);
535 int data_order1 = nummodes[mode_offset + 1];
536 int order1 = m_base[1]->GetNumModes();
537 int fillorder1 = min(order1, data_order1);
538 int data_order2 = nummodes[mode_offset + 2];
539 int order2 = m_base[2]->GetNumModes();
540 int fillorder2 = min(order2, data_order2);
541
542 switch (m_base[0]->GetBasisType())
543 {
545 {
546 int i, j;
547 int cnt = 0;
548 int cnt1 = 0;
549
551 "Extraction routine not set up for this basis");
553 "Extraction routine not set up for this basis");
554
555 Vmath::Zero(m_ncoeffs, coeffs, 1);
556 for (j = 0; j < fillorder0; ++j)
557 {
558 for (i = 0; i < fillorder1; ++i)
559 {
560 Vmath::Vcopy(fillorder2 - j, &data[cnt], 1, &coeffs[cnt1],
561 1);
562 cnt += data_order2 - j;
563 cnt1 += order2 - j;
564 }
565
566 // count out data for j iteration
567 for (i = fillorder1; i < data_order1; ++i)
568 {
569 cnt += data_order2 - j;
570 }
571
572 for (i = fillorder1; i < order1; ++i)
573 {
574 cnt1 += order2 - j;
575 }
576 }
577 }
578 break;
579 default:
580 ASSERTL0(false, "basis is either not set up or not "
581 "hierarchicial");
582 }
583}
584
585void PrismExp::v_GetTracePhysMap(const int face, Array<OneD, int> &outarray)
586{
587 int nquad0 = m_base[0]->GetNumPoints();
588 int nquad1 = m_base[1]->GetNumPoints();
589 int nquad2 = m_base[2]->GetNumPoints();
590 int nq0 = 0;
591 int nq1 = 0;
592
593 switch (face)
594 {
595 case 0:
596 nq0 = nquad0;
597 nq1 = nquad1;
598 if (outarray.size() != nq0 * nq1)
599 {
600 outarray = Array<OneD, int>(nq0 * nq1);
601 }
602
603 // Directions A and B positive
604 for (int i = 0; i < nquad0 * nquad1; ++i)
605 {
606 outarray[i] = i;
607 }
608 break;
609 case 1:
610
611 nq0 = nquad0;
612 nq1 = nquad2;
613 if (outarray.size() != nq0 * nq1)
614 {
615 outarray = Array<OneD, int>(nq0 * nq1);
616 }
617
618 // Direction A and B positive
619 for (int k = 0; k < nquad2; k++)
620 {
621 for (int i = 0; i < nquad0; ++i)
622 {
623 outarray[k * nquad0 + i] = (nquad0 * nquad1 * k) + i;
624 }
625 }
626
627 break;
628 case 2:
629
630 nq0 = nquad1;
631 nq1 = nquad2;
632 if (outarray.size() != nq0 * nq1)
633 {
634 outarray = Array<OneD, int>(nq0 * nq1);
635 }
636
637 // Directions A and B positive
638 for (int j = 0; j < nquad1 * nquad2; ++j)
639 {
640 outarray[j] = nquad0 - 1 + j * nquad0;
641 }
642 break;
643 case 3:
644 nq0 = nquad0;
645 nq1 = nquad2;
646 if (outarray.size() != nq0 * nq1)
647 {
648 outarray = Array<OneD, int>(nq0 * nq1);
649 }
650
651 // Direction A and B positive
652 for (int k = 0; k < nquad2; k++)
653 {
654 for (int i = 0; i < nquad0; ++i)
655 {
656 outarray[k * nquad0 + i] =
657 nquad0 * (nquad1 - 1) + (nquad0 * nquad1 * k) + i;
658 }
659 }
660 break;
661 case 4:
662
663 nq0 = nquad1;
664 nq1 = nquad2;
665 if (outarray.size() != nq0 * nq1)
666 {
667 outarray = Array<OneD, int>(nq0 * nq1);
668 }
669
670 // Directions A and B positive
671 for (int j = 0; j < nquad1 * nquad2; ++j)
672 {
673 outarray[j] = j * nquad0;
674 }
675 break;
676 default:
677 ASSERTL0(false, "face value (> 4) is out of range");
678 break;
679 }
680}
681
682/** \brief Get the normals along specficied face
683 * Get the face normals interplated to a points0 x points 0
684 * type distribution
685 **/
687{
688 const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
690
692 for (int i = 0; i < ptsKeys.size(); ++i)
693 {
694 // Need at least 2 points for computing normals
695 if (ptsKeys[i].GetNumPoints() == 1)
696 {
697 LibUtilities::PointsKey pKey(2, ptsKeys[i].GetPointsType());
698 ptsKeys[i] = pKey;
699 }
700 }
701
702 SpatialDomains::GeomType type = geomFactors->GetGtype();
704 geomFactors->GetDerivFactors(ptsKeys);
705 const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
706
707 int nq0 = ptsKeys[0].GetNumPoints();
708 int nq1 = ptsKeys[1].GetNumPoints();
709 int nq2 = ptsKeys[2].GetNumPoints();
710 int nq01 = nq0 * nq1;
711 int nqtot;
712
713 LibUtilities::BasisKey tobasis0 = GetTraceBasisKey(face, 0);
714 LibUtilities::BasisKey tobasis1 = GetTraceBasisKey(face, 1);
715
716 // Number of quadrature points in face expansion.
717 int nq_face = tobasis0.GetNumPoints() * tobasis1.GetNumPoints();
718
719 int vCoordDim = GetCoordim();
720 int i;
721
724 for (i = 0; i < vCoordDim; ++i)
725 {
726 normal[i] = Array<OneD, NekDouble>(nq_face);
727 }
728
729 size_t nqb = nq_face;
730 size_t nbnd = face;
733
734 // Regular geometry case
735 if (type == SpatialDomains::eRegular ||
737 {
738 NekDouble fac;
739 // Set up normals
740 switch (face)
741 {
742 case 0:
743 {
744 for (i = 0; i < vCoordDim; ++i)
745 {
746 normal[i][0] = -df[3 * i + 2][0];
747 ;
748 }
749 break;
750 }
751 case 1:
752 {
753 for (i = 0; i < vCoordDim; ++i)
754 {
755 normal[i][0] = -df[3 * i + 1][0];
756 }
757 break;
758 }
759 case 2:
760 {
761 for (i = 0; i < vCoordDim; ++i)
762 {
763 normal[i][0] = df[3 * i][0] + df[3 * i + 2][0];
764 }
765 break;
766 }
767 case 3:
768 {
769 for (i = 0; i < vCoordDim; ++i)
770 {
771 normal[i][0] = df[3 * i + 1][0];
772 }
773 break;
774 }
775 case 4:
776 {
777 for (i = 0; i < vCoordDim; ++i)
778 {
779 normal[i][0] = -df[3 * i][0];
780 }
781 break;
782 }
783 default:
784 ASSERTL0(false, "face is out of range (face < 4)");
785 }
786
787 // Normalise resulting vector.
788 fac = 0.0;
789 for (i = 0; i < vCoordDim; ++i)
790 {
791 fac += normal[i][0] * normal[i][0];
792 }
793 fac = 1.0 / sqrt(fac);
794
795 Vmath::Fill(nqb, fac, length, 1);
796
797 for (i = 0; i < vCoordDim; ++i)
798 {
799 Vmath::Fill(nq_face, fac * normal[i][0], normal[i], 1);
800 }
801 }
802 else
803 {
804 // Set up deformed normals.
805 int j, k;
806
807 // Determine number of quadrature points on the face of 3D elmt
808 if (face == 0)
809 {
810 nqtot = nq0 * nq1;
811 }
812 else if (face == 1 || face == 3)
813 {
814 nqtot = nq0 * nq2;
815 }
816 else
817 {
818 nqtot = nq1 * nq2;
819 }
820
823
824 Array<OneD, NekDouble> faceJac(nqtot);
825 Array<OneD, NekDouble> normals(vCoordDim * nqtot, 0.0);
826
827 // Extract Jacobian along face and recover local derivatives
828 // (dx/dr) for polynomial interpolation by multiplying m_gmat by
829 // jacobian
830 switch (face)
831 {
832 case 0:
833 {
834 for (j = 0; j < nq01; ++j)
835 {
836 normals[j] = -df[2][j] * jac[j];
837 normals[nqtot + j] = -df[5][j] * jac[j];
838 normals[2 * nqtot + j] = -df[8][j] * jac[j];
839 faceJac[j] = jac[j];
840 }
841
842 points0 = ptsKeys[0];
843 points1 = ptsKeys[1];
844 break;
845 }
846
847 case 1:
848 {
849 for (j = 0; j < nq0; ++j)
850 {
851 for (k = 0; k < nq2; ++k)
852 {
853 int tmp = j + nq01 * k;
854 normals[j + k * nq0] = -df[1][tmp] * jac[tmp];
855 normals[nqtot + j + k * nq0] = -df[4][tmp] * jac[tmp];
856 normals[2 * nqtot + j + k * nq0] =
857 -df[7][tmp] * jac[tmp];
858 faceJac[j + k * nq0] = jac[tmp];
859 }
860 }
861
862 points0 = ptsKeys[0];
863 points1 = ptsKeys[2];
864 break;
865 }
866
867 case 2:
868 {
869 for (j = 0; j < nq1; ++j)
870 {
871 for (k = 0; k < nq2; ++k)
872 {
873 int tmp = nq0 - 1 + nq0 * j + nq01 * k;
874 normals[j + k * nq1] =
875 (df[0][tmp] + df[2][tmp]) * jac[tmp];
876 normals[nqtot + j + k * nq1] =
877 (df[3][tmp] + df[5][tmp]) * jac[tmp];
878 normals[2 * nqtot + j + k * nq1] =
879 (df[6][tmp] + df[8][tmp]) * jac[tmp];
880 faceJac[j + k * nq1] = jac[tmp];
881 }
882 }
883
884 points0 = ptsKeys[1];
885 points1 = ptsKeys[2];
886 break;
887 }
888
889 case 3:
890 {
891 for (j = 0; j < nq0; ++j)
892 {
893 for (k = 0; k < nq2; ++k)
894 {
895 int tmp = nq0 * (nq1 - 1) + j + nq01 * k;
896 normals[j + k * nq0] = df[1][tmp] * jac[tmp];
897 normals[nqtot + j + k * nq0] = df[4][tmp] * jac[tmp];
898 normals[2 * nqtot + j + k * nq0] =
899 df[7][tmp] * jac[tmp];
900 faceJac[j + k * nq0] = jac[tmp];
901 }
902 }
903
904 points0 = ptsKeys[0];
905 points1 = ptsKeys[2];
906 break;
907 }
908
909 case 4:
910 {
911 for (j = 0; j < nq1; ++j)
912 {
913 for (k = 0; k < nq2; ++k)
914 {
915 int tmp = j * nq0 + nq01 * k;
916 normals[j + k * nq1] = -df[0][tmp] * jac[tmp];
917 normals[nqtot + j + k * nq1] = -df[3][tmp] * jac[tmp];
918 normals[2 * nqtot + j + k * nq1] =
919 -df[6][tmp] * jac[tmp];
920 faceJac[j + k * nq1] = jac[tmp];
921 }
922 }
923
924 points0 = ptsKeys[1];
925 points1 = ptsKeys[2];
926 break;
927 }
928
929 default:
930 ASSERTL0(false, "face is out of range (face < 4)");
931 }
932
933 Array<OneD, NekDouble> work(nq_face, 0.0);
934 // Interpolate Jacobian and invert
935 LibUtilities::Interp2D(points0, points1, faceJac,
936 tobasis0.GetPointsKey(), tobasis1.GetPointsKey(),
937 work);
938 Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
939
940 // Interpolate normal and multiply by inverse Jacobian.
941 for (i = 0; i < vCoordDim; ++i)
942 {
943 LibUtilities::Interp2D(points0, points1, &normals[i * nqtot],
944 tobasis0.GetPointsKey(),
945 tobasis1.GetPointsKey(), &normal[i][0]);
946 Vmath::Vmul(nq_face, work, 1, normal[i], 1, normal[i], 1);
947 }
948
949 // Normalise to obtain unit normals.
950 Vmath::Zero(nq_face, work, 1);
951 for (i = 0; i < GetCoordim(); ++i)
952 {
953 Vmath::Vvtvp(nq_face, normal[i], 1, normal[i], 1, work, 1, work, 1);
954 }
955
956 Vmath::Vsqrt(nq_face, work, 1, work, 1);
957 Vmath::Sdiv(nq_face, 1.0, work, 1, work, 1);
958
959 Vmath::Vcopy(nqb, work, 1, length, 1);
960
961 for (i = 0; i < GetCoordim(); ++i)
962 {
963 Vmath::Vmul(nq_face, normal[i], 1, work, 1, normal[i], 1);
964 }
965 }
966}
967
969 Array<OneD, NekDouble> &outarray,
970 const StdRegions::StdMatrixKey &mkey)
971{
972 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
973}
974
976 Array<OneD, NekDouble> &outarray,
977 const StdRegions::StdMatrixKey &mkey)
978{
979 PrismExp::LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
980}
981
982void PrismExp::v_LaplacianMatrixOp(const int k1, const int k2,
983 const Array<OneD, const NekDouble> &inarray,
984 Array<OneD, NekDouble> &outarray,
985 const StdRegions::StdMatrixKey &mkey)
986{
987 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
988}
989
991 Array<OneD, NekDouble> &outarray,
992 const StdRegions::StdMatrixKey &mkey)
993{
994 PrismExp::v_HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
995}
996
998 const StdRegions::StdMatrixKey &mkey)
999{
1000 int nq = GetTotPoints();
1001
1002 // Calculate sqrt of the Jacobian
1004 Array<OneD, NekDouble> sqrt_jac(nq);
1005 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1006 {
1007 Vmath::Vsqrt(nq, jac, 1, sqrt_jac, 1);
1008 }
1009 else
1010 {
1011 Vmath::Fill(nq, sqrt(jac[0]), sqrt_jac, 1);
1012 }
1013
1014 // Multiply array by sqrt(Jac)
1015 Vmath::Vmul(nq, sqrt_jac, 1, array, 1, array, 1);
1016
1017 // Apply std region filter
1018 StdPrismExp::v_SVVLaplacianFilter(array, mkey);
1019
1020 // Divide by sqrt(Jac)
1021 Vmath::Vdiv(nq, array, 1, sqrt_jac, 1, array, 1);
1022}
1023
1024//---------------------------------------
1025// Matrix creation functions
1026//---------------------------------------
1027
1029{
1030 DNekMatSharedPtr returnval;
1031
1032 switch (mkey.GetMatrixType())
1033 {
1041 returnval = Expansion3D::v_GenMatrix(mkey);
1042 break;
1043 default:
1044 returnval = StdPrismExp::v_GenMatrix(mkey);
1045 break;
1046 }
1047
1048 return returnval;
1049}
1050
1052 const StdRegions::StdMatrixKey &mkey)
1053{
1054 LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
1055 LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
1056 LibUtilities::BasisKey bkey2 = m_base[2]->GetBasisKey();
1059
1060 return tmp->GetStdMatrix(mkey);
1061}
1062
1067
1069{
1070 m_matrixManager.DeleteObject(mkey);
1071}
1072
1078
1080{
1081 m_staticCondMatrixManager.DeleteObject(mkey);
1082}
1083
1084/**
1085 * @brief Calculate the Laplacian multiplication in a matrix-free
1086 * manner.
1087 *
1088 * This function is the kernel of the Laplacian matrix-free operator,
1089 * and is used in #v_HelmholtzMatrixOp_MatFree to determine the effect
1090 * of the Helmholtz operator in a similar fashion.
1091 *
1092 * The majority of the calculation is precisely the same as in the
1093 * hexahedral expansion; however the collapsed co-ordinate system must
1094 * be taken into account when constructing the geometric factors. How
1095 * this is done is detailed more exactly in the tetrahedral expansion.
1096 * On entry to this function, the input #inarray must be in its
1097 * backwards-transformed state (i.e. \f$\mathbf{u} =
1098 * \mathbf{B}\hat{\mathbf{u}}\f$). The output is in coefficient space.
1099 *
1100 * @see %TetExp::v_HelmholtzMatrixOp_MatFree
1101 */
1103 const Array<OneD, const NekDouble> &inarray,
1105{
1106 int nquad0 = m_base[0]->GetNumPoints();
1107 int nquad1 = m_base[1]->GetNumPoints();
1108 int nquad2 = m_base[2]->GetNumPoints();
1109 int nqtot = nquad0 * nquad1 * nquad2;
1110 int i;
1111
1112 // Set up temporary storage.
1113 Array<OneD, NekDouble> alloc(11 * nqtot, 0.0);
1114 Array<OneD, NekDouble> wsp1(alloc); // TensorDeriv 1
1115 Array<OneD, NekDouble> wsp2(alloc + 1 * nqtot); // TensorDeriv 2
1116 Array<OneD, NekDouble> wsp3(alloc + 2 * nqtot); // TensorDeriv 3
1117 Array<OneD, NekDouble> g0(alloc + 3 * nqtot); // g0
1118 Array<OneD, NekDouble> g1(alloc + 4 * nqtot); // g1
1119 Array<OneD, NekDouble> g2(alloc + 5 * nqtot); // g2
1120 Array<OneD, NekDouble> g3(alloc + 6 * nqtot); // g3
1121 Array<OneD, NekDouble> g4(alloc + 7 * nqtot); // g4
1122 Array<OneD, NekDouble> g5(alloc + 8 * nqtot); // g5
1123 Array<OneD, NekDouble> h0(alloc + 3 * nqtot); // h0 == g0
1124 Array<OneD, NekDouble> h1(alloc + 6 * nqtot); // h1 == g3
1125 Array<OneD, NekDouble> wsp4(alloc + 4 * nqtot); // wsp4 == g1
1126 Array<OneD, NekDouble> wsp5(alloc + 5 * nqtot); // wsp5 == g2
1127 Array<OneD, NekDouble> wsp6(alloc + 8 * nqtot); // wsp6 == g5
1128 Array<OneD, NekDouble> wsp7(alloc + 3 * nqtot); // wsp7 == g0
1129 Array<OneD, NekDouble> wsp8(alloc + 9 * nqtot); // wsp8
1130 Array<OneD, NekDouble> wsp9(alloc + 10 * nqtot); // wsp9
1131
1132 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
1133 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
1134 const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
1135 const Array<OneD, const NekDouble> &dbase0 = m_base[0]->GetDbdata();
1136 const Array<OneD, const NekDouble> &dbase1 = m_base[1]->GetDbdata();
1137 const Array<OneD, const NekDouble> &dbase2 = m_base[2]->GetDbdata();
1138
1139 // Step 1. LAPLACIAN MATRIX OPERATION
1140 // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
1141 // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
1142 // wsp3 = du_dxi3 = D_xi3 * wsp0 = D_xi3 * u
1143 StdExpansion3D::PhysTensorDeriv(inarray, wsp1, wsp2, wsp3);
1144
1146 m_metricinfo->GetDerivFactors(GetPointsKeys());
1147 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
1148 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
1149
1150 // Step 2. Calculate the metric terms of the collapsed
1151 // coordinate transformation (Spencer's book P152)
1152 for (i = 0; i < nquad2; ++i)
1153 {
1154 Vmath::Fill(nquad0 * nquad1, 2.0 / (1.0 - z2[i]),
1155 &h0[0] + i * nquad0 * nquad1, 1);
1156 Vmath::Fill(nquad0 * nquad1, 2.0 / (1.0 - z2[i]),
1157 &h1[0] + i * nquad0 * nquad1, 1);
1158 }
1159 for (i = 0; i < nquad0; i++)
1160 {
1161 Blas::Dscal(nquad1 * nquad2, 0.5 * (1 + z0[i]), &h1[0] + i, nquad0);
1162 }
1163
1164 // Step 3. Construct combined metric terms for physical space to
1165 // collapsed coordinate system. Order of construction optimised
1166 // to minimise temporary storage
1167 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1168 {
1169 // wsp4 = d eta_1/d x_1
1170 Vmath::Vvtvvtp(nqtot, &df[0][0], 1, &h0[0], 1, &df[2][0], 1, &h1[0], 1,
1171 &wsp4[0], 1);
1172 // wsp5 = d eta_2/d x_1
1173 Vmath::Vvtvvtp(nqtot, &df[3][0], 1, &h0[0], 1, &df[5][0], 1, &h1[0], 1,
1174 &wsp5[0], 1);
1175 // wsp6 = d eta_3/d x_1d
1176 Vmath::Vvtvvtp(nqtot, &df[6][0], 1, &h0[0], 1, &df[8][0], 1, &h1[0], 1,
1177 &wsp6[0], 1);
1178
1179 // g0 (overwrites h0)
1180 Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0],
1181 1, &g0[0], 1);
1182 Vmath::Vvtvp(nqtot, &wsp6[0], 1, &wsp6[0], 1, &g0[0], 1, &g0[0], 1);
1183
1184 // g3 (overwrites h1)
1185 Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &wsp4[0], 1, &df[4][0], 1, &wsp5[0],
1186 1, &g3[0], 1);
1187 Vmath::Vvtvp(nqtot, &df[7][0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1188
1189 // g4
1190 Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp4[0], 1, &df[5][0], 1, &wsp5[0],
1191 1, &g4[0], 1);
1192 Vmath::Vvtvp(nqtot, &df[8][0], 1, &wsp6[0], 1, &g4[0], 1, &g4[0], 1);
1193
1194 // Overwrite wsp4/5/6 with g1/2/5
1195 // g1
1196 Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &df[1][0], 1, &df[4][0], 1,
1197 &df[4][0], 1, &g1[0], 1);
1198 Vmath::Vvtvp(nqtot, &df[7][0], 1, &df[7][0], 1, &g1[0], 1, &g1[0], 1);
1199
1200 // g2
1201 Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &df[2][0], 1, &df[5][0], 1,
1202 &df[5][0], 1, &g2[0], 1);
1203 Vmath::Vvtvp(nqtot, &df[8][0], 1, &df[8][0], 1, &g2[0], 1, &g2[0], 1);
1204
1205 // g5
1206 Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &df[2][0], 1, &df[4][0], 1,
1207 &df[5][0], 1, &g5[0], 1);
1208 Vmath::Vvtvp(nqtot, &df[7][0], 1, &df[8][0], 1, &g5[0], 1, &g5[0], 1);
1209 }
1210 else
1211 {
1212 // wsp4 = d eta_1/d x_1
1213 Vmath::Svtsvtp(nqtot, df[0][0], &h0[0], 1, df[2][0], &h1[0], 1,
1214 &wsp4[0], 1);
1215 // wsp5 = d eta_2/d x_1
1216 Vmath::Svtsvtp(nqtot, df[3][0], &h0[0], 1, df[5][0], &h1[0], 1,
1217 &wsp5[0], 1);
1218 // wsp6 = d eta_3/d x_1
1219 Vmath::Svtsvtp(nqtot, df[6][0], &h0[0], 1, df[8][0], &h1[0], 1,
1220 &wsp6[0], 1);
1221
1222 // g0 (overwrites h0)
1223 Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0],
1224 1, &g0[0], 1);
1225 Vmath::Vvtvp(nqtot, &wsp6[0], 1, &wsp6[0], 1, &g0[0], 1, &g0[0], 1);
1226
1227 // g3 (overwrites h1)
1228 Vmath::Svtsvtp(nqtot, df[1][0], &wsp4[0], 1, df[4][0], &wsp5[0], 1,
1229 &g3[0], 1);
1230 Vmath::Svtvp(nqtot, df[7][0], &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1231
1232 // g4
1233 Vmath::Svtsvtp(nqtot, df[2][0], &wsp4[0], 1, df[5][0], &wsp5[0], 1,
1234 &g4[0], 1);
1235 Vmath::Svtvp(nqtot, df[8][0], &wsp6[0], 1, &g4[0], 1, &g4[0], 1);
1236
1237 // Overwrite wsp4/5/6 with g1/2/5
1238 // g1
1239 Vmath::Fill(nqtot,
1240 df[1][0] * df[1][0] + df[4][0] * df[4][0] +
1241 df[7][0] * df[7][0],
1242 &g1[0], 1);
1243
1244 // g2
1245 Vmath::Fill(nqtot,
1246 df[2][0] * df[2][0] + df[5][0] * df[5][0] +
1247 df[8][0] * df[8][0],
1248 &g2[0], 1);
1249
1250 // g5
1251 Vmath::Fill(nqtot,
1252 df[1][0] * df[2][0] + df[4][0] * df[5][0] +
1253 df[7][0] * df[8][0],
1254 &g5[0], 1);
1255 }
1256 // Compute component derivatives into wsp7, 8, 9 (wsp7 overwrites
1257 // g0).
1258 Vmath::Vvtvvtp(nqtot, &g0[0], 1, &wsp1[0], 1, &g3[0], 1, &wsp2[0], 1,
1259 &wsp7[0], 1);
1260 Vmath::Vvtvp(nqtot, &g4[0], 1, &wsp3[0], 1, &wsp7[0], 1, &wsp7[0], 1);
1261 Vmath::Vvtvvtp(nqtot, &g1[0], 1, &wsp2[0], 1, &g3[0], 1, &wsp1[0], 1,
1262 &wsp8[0], 1);
1263 Vmath::Vvtvp(nqtot, &g5[0], 1, &wsp3[0], 1, &wsp8[0], 1, &wsp8[0], 1);
1264 Vmath::Vvtvvtp(nqtot, &g2[0], 1, &wsp3[0], 1, &g4[0], 1, &wsp1[0], 1,
1265 &wsp9[0], 1);
1266 Vmath::Vvtvp(nqtot, &g5[0], 1, &wsp2[0], 1, &wsp9[0], 1, &wsp9[0], 1);
1267
1268 // Step 4.
1269 // Multiply by quadrature metric
1270 MultiplyByQuadratureMetric(wsp7, wsp7);
1271 MultiplyByQuadratureMetric(wsp8, wsp8);
1272 MultiplyByQuadratureMetric(wsp9, wsp9);
1273
1274 // Perform inner product w.r.t derivative bases.
1275 IProductWRTBase_SumFacKernel(dbase0, base1, base2, wsp7, wsp1, wsp, false,
1276 true, true);
1277 IProductWRTBase_SumFacKernel(base0, dbase1, base2, wsp8, wsp2, wsp, true,
1278 false, true);
1279 IProductWRTBase_SumFacKernel(base0, base1, dbase2, wsp9, outarray, wsp,
1280 true, true, false);
1281
1282 // Step 5.
1283 // Sum contributions from wsp1, wsp2 and outarray.
1284 Vmath::Vadd(m_ncoeffs, wsp1.data(), 1, outarray.data(), 1, outarray.data(),
1285 1);
1286 Vmath::Vadd(m_ncoeffs, wsp2.data(), 1, outarray.data(), 1, outarray.data(),
1287 1);
1288}
1289
1291 Array<OneD, int> &conn, [[maybe_unused]] bool oldstandard)
1292{
1293 int np0 = m_base[0]->GetNumPoints();
1294 int np1 = m_base[1]->GetNumPoints();
1295 int np2 = m_base[2]->GetNumPoints();
1296 int np = max(np0, max(np1, np2));
1297 Array<OneD, int> prismpt(6);
1298 bool standard = true;
1299
1300 int vid0 = m_geom->GetVid(0);
1301 int vid1 = m_geom->GetVid(1);
1302 int vid2 = m_geom->GetVid(4);
1303 int rotate = 0;
1304
1305 // sort out prism rotation according to
1306 if ((vid2 < vid1) && (vid2 < vid0)) // top triangle vertex is lowest id
1307 {
1308 rotate = 0;
1309 if (vid0 > vid1)
1310 {
1311 standard = false; // reverse base direction
1312 }
1313 }
1314 else if ((vid1 < vid2) && (vid1 < vid0))
1315 {
1316 rotate = 1;
1317 if (vid2 > vid0)
1318 {
1319 standard = false; // reverse base direction
1320 }
1321 }
1322 else if ((vid0 < vid2) && (vid0 < vid1))
1323 {
1324 rotate = 2;
1325 if (vid1 > vid2)
1326 {
1327 standard = false; // reverse base direction
1328 }
1329 }
1330
1331 conn = Array<OneD, int>(12 * (np - 1) * (np - 1) * (np - 1));
1332
1333 int row = 0;
1334 int rowp1 = 0;
1335 int plane = 0;
1336 int row1 = 0;
1337 int row1p1 = 0;
1338 int planep1 = 0;
1339 int cnt = 0;
1340
1341 Array<OneD, int> rot(3);
1342
1343 rot[0] = (0 + rotate) % 3;
1344 rot[1] = (1 + rotate) % 3;
1345 rot[2] = (2 + rotate) % 3;
1346
1347 // lower diagonal along 1-3 on base
1348 for (int i = 0; i < np - 1; ++i)
1349 {
1350 planep1 += (np - i) * np;
1351 row = 0; // current plane row offset
1352 rowp1 = 0; // current plane row plus one offset
1353 row1 = 0; // next plane row offset
1354 row1p1 = 0; // nex plane row plus one offset
1355 if (standard == false)
1356 {
1357 for (int j = 0; j < np - 1; ++j)
1358 {
1359 rowp1 += np - i;
1360 row1p1 += np - i - 1;
1361 for (int k = 0; k < np - i - 2; ++k)
1362 {
1363 // bottom prism block
1364 prismpt[rot[0]] = plane + row + k;
1365 prismpt[rot[1]] = plane + row + k + 1;
1366 prismpt[rot[2]] = planep1 + row1 + k;
1367
1368 prismpt[3 + rot[0]] = plane + rowp1 + k;
1369 prismpt[3 + rot[1]] = plane + rowp1 + k + 1;
1370 prismpt[3 + rot[2]] = planep1 + row1p1 + k;
1371
1372 conn[cnt++] = prismpt[0];
1373 conn[cnt++] = prismpt[1];
1374 conn[cnt++] = prismpt[3];
1375 conn[cnt++] = prismpt[2];
1376
1377 conn[cnt++] = prismpt[5];
1378 conn[cnt++] = prismpt[2];
1379 conn[cnt++] = prismpt[3];
1380 conn[cnt++] = prismpt[4];
1381
1382 conn[cnt++] = prismpt[3];
1383 conn[cnt++] = prismpt[1];
1384 conn[cnt++] = prismpt[4];
1385 conn[cnt++] = prismpt[2];
1386
1387 // upper prism block.
1388 prismpt[rot[0]] = planep1 + row1 + k + 1;
1389 prismpt[rot[1]] = planep1 + row1 + k;
1390 prismpt[rot[2]] = plane + row + k + 1;
1391
1392 prismpt[3 + rot[0]] = planep1 + row1p1 + k + 1;
1393 prismpt[3 + rot[1]] = planep1 + row1p1 + k;
1394 prismpt[3 + rot[2]] = plane + rowp1 + k + 1;
1395
1396 conn[cnt++] = prismpt[0];
1397 conn[cnt++] = prismpt[1];
1398 conn[cnt++] = prismpt[2];
1399 conn[cnt++] = prismpt[5];
1400
1401 conn[cnt++] = prismpt[5];
1402 conn[cnt++] = prismpt[0];
1403 conn[cnt++] = prismpt[4];
1404 conn[cnt++] = prismpt[1];
1405
1406 conn[cnt++] = prismpt[3];
1407 conn[cnt++] = prismpt[4];
1408 conn[cnt++] = prismpt[0];
1409 conn[cnt++] = prismpt[5];
1410 }
1411
1412 // bottom prism block
1413 prismpt[rot[0]] = plane + row + np - i - 2;
1414 prismpt[rot[1]] = plane + row + np - i - 1;
1415 prismpt[rot[2]] = planep1 + row1 + np - i - 2;
1416
1417 prismpt[3 + rot[0]] = plane + rowp1 + np - i - 2;
1418 prismpt[3 + rot[1]] = plane + rowp1 + np - i - 1;
1419 prismpt[3 + rot[2]] = planep1 + row1p1 + np - i - 2;
1420
1421 conn[cnt++] = prismpt[0];
1422 conn[cnt++] = prismpt[1];
1423 conn[cnt++] = prismpt[3];
1424 conn[cnt++] = prismpt[2];
1425
1426 conn[cnt++] = prismpt[5];
1427 conn[cnt++] = prismpt[2];
1428 conn[cnt++] = prismpt[3];
1429 conn[cnt++] = prismpt[4];
1430
1431 conn[cnt++] = prismpt[3];
1432 conn[cnt++] = prismpt[1];
1433 conn[cnt++] = prismpt[4];
1434 conn[cnt++] = prismpt[2];
1435
1436 row += np - i;
1437 row1 += np - i - 1;
1438 }
1439 }
1440 else
1441 { // lower diagonal along 0-4 on base
1442 for (int j = 0; j < np - 1; ++j)
1443 {
1444 rowp1 += np - i;
1445 row1p1 += np - i - 1;
1446 for (int k = 0; k < np - i - 2; ++k)
1447 {
1448 // bottom prism block
1449 prismpt[rot[0]] = plane + row + k;
1450 prismpt[rot[1]] = plane + row + k + 1;
1451 prismpt[rot[2]] = planep1 + row1 + k;
1452
1453 prismpt[3 + rot[0]] = plane + rowp1 + k;
1454 prismpt[3 + rot[1]] = plane + rowp1 + k + 1;
1455 prismpt[3 + rot[2]] = planep1 + row1p1 + k;
1456
1457 conn[cnt++] = prismpt[0];
1458 conn[cnt++] = prismpt[1];
1459 conn[cnt++] = prismpt[4];
1460 conn[cnt++] = prismpt[2];
1461
1462 conn[cnt++] = prismpt[4];
1463 conn[cnt++] = prismpt[3];
1464 conn[cnt++] = prismpt[0];
1465 conn[cnt++] = prismpt[2];
1466
1467 conn[cnt++] = prismpt[3];
1468 conn[cnt++] = prismpt[4];
1469 conn[cnt++] = prismpt[5];
1470 conn[cnt++] = prismpt[2];
1471
1472 // upper prism block.
1473 prismpt[rot[0]] = planep1 + row1 + k + 1;
1474 prismpt[rot[1]] = planep1 + row1 + k;
1475 prismpt[rot[2]] = plane + row + k + 1;
1476
1477 prismpt[3 + rot[0]] = planep1 + row1p1 + k + 1;
1478 prismpt[3 + rot[1]] = planep1 + row1p1 + k;
1479 prismpt[3 + rot[2]] = plane + rowp1 + k + 1;
1480
1481 conn[cnt++] = prismpt[0];
1482 conn[cnt++] = prismpt[2];
1483 conn[cnt++] = prismpt[1];
1484 conn[cnt++] = prismpt[5];
1485
1486 conn[cnt++] = prismpt[3];
1487 conn[cnt++] = prismpt[5];
1488 conn[cnt++] = prismpt[0];
1489 conn[cnt++] = prismpt[1];
1490
1491 conn[cnt++] = prismpt[5];
1492 conn[cnt++] = prismpt[3];
1493 conn[cnt++] = prismpt[4];
1494 conn[cnt++] = prismpt[1];
1495 }
1496
1497 // bottom prism block
1498 prismpt[rot[0]] = plane + row + np - i - 2;
1499 prismpt[rot[1]] = plane + row + np - i - 1;
1500 prismpt[rot[2]] = planep1 + row1 + np - i - 2;
1501
1502 prismpt[3 + rot[0]] = plane + rowp1 + np - i - 2;
1503 prismpt[3 + rot[1]] = plane + rowp1 + np - i - 1;
1504 prismpt[3 + rot[2]] = planep1 + row1p1 + np - i - 2;
1505
1506 conn[cnt++] = prismpt[0];
1507 conn[cnt++] = prismpt[1];
1508 conn[cnt++] = prismpt[4];
1509 conn[cnt++] = prismpt[2];
1510
1511 conn[cnt++] = prismpt[4];
1512 conn[cnt++] = prismpt[3];
1513 conn[cnt++] = prismpt[0];
1514 conn[cnt++] = prismpt[2];
1515
1516 conn[cnt++] = prismpt[3];
1517 conn[cnt++] = prismpt[4];
1518 conn[cnt++] = prismpt[5];
1519 conn[cnt++] = prismpt[2];
1520
1521 row += np - i;
1522 row1 += np - i - 1;
1523 }
1524 }
1525 plane += (np - i) * np;
1526 }
1527}
1528
1529/** @brief: This method gets all of the factors which are
1530 required as part of the Gradient Jump Penalty
1531 stabilisation and involves the product of the normal and
1532 geometric factors along the element trace.
1533*/
1535 Array<OneD, Array<OneD, NekDouble>> &d0factors,
1536 Array<OneD, Array<OneD, NekDouble>> &d1factors,
1537 Array<OneD, Array<OneD, NekDouble>> &d2factors)
1538{
1539 int nquad0 = GetNumPoints(0);
1540 int nquad1 = GetNumPoints(1);
1541 int nquad2 = GetNumPoints(2);
1542
1544 m_metricinfo->GetDerivFactors(GetPointsKeys());
1545
1546 if (d0factors.size() != 5)
1547 {
1548 d0factors = Array<OneD, Array<OneD, NekDouble>>(5);
1549 d1factors = Array<OneD, Array<OneD, NekDouble>>(5);
1550 d2factors = Array<OneD, Array<OneD, NekDouble>>(5);
1551 }
1552
1553 if (d0factors[0].size() != nquad0 * nquad1)
1554 {
1555 d0factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1556 d1factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1557 d2factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1558 }
1559
1560 if (d0factors[1].size() != nquad0 * nquad2)
1561 {
1562 d0factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1563 d0factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1564 d1factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1565 d1factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1566 d2factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1567 d2factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1568 }
1569
1570 if (d0factors[2].size() != nquad1 * nquad2)
1571 {
1572 d0factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1573 d0factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1574 d1factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1575 d1factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1576 d2factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1577 d2factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1578 }
1579
1580 // Outwards normals
1582 GetTraceNormal(0);
1584 GetTraceNormal(1);
1586 GetTraceNormal(2);
1588 GetTraceNormal(3);
1590 GetTraceNormal(4);
1591
1592 int ncoords = normal_0.size();
1593
1594 // first gather together standard cartesian inner products
1595 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1596 {
1597 // face 0
1598 for (int i = 0; i < nquad0 * nquad1; ++i)
1599 {
1600 d0factors[0][i] = df[0][i] * normal_0[0][i];
1601 d1factors[0][i] = df[1][i] * normal_0[0][i];
1602 d2factors[0][i] = df[2][i] * normal_0[0][i];
1603 }
1604
1605 for (int n = 1; n < ncoords; ++n)
1606 {
1607 for (int i = 0; i < nquad0 * nquad1; ++i)
1608 {
1609 d0factors[0][i] += df[3 * n][i] * normal_0[n][i];
1610 d1factors[0][i] += df[3 * n + 1][i] * normal_0[n][i];
1611 d2factors[0][i] += df[3 * n + 2][i] * normal_0[n][i];
1612 }
1613 }
1614
1615 // faces 1 and 3
1616 for (int j = 0; j < nquad2; ++j)
1617 {
1618 for (int i = 0; i < nquad0; ++i)
1619 {
1620 d0factors[1][j * nquad0 + i] = df[0][j * nquad0 * nquad1 + i] *
1621 normal_1[0][j * nquad0 + i];
1622 d1factors[1][j * nquad0 + i] = df[1][j * nquad0 * nquad1 + i] *
1623 normal_1[0][j * nquad0 + i];
1624 d2factors[1][j * nquad0 + i] = df[2][j * nquad0 * nquad1 + i] *
1625 normal_1[0][j * nquad0 + i];
1626
1627 d0factors[3][j * nquad0 + i] =
1628 df[0][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1629 normal_3[0][j * nquad0 + i];
1630 d1factors[3][j * nquad0 + i] =
1631 df[1][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1632 normal_3[0][j * nquad0 + i];
1633 d2factors[3][j * nquad0 + i] =
1634 df[2][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1635 normal_3[0][j * nquad0 + i];
1636 }
1637 }
1638
1639 for (int n = 1; n < ncoords; ++n)
1640 {
1641 for (int j = 0; j < nquad2; ++j)
1642 {
1643 for (int i = 0; i < nquad0; ++i)
1644 {
1645 d0factors[1][j * nquad0 + i] +=
1646 df[3 * n][j * nquad0 * nquad1 + i] *
1647 normal_1[n][j * nquad0 + i];
1648 d1factors[1][j * nquad0 + i] +=
1649 df[3 * n + 1][j * nquad0 * nquad1 + i] *
1650 normal_1[n][j * nquad0 + i];
1651 d2factors[1][j * nquad0 + i] +=
1652 df[3 * n + 2][j * nquad0 * nquad1 + i] *
1653 normal_1[n][j * nquad0 + i];
1654
1655 d0factors[3][j * nquad0 + i] +=
1656 df[3 * n][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1657 normal_3[n][j * nquad0 + i];
1658 d1factors[3][j * nquad0 + i] +=
1659 df[3 * n + 1][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1660 normal_3[n][j * nquad0 + i];
1661 d2factors[3][j * nquad0 + i] +=
1662 df[3 * n + 2][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1663 normal_3[n][j * nquad0 + i];
1664 }
1665 }
1666 }
1667
1668 // faces 2 and 4
1669 for (int j = 0; j < nquad2; ++j)
1670 {
1671 for (int i = 0; i < nquad1; ++i)
1672 {
1673 d0factors[2][j * nquad1 + i] =
1674 df[0][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1675 normal_2[0][j * nquad1 + i];
1676 d1factors[2][j * nquad1 + i] =
1677 df[1][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1678 normal_2[0][j * nquad1 + i];
1679 d2factors[2][j * nquad1 + i] =
1680 df[2][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1681 normal_2[0][j * nquad1 + i];
1682
1683 d0factors[4][j * nquad1 + i] =
1684 df[0][j * nquad0 * nquad1 + i * nquad0] *
1685 normal_4[0][j * nquad1 + i];
1686 d1factors[4][j * nquad1 + i] =
1687 df[1][j * nquad0 * nquad1 + i * nquad0] *
1688 normal_4[0][j * nquad1 + i];
1689 d2factors[4][j * nquad1 + i] =
1690 df[2][j * nquad0 * nquad1 + i * nquad0] *
1691 normal_4[0][j * nquad1 + i];
1692 }
1693 }
1694
1695 for (int n = 1; n < ncoords; ++n)
1696 {
1697 for (int j = 0; j < nquad2; ++j)
1698 {
1699 for (int i = 0; i < nquad1; ++i)
1700 {
1701 d0factors[2][j * nquad1 + i] +=
1702 df[3 * n][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1703 normal_2[n][j * nquad1 + i];
1704 d1factors[2][j * nquad1 + i] +=
1705 df[3 * n + 1]
1706 [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1707 normal_2[n][j * nquad1 + i];
1708 d2factors[2][j * nquad1 + i] +=
1709 df[3 * n + 2]
1710 [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1711 normal_2[n][j * nquad1 + i];
1712
1713 d0factors[4][j * nquad1 + i] +=
1714 df[3 * n][j * nquad0 * nquad1 + i * nquad0] *
1715 normal_4[n][j * nquad1 + i];
1716 d1factors[4][j * nquad1 + i] +=
1717 df[3 * n + 1][j * nquad0 * nquad1 + i * nquad0] *
1718 normal_4[n][j * nquad1 + i];
1719 d2factors[4][j * nquad1 + i] +=
1720 df[3 * n + 2][j * nquad0 * nquad1 + i * nquad0] *
1721 normal_4[n][j * nquad1 + i];
1722 }
1723 }
1724 }
1725 }
1726 else
1727 {
1728 // Face 0
1729 for (int i = 0; i < nquad0 * nquad1; ++i)
1730 {
1731 d0factors[0][i] = df[0][0] * normal_0[0][i];
1732 d1factors[0][i] = df[1][0] * normal_0[0][i];
1733 d2factors[0][i] = df[2][0] * normal_0[0][i];
1734 }
1735
1736 for (int n = 1; n < ncoords; ++n)
1737 {
1738 for (int i = 0; i < nquad0 * nquad1; ++i)
1739 {
1740 d0factors[0][i] += df[3 * n][0] * normal_0[n][i];
1741 d1factors[0][i] += df[3 * n + 1][0] * normal_0[n][i];
1742 d2factors[0][i] += df[3 * n + 2][0] * normal_0[n][i];
1743 }
1744 }
1745
1746 // faces 1 and 3
1747 for (int i = 0; i < nquad0 * nquad2; ++i)
1748 {
1749 d0factors[1][i] = df[0][0] * normal_1[0][i];
1750 d0factors[3][i] = df[0][0] * normal_3[0][i];
1751
1752 d1factors[1][i] = df[1][0] * normal_1[0][i];
1753 d1factors[3][i] = df[1][0] * normal_3[0][i];
1754
1755 d2factors[1][i] = df[2][0] * normal_1[0][i];
1756 d2factors[3][i] = df[2][0] * normal_3[0][i];
1757 }
1758
1759 for (int n = 1; n < ncoords; ++n)
1760 {
1761 for (int i = 0; i < nquad0 * nquad2; ++i)
1762 {
1763 d0factors[1][i] += df[3 * n][0] * normal_1[n][i];
1764 d0factors[3][i] += df[3 * n][0] * normal_3[n][i];
1765
1766 d1factors[1][i] += df[3 * n + 1][0] * normal_1[n][i];
1767 d1factors[3][i] += df[3 * n + 1][0] * normal_3[n][i];
1768
1769 d2factors[1][i] += df[3 * n + 2][0] * normal_1[n][i];
1770 d2factors[3][i] += df[3 * n + 2][0] * normal_3[n][i];
1771 }
1772 }
1773
1774 // faces 2 and 4
1775 for (int i = 0; i < nquad1 * nquad2; ++i)
1776 {
1777 d0factors[2][i] = df[0][0] * normal_2[0][i];
1778 d0factors[4][i] = df[0][0] * normal_4[0][i];
1779
1780 d1factors[2][i] = df[1][0] * normal_2[0][i];
1781 d1factors[4][i] = df[1][0] * normal_4[0][i];
1782
1783 d2factors[2][i] = df[2][0] * normal_2[0][i];
1784 d2factors[4][i] = df[2][0] * normal_4[0][i];
1785 }
1786
1787 for (int n = 1; n < ncoords; ++n)
1788 {
1789 for (int i = 0; i < nquad1 * nquad2; ++i)
1790 {
1791 d0factors[2][i] += df[3 * n][0] * normal_2[n][i];
1792 d0factors[4][i] += df[3 * n][0] * normal_4[n][i];
1793
1794 d1factors[2][i] += df[3 * n + 1][0] * normal_2[n][i];
1795 d1factors[4][i] += df[3 * n + 1][0] * normal_4[n][i];
1796
1797 d2factors[2][i] += df[3 * n + 2][0] * normal_2[n][i];
1798 d2factors[4][i] += df[3 * n + 2][0] * normal_4[n][i];
1799 }
1800 }
1801 }
1802}
1803} // namespace Nektar::LocalRegions
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
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:282
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:292
SpatialDomains::Geometry * GetGeom() const
SpatialDomains::Geometry * m_geom
Definition Expansion.h:279
void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition Expansion.h:280
const NormalVector & GetTraceNormal(const int id)
PrismExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, SpatialDomains::Geometry3D *geom)
Constructor using BasisKey class for quadrature points and order definition.
Definition PrismExp.cpp:45
DNekMatSharedPtr v_CreateStdMatrix(const StdRegions::StdMatrixKey &mkey) override
DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const override
Definition PrismExp.cpp:440
void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition PrismExp.cpp:975
void v_DropLocStaticCondMatrix(const MatrixKey &mkey) override
void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition PrismExp.cpp:990
NekDouble v_PhysEvalFirstDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
Definition PrismExp.cpp:513
void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true) override
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 PrismExp.cpp:261
DNekScalMatSharedPtr v_GetLocMatrix(const MatrixKey &mkey) override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculates the inner product .
Definition PrismExp.cpp:327
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 PrismExp.cpp:123
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition PrismExp.h:200
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition PrismExp.cpp:267
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 PrismExp.cpp:528
void v_GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
Get the coordinates #coords at the local coordinates #Lcoords.
Definition PrismExp.cpp:464
void v_NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors, Array< OneD, Array< OneD, NekDouble > > &d2factors) override
: This method gets all of the factors which are required as part of the Gradient Jump Penalty stabili...
void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
Definition PrismExp.cpp:378
StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const override
Definition PrismExp.cpp:447
void v_ComputeTraceNormal(const int face) override
Get the normals along specficied face Get the face normals interplated to a points0 x points 0 type d...
Definition PrismExp.cpp:686
void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition PrismExp.cpp:968
void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp) override
Calculate the Laplacian multiplication in a matrix-free manner.
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition PrismExp.cpp:334
NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
Definition PrismExp.cpp:493
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 PrismExp.cpp:209
void v_GetTracePhysMap(const int face, Array< OneD, int > &outarray) override
Definition PrismExp.cpp:585
void v_DropLocMatrix(const MatrixKey &mkey) override
NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals) override
Definition PrismExp.cpp:501
void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey) override
Definition PrismExp.cpp:997
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition PrismExp.h:198
DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const MatrixKey &mkey) override
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrate the physical point list inarray over prismatic region and return the value.
Definition PrismExp.cpp:96
void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition PrismExp.cpp:481
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
3D geometry information
Definition Geometry3D.h:50
NekDouble GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.
Definition Geometry.h:558
NekDouble GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
Definition Geometry.h:548
int GetVid(int i) const
Returns global id of vertex i of this object.
Definition Geometry.h:353
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition Geometry.h:279
void FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition Geometry.h:460
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition Geometry.h:306
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.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
const LibUtilities::PointsKeyVector GetPointsKeys() const
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
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
@ eModified_B
Principle Modified Functions .
Definition BasisType.h:49
@ eModified_A
Principle Modified Functions .
Definition BasisType.h:48
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition GeomFactors.h:58
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< StdPrismExp > StdPrismExpSharedPtr
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
std::shared_ptr< DNekMat > DNekMatSharedPtr
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 > max(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:305
scalarT< T > min(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:300
scalarT< T > sqrt(scalarT< T > in)
Definition scalar.hpp:290