Nektar++
HexExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: HexExp.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: Methods for Hex expansion in local regoins
32//
33///////////////////////////////////////////////////////////////////////////////
34
37#include <LocalRegions/HexExp.h>
39
40using namespace std;
41
42namespace Nektar
43{
44namespace LocalRegions
45{
46/**
47 * @class HexExp
48 * Defines a hexahedral local expansion.
49 */
50
51/**
52 * \brief Constructor using BasisKey class for quadrature points and
53 * order definition
54 *
55 * @param Ba Basis key for first coordinate.
56 * @param Bb Basis key for second coordinate.
57 * @param Bc Basis key for third coordinate.
58 */
60 const LibUtilities::BasisKey &Bb,
61 const LibUtilities::BasisKey &Bc,
63 : StdExpansion(Ba.GetNumModes() * Bb.GetNumModes() * Bc.GetNumModes(), 3,
64 Ba, Bb, Bc),
65 StdExpansion3D(Ba.GetNumModes() * Bb.GetNumModes() * Bc.GetNumModes(), Ba,
66 Bb, Bc),
67 StdHexExp(Ba, Bb, Bc), Expansion(geom), Expansion3D(geom),
68 m_matrixManager(
69 std::bind(&Expansion3D::CreateMatrix, this, std::placeholders::_1),
70 std::string("HexExpMatrix")),
71 m_staticCondMatrixManager(std::bind(&Expansion::CreateStaticCondMatrix,
72 this, std::placeholders::_1),
73 std::string("HexExpStaticCondMatrix"))
74{
75}
76
77/**
78 * \brief Copy Constructor
79 *
80 * @param T HexExp to copy.
81 */
83 : StdExpansion(T), StdExpansion3D(T), StdHexExp(T), Expansion(T),
84 Expansion3D(T), m_matrixManager(T.m_matrixManager),
85 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
86{
87}
88
89//-----------------------------
90// Integration Methods
91//-----------------------------
92/**
93 * \brief Integrate the physical point list \a inarray over region
94 *
95 * @param inarray definition of function to be returned at
96 * quadrature points of expansion.
97 * @returns \f$\int^1_{-1}\int^1_{-1} \int^1_{-1}
98 * u(\eta_1, \eta_2, \eta_3) J[i,j,k] d \eta_1 d \eta_2 d \eta_3 \f$
99 * where \f$inarray[i,j,k] = u(\eta_{1i},\eta_{2j},\eta_{3k}) \f$
100 * and \f$ J[i,j,k] \f$ is the Jacobian evaluated at the quadrature
101 * point.
102 */
104{
105 int nquad0 = m_base[0]->GetNumPoints();
106 int nquad1 = m_base[1]->GetNumPoints();
107 int nquad2 = m_base[2]->GetNumPoints();
109 NekDouble returnVal;
110 Array<OneD, NekDouble> tmp(nquad0 * nquad1 * nquad2);
111
112 // multiply inarray with Jacobian
113
114 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
115 {
116 Vmath::Vmul(nquad0 * nquad1 * nquad2, &jac[0], 1,
117 (NekDouble *)&inarray[0], 1, &tmp[0], 1);
118 }
119 else
120 {
121 Vmath::Smul(nquad0 * nquad1 * nquad2, (NekDouble)jac[0],
122 (NekDouble *)&inarray[0], 1, &tmp[0], 1);
123 }
124
125 // call StdHexExp version;
126 returnVal = StdHexExp::v_Integral(tmp);
127
128 return returnVal;
129}
130
131//-----------------------------
132// Differentiation Methods
133//-----------------------------
134/**
135 * \brief Calculate the derivative of the physical points
136 *
137 * For Hexahedral region can use the Tensor_Deriv function defined
138 * under StdExpansion.
139 * @param inarray Input array
140 * @param out_d0 Derivative of \a inarray in first direction.
141 * @param out_d1 Derivative of \a inarray in second direction.
142 * @param out_d2 Derivative of \a inarray in third direction.
143 */
148{
149 int nquad0 = m_base[0]->GetNumPoints();
150 int nquad1 = m_base[1]->GetNumPoints();
151 int nquad2 = m_base[2]->GetNumPoints();
152 int ntot = nquad0 * nquad1 * nquad2;
153
155 m_metricinfo->GetDerivFactors(GetPointsKeys());
159
160 StdHexExp::v_PhysDeriv(inarray, Diff0, Diff1, Diff2);
161
162 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
163 {
164 if (out_d0.size())
165 {
166 Vmath::Vmul(ntot, &df[0][0], 1, &Diff0[0], 1, &out_d0[0], 1);
167 Vmath::Vvtvp(ntot, &df[1][0], 1, &Diff1[0], 1, &out_d0[0], 1,
168 &out_d0[0], 1);
169 Vmath::Vvtvp(ntot, &df[2][0], 1, &Diff2[0], 1, &out_d0[0], 1,
170 &out_d0[0], 1);
171 }
172
173 if (out_d1.size())
174 {
175 Vmath::Vmul(ntot, &df[3][0], 1, &Diff0[0], 1, &out_d1[0], 1);
176 Vmath::Vvtvp(ntot, &df[4][0], 1, &Diff1[0], 1, &out_d1[0], 1,
177 &out_d1[0], 1);
178 Vmath::Vvtvp(ntot, &df[5][0], 1, &Diff2[0], 1, &out_d1[0], 1,
179 &out_d1[0], 1);
180 }
181
182 if (out_d2.size())
183 {
184 Vmath::Vmul(ntot, &df[6][0], 1, &Diff0[0], 1, &out_d2[0], 1);
185 Vmath::Vvtvp(ntot, &df[7][0], 1, &Diff1[0], 1, &out_d2[0], 1,
186 &out_d2[0], 1);
187 Vmath::Vvtvp(ntot, &df[8][0], 1, &Diff2[0], 1, &out_d2[0], 1,
188 &out_d2[0], 1);
189 }
190 }
191 else // regular geometry
192 {
193 if (out_d0.size())
194 {
195 Vmath::Smul(ntot, df[0][0], &Diff0[0], 1, &out_d0[0], 1);
196 Blas::Daxpy(ntot, df[1][0], &Diff1[0], 1, &out_d0[0], 1);
197 Blas::Daxpy(ntot, df[2][0], &Diff2[0], 1, &out_d0[0], 1);
198 }
199
200 if (out_d1.size())
201 {
202 Vmath::Smul(ntot, df[3][0], &Diff0[0], 1, &out_d1[0], 1);
203 Blas::Daxpy(ntot, df[4][0], &Diff1[0], 1, &out_d1[0], 1);
204 Blas::Daxpy(ntot, df[5][0], &Diff2[0], 1, &out_d1[0], 1);
205 }
206
207 if (out_d2.size())
208 {
209 Vmath::Smul(ntot, df[6][0], &Diff0[0], 1, &out_d2[0], 1);
210 Blas::Daxpy(ntot, df[7][0], &Diff1[0], 1, &out_d2[0], 1);
211 Blas::Daxpy(ntot, df[8][0], &Diff2[0], 1, &out_d2[0], 1);
212 }
213 }
214}
215
216/**
217 * \brief Calculate the derivative of the physical points in a single
218 * direction.
219 *
220 * @param dir Direction in which to compute derivative.
221 * Valid values are 0, 1, 2.
222 * @param inarray Input array.
223 * @param outarray Output array.
224 */
225void HexExp::v_PhysDeriv(const int dir,
226 const Array<OneD, const NekDouble> &inarray,
227 Array<OneD, NekDouble> &outarray)
228{
229 switch (dir)
230 {
231 case 0:
232 {
233 PhysDeriv(inarray, outarray, NullNekDouble1DArray,
235 }
236 break;
237 case 1:
238 {
239 PhysDeriv(inarray, NullNekDouble1DArray, outarray,
241 }
242 break;
243 case 2:
244 {
246 outarray);
247 }
248 break;
249 default:
250 {
251 ASSERTL1(false, "input dir is out of range");
252 }
253 break;
254 }
255}
256
258 const Array<OneD, const NekDouble> &inarray,
259 const Array<OneD, const NekDouble> &direction,
260 Array<OneD, NekDouble> &outarray)
261{
262
263 int shapedim = 3;
264 int nquad0 = m_base[0]->GetNumPoints();
265 int nquad1 = m_base[1]->GetNumPoints();
266 int nquad2 = m_base[2]->GetNumPoints();
267 int ntot = nquad0 * nquad1 * nquad2;
268
270 m_metricinfo->GetDerivFactors(GetPointsKeys());
274
275 StdHexExp::v_PhysDeriv(inarray, Diff0, Diff1, Diff2);
276
277 Array<OneD, Array<OneD, NekDouble>> dfdir(shapedim);
278 Expansion::ComputeGmatcdotMF(df, direction, dfdir);
279
280 Vmath::Vmul(ntot, &dfdir[0][0], 1, &Diff0[0], 1, &outarray[0], 1);
281 Vmath::Vvtvp(ntot, &dfdir[1][0], 1, &Diff1[0], 1, &outarray[0], 1,
282 &outarray[0], 1);
283 Vmath::Vvtvp(ntot, &dfdir[2][0], 1, &Diff2[0], 1, &outarray[0], 1,
284 &outarray[0], 1);
285}
286
287//-----------------------------
288// Transforms
289//-----------------------------
290
291/**
292 * \brief Forward transform from physical quadrature space stored in \a
293 * inarray and evaluate the expansion coefficients and store in
294 * \a (this)->_coeffs
295 *
296 * @param inarray Input array
297 * @param outarray Output array
298 */
300 Array<OneD, NekDouble> &outarray)
301{
302 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
303 m_base[2]->Collocation())
304 {
305 Vmath::Vcopy(GetNcoeffs(), &inarray[0], 1, &outarray[0], 1);
306 }
307 else
308 {
309 IProductWRTBase(inarray, outarray);
310
311 // get Mass matrix inverse
313 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
314
315 // copy inarray in case inarray == outarray
316 DNekVec in(m_ncoeffs, outarray);
317 DNekVec out(m_ncoeffs, outarray, eWrapper);
318
319 out = (*matsys) * in;
320 }
321}
322
323//-----------------------------
324// Inner product functions
325//-----------------------------
326
327/**
328 * \brief Calculate the inner product of inarray with respect to the
329 * elements basis.
330 *
331 * @param inarray Input array of physical space data.
332 * @param outarray Output array of data.
333 */
335 Array<OneD, NekDouble> &outarray)
336{
337 HexExp::v_IProductWRTBase_SumFac(inarray, outarray);
338}
339
340/**
341 * \brief Calculate the inner product of inarray with respect to the
342 * given basis B = base0 * base1 * base2.
343 *
344 * \f$ \begin{array}{rcl} I_{pqr} = (\phi_{pqr}, u)_{\delta}
345 * & = & \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \sum_{k=0}^{nq_2}
346 * \psi_{p}^{a} (\xi_{1i}) \psi_{q}^{a} (\xi_{2j}) \psi_{r}^{a}
347 * (\xi_{3k}) w_i w_j w_k u(\xi_{1,i} \xi_{2,j} \xi_{3,k})
348 * J_{i,j,k}\\ & = & \sum_{i=0}^{nq_0} \psi_p^a(\xi_{1,i})
349 * \sum_{j=0}^{nq_1} \psi_{q}^a(\xi_{2,j}) \sum_{k=0}^{nq_2}
350 * \psi_{r}^a u(\xi_{1i},\xi_{2j},\xi_{3k})
351 * J_{i,j,k} \end{array} \f$ \n
352 * where
353 * \f$ \phi_{pqr} (\xi_1 , \xi_2 , \xi_3)
354 * = \psi_p^a ( \xi_1) \psi_{q}^a (\xi_2) \psi_{r}^a (\xi_3) \f$ \n
355 * which can be implemented as \n
356 * \f$f_{r} (\xi_{3k})
357 * = \sum_{k=0}^{nq_3} \psi_{r}^a u(\xi_{1i},\xi_{2j},\xi_{3k})
358 * J_{i,j,k} = {\bf B_3 U} \f$ \n
359 * \f$ g_{q} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{q}^a (\xi_{2j})
360 * f_{r} (\xi_{3k}) = {\bf B_2 F} \f$ \n
361 * \f$ (\phi_{pqr}, u)_{\delta}
362 * = \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{q} (\xi_{3k})
363 * = {\bf B_1 G} \f$
364 *
365 * @param base0 Basis to integrate wrt in first dimension.
366 * @param base1 Basis to integrate wrt in second dimension.
367 * @param base2 Basis to integrate wrt in third dimension.
368 * @param inarray Input array.
369 * @param outarray Output array.
370 * @param coll_check (not used)
371 */
373 const Array<OneD, const NekDouble> &inarray,
374 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
375{
376 int nquad0 = m_base[0]->GetNumPoints();
377 int nquad1 = m_base[1]->GetNumPoints();
378 int nquad2 = m_base[2]->GetNumPoints();
379 int order0 = m_base[0]->GetNumModes();
380 int order1 = m_base[1]->GetNumModes();
381
382 Array<OneD, NekDouble> wsp(nquad0 * nquad1 * (nquad2 + order0) +
383 order0 * order1 * nquad2);
384
385 if (multiplybyweights)
386 {
387 Array<OneD, NekDouble> tmp(inarray.size());
388
389 MultiplyByQuadratureMetric(inarray, tmp);
391 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
392 tmp, outarray, wsp, true, true, true);
393 }
394 else
395 {
397 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
398 inarray, outarray, wsp, true, true, true);
399 }
400}
401
403 const Array<OneD, const NekDouble> &inarray,
404 Array<OneD, NekDouble> &outarray)
405{
406 HexExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
407}
408
409/**
410 * @brief Calculates the inner product \f$ I_{pqr} = (u,
411 * \partial_{x_i} \phi_{pqr}) \f$.
412 *
413 * The derivative of the basis functions is performed using the chain
414 * rule in order to incorporate the geometric factors. Assuming that
415 * the basis functions are a tensor product
416 * \f$\phi_{pqr}(\xi_1,\xi_2,\xi_3) =
417 * \phi_1(\xi_1)\phi_2(\xi_2)\phi_3(\xi_3)\f$, in the hexahedral
418 * element, this is straightforward and yields the result
419 *
420 * \f[
421 * I_{pqr} = \sum_{k=1}^3 \left(u, \frac{\partial u}{\partial \xi_k}
422 * \frac{\partial \xi_k}{\partial x_i}\right)
423 * \f]
424 *
425 * @param dir Direction in which to take the derivative.
426 * @param inarray The function \f$ u \f$.
427 * @param outarray Value of the inner product.
428 */
430 const int dir, const Array<OneD, const NekDouble> &inarray,
431 Array<OneD, NekDouble> &outarray)
432{
433 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
434
435 const int nq0 = m_base[0]->GetNumPoints();
436 const int nq1 = m_base[1]->GetNumPoints();
437 const int nq2 = m_base[2]->GetNumPoints();
438 const int nq = nq0 * nq1 * nq2;
439 const int nm0 = m_base[0]->GetNumModes();
440 const int nm1 = m_base[1]->GetNumModes();
441
442 Array<OneD, NekDouble> alloc(4 * nq + m_ncoeffs + nm0 * nq2 * (nq1 + nm1));
443 Array<OneD, NekDouble> tmp1(alloc); // Quad metric
444 Array<OneD, NekDouble> tmp2(alloc + nq); // Dir1 metric
445 Array<OneD, NekDouble> tmp3(alloc + 2 * nq); // Dir2 metric
446 Array<OneD, NekDouble> tmp4(alloc + 3 * nq); // Dir3 metric
447 Array<OneD, NekDouble> tmp5(alloc + 4 * nq); // iprod tmp
448 Array<OneD, NekDouble> wsp(tmp5 + m_ncoeffs); // Wsp
449
450 MultiplyByQuadratureMetric(inarray, tmp1);
451
453 tmp2D[0] = tmp2;
454 tmp2D[1] = tmp3;
455 tmp2D[2] = tmp4;
456
457 HexExp::v_AlignVectorToCollapsedDir(dir, tmp1, tmp2D);
458
459 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
460 m_base[2]->GetBdata(), tmp2, outarray, wsp,
461 false, true, true);
462
463 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
464 m_base[2]->GetBdata(), tmp3, tmp5, wsp, true,
465 false, true);
466 Vmath::Vadd(m_ncoeffs, tmp5, 1, outarray, 1, outarray, 1);
467
468 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
469 m_base[2]->GetDbdata(), tmp4, tmp5, wsp, true,
470 true, false);
471 Vmath::Vadd(m_ncoeffs, tmp5, 1, outarray, 1, outarray, 1);
472}
473
475 const int dir, const Array<OneD, const NekDouble> &inarray,
477{
478 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
479
480 const int nq0 = m_base[0]->GetNumPoints();
481 const int nq1 = m_base[1]->GetNumPoints();
482 const int nq2 = m_base[2]->GetNumPoints();
483 const int nq = nq0 * nq1 * nq2;
484
486 m_metricinfo->GetDerivFactors(GetPointsKeys());
487
488 Array<OneD, NekDouble> tmp1(nq); // Quad metric
489
490 Array<OneD, NekDouble> tmp2 = outarray[0]; // Dir1 metric
491 Array<OneD, NekDouble> tmp3 = outarray[1]; // Dir2 metric
492 Array<OneD, NekDouble> tmp4 = outarray[2];
493
494 Vmath::Vcopy(nq, inarray, 1, tmp1, 1); // Dir3 metric
495
496 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
497 {
498 Vmath::Vmul(nq, &df[3 * dir][0], 1, tmp1.get(), 1, tmp2.get(), 1);
499 Vmath::Vmul(nq, &df[3 * dir + 1][0], 1, tmp1.get(), 1, tmp3.get(), 1);
500 Vmath::Vmul(nq, &df[3 * dir + 2][0], 1, tmp1.get(), 1, tmp4.get(), 1);
501 }
502 else
503 {
504 Vmath::Smul(nq, df[3 * dir][0], tmp1.get(), 1, tmp2.get(), 1);
505 Vmath::Smul(nq, df[3 * dir + 1][0], tmp1.get(), 1, tmp3.get(), 1);
506 Vmath::Smul(nq, df[3 * dir + 2][0], tmp1.get(), 1, tmp4.get(), 1);
507 }
508}
509
510/**
511 *
512 * @param dir Vector direction in which to take the derivative.
513 * @param inarray The function \f$ u \f$.
514 * @param outarray Value of the inner product.
515 */
517 const Array<OneD, const NekDouble> &direction,
518 const Array<OneD, const NekDouble> &inarray,
519 Array<OneD, NekDouble> &outarray)
520{
521 int shapedim = 3;
522 const int nq0 = m_base[0]->GetNumPoints();
523 const int nq1 = m_base[1]->GetNumPoints();
524 const int nq2 = m_base[2]->GetNumPoints();
525 const int nq = nq0 * nq1 * nq2;
526 const int nm0 = m_base[0]->GetNumModes();
527 const int nm1 = m_base[1]->GetNumModes();
528
530 m_metricinfo->GetDerivFactors(GetPointsKeys());
531
532 Array<OneD, NekDouble> alloc(4 * nq + m_ncoeffs + nm0 * nq2 * (nq1 + nm1));
533 Array<OneD, NekDouble> tmp1(alloc); // Quad metric
534 Array<OneD, NekDouble> tmp2(alloc + nq); // Dir1 metric
535 Array<OneD, NekDouble> tmp3(alloc + 2 * nq); // Dir2 metric
536 Array<OneD, NekDouble> tmp4(alloc + 3 * nq); // Dir3 metric
537 Array<OneD, NekDouble> tmp5(alloc + 4 * nq); // iprod tmp
538 Array<OneD, NekDouble> wsp(tmp5 + m_ncoeffs); // Wsp
539
540 MultiplyByQuadratureMetric(inarray, tmp1);
541
542 Array<OneD, Array<OneD, NekDouble>> dfdir(shapedim);
543 Expansion::ComputeGmatcdotMF(df, direction, dfdir);
544
545 Vmath::Vmul(nq, &dfdir[0][0], 1, tmp1.get(), 1, tmp2.get(), 1);
546 Vmath::Vmul(nq, &dfdir[1][0], 1, tmp1.get(), 1, tmp3.get(), 1);
547 Vmath::Vmul(nq, &dfdir[2][0], 1, tmp1.get(), 1, tmp4.get(), 1);
548
549 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
550 m_base[2]->GetBdata(), tmp2, outarray, wsp,
551 false, true, true);
552
553 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
554 m_base[2]->GetBdata(), tmp3, tmp5, wsp, true,
555 false, true);
556
557 Vmath::Vadd(m_ncoeffs, tmp5, 1, outarray, 1, outarray, 1);
558
559 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
560 m_base[2]->GetDbdata(), tmp4, tmp5, wsp, true,
561 true, false);
562
563 Vmath::Vadd(m_ncoeffs, tmp5, 1, outarray, 1, outarray, 1);
564}
565
566//-----------------------------
567// Evaluation functions
568//-----------------------------
569
570/**
571 * Given the local cartesian coordinate \a Lcoord evaluate the
572 * value of physvals at this point by calling through to the
573 * StdExpansion method
574 */
576 const Array<OneD, const NekDouble> &Lcoord,
577 const Array<OneD, const NekDouble> &physvals)
578{
579 // Evaluate point in local coordinates.
580 return StdExpansion3D::v_PhysEvaluate(Lcoord, physvals);
581}
582
584 const Array<OneD, const NekDouble> &physvals)
585{
587
588 ASSERTL0(m_geom, "m_geom not defined");
589 m_geom->GetLocCoords(coord, Lcoord);
590 return StdExpansion3D::v_PhysEvaluate(Lcoord, physvals);
591}
592
594 const Array<OneD, const NekDouble> &inarray,
595 std::array<NekDouble, 3> &firstOrderDerivs)
596{
597 Array<OneD, NekDouble> Lcoord(3);
598 ASSERTL0(m_geom, "m_geom not defined");
599 m_geom->GetLocCoords(coord, Lcoord);
600 return StdHexExp::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
601}
602
604{
606 m_base[0]->GetBasisKey(), m_base[1]->GetBasisKey(),
607 m_base[2]->GetBasisKey());
608}
609
611{
613 m_base[0]->GetPointsKey());
615 m_base[1]->GetPointsKey());
617 m_base[2]->GetPointsKey());
618
620 bkey2);
621}
622
623/**
624 * \brief Retrieves the physical coordinates of a given set of
625 * reference coordinates.
626 *
627 * @param Lcoords Local coordinates in reference space.
628 * @param coords Corresponding coordinates in physical space.
629 */
632{
633 int i;
634
635 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[0] <= 1.0 && Lcoords[1] >= -1.0 &&
636 Lcoords[1] <= 1.0 && Lcoords[2] >= -1.0 && Lcoords[2] <= 1.0,
637 "Local coordinates are not in region [-1,1]");
638
639 m_geom->FillGeom();
640
641 for (i = 0; i < m_geom->GetCoordim(); ++i)
642 {
643 coords[i] = m_geom->GetCoord(i, Lcoords);
644 }
645}
646
648 Array<OneD, NekDouble> &coords_1,
649 Array<OneD, NekDouble> &coords_2)
650{
651 Expansion::v_GetCoords(coords_0, coords_1, coords_2);
652}
653
654//-----------------------------
655// Helper functions
656//-----------------------------
657
658/// Return the region shape using the enum-list of ShapeType
660{
662}
663
665 const NekDouble *data, const std::vector<unsigned int> &nummodes,
666 const int mode_offset, NekDouble *coeffs,
667 std::vector<LibUtilities::BasisType> &fromType)
668{
669 int data_order0 = nummodes[mode_offset];
670 int fillorder0 = min(m_base[0]->GetNumModes(), data_order0);
671 int data_order1 = nummodes[mode_offset + 1];
672 int order1 = m_base[1]->GetNumModes();
673 int fillorder1 = min(order1, data_order1);
674 int data_order2 = nummodes[mode_offset + 2];
675 int order2 = m_base[2]->GetNumModes();
676 int fillorder2 = min(order2, data_order2);
677
678 // Check if same basis
679 if (fromType[0] != m_base[0]->GetBasisType() ||
680 fromType[1] != m_base[1]->GetBasisType() ||
681 fromType[2] != m_base[2]->GetBasisType())
682 {
683 // Construct a hex with the appropriate basis type at our
684 // quadrature points, and one more to do a forwards
685 // transform. We can then copy the output to coeffs.
687 LibUtilities::BasisKey(fromType[0], data_order0,
688 m_base[0]->GetPointsKey()),
689 LibUtilities::BasisKey(fromType[1], data_order1,
690 m_base[1]->GetPointsKey()),
691 LibUtilities::BasisKey(fromType[2], data_order2,
692 m_base[2]->GetPointsKey()));
693 StdRegions::StdHexExp tmpHex2(m_base[0]->GetBasisKey(),
694 m_base[1]->GetBasisKey(),
695 m_base[2]->GetBasisKey());
696
697 Array<OneD, const NekDouble> tmpData(tmpHex.GetNcoeffs(), data);
698 Array<OneD, NekDouble> tmpBwd(tmpHex2.GetTotPoints());
699 Array<OneD, NekDouble> tmpOut(tmpHex2.GetNcoeffs());
700
701 tmpHex.BwdTrans(tmpData, tmpBwd);
702 tmpHex2.FwdTrans(tmpBwd, tmpOut);
703 Vmath::Vcopy(tmpOut.size(), &tmpOut[0], 1, coeffs, 1);
704
705 return;
706 }
707
708 switch (m_base[0]->GetBasisType())
709 {
711 {
712 int i, j;
713 int cnt = 0;
714 int cnt1 = 0;
715
717 "Extraction routine not set up for this basis");
719 "Extraction routine not set up for this basis");
720
721 Vmath::Zero(m_ncoeffs, coeffs, 1);
722 for (j = 0; j < fillorder0; ++j)
723 {
724 for (i = 0; i < fillorder1; ++i)
725 {
726 Vmath::Vcopy(fillorder2, &data[cnt], 1, &coeffs[cnt1], 1);
727 cnt += data_order2;
728 cnt1 += order2;
729 }
730
731 // count out data for j iteration
732 for (i = fillorder1; i < data_order1; ++i)
733 {
734 cnt += data_order2;
735 }
736
737 for (i = fillorder1; i < order1; ++i)
738 {
739 cnt1 += order2;
740 }
741 }
742 break;
743 }
745 {
746 LibUtilities::PointsKey p0(nummodes[0],
748 LibUtilities::PointsKey p1(nummodes[1],
750 LibUtilities::PointsKey p2(nummodes[2],
752 LibUtilities::PointsKey t0(m_base[0]->GetNumModes(),
754 LibUtilities::PointsKey t1(m_base[1]->GetNumModes(),
756 LibUtilities::PointsKey t2(m_base[2]->GetNumModes(),
758 LibUtilities::Interp3D(p0, p1, p2, data, t0, t1, t2, coeffs);
759 }
760 break;
761 default:
762 ASSERTL0(false, "basis is either not set up or not "
763 "hierarchicial");
764 }
765}
766
767void HexExp::v_GetTracePhysMap(const int face, Array<OneD, int> &outarray)
768{
769 int nquad0 = m_base[0]->GetNumPoints();
770 int nquad1 = m_base[1]->GetNumPoints();
771 int nquad2 = m_base[2]->GetNumPoints();
772
773 int nq0 = 0;
774 int nq1 = 0;
775
776 switch (face)
777 {
778 case 0:
779 nq0 = nquad0;
780 nq1 = nquad1;
781
782 // Directions A and B positive
783 if (outarray.size() != nq0 * nq1)
784 {
785 outarray = Array<OneD, int>(nq0 * nq1);
786 }
787
788 for (int i = 0; i < nquad0 * nquad1; ++i)
789 {
790 outarray[i] = i;
791 }
792
793 break;
794 case 1:
795 nq0 = nquad0;
796 nq1 = nquad2;
797 // Direction A and B positive
798 if (outarray.size() != nq0 * nq1)
799 {
800 outarray = Array<OneD, int>(nq0 * nq1);
801 }
802
803 // Direction A and B positive
804 for (int k = 0; k < nquad2; k++)
805 {
806 for (int i = 0; i < nquad0; ++i)
807 {
808 outarray[k * nquad0 + i] = nquad0 * nquad1 * k + i;
809 }
810 }
811 break;
812 case 2:
813 nq0 = nquad1;
814 nq1 = nquad2;
815
816 // Direction A and B positive
817 if (outarray.size() != nq0 * nq1)
818 {
819 outarray = Array<OneD, int>(nq0 * nq1);
820 }
821
822 for (int i = 0; i < nquad1 * nquad2; i++)
823 {
824 outarray[i] = nquad0 - 1 + i * nquad0;
825 }
826 break;
827 case 3:
828 nq0 = nquad0;
829 nq1 = nquad2;
830
831 // Direction A and B positive
832 if (outarray.size() != nq0 * nq1)
833 {
834 outarray = Array<OneD, int>(nq0 * nq1);
835 }
836
837 for (int k = 0; k < nquad2; k++)
838 {
839 for (int i = 0; i < nquad0; i++)
840 {
841 outarray[k * nquad0 + i] =
842 (nquad0 * (nquad1 - 1)) + (k * nquad0 * nquad1) + i;
843 }
844 }
845 break;
846 case 4:
847 nq0 = nquad1;
848 nq1 = nquad2;
849
850 // Direction A and B positive
851 if (outarray.size() != nq0 * nq1)
852 {
853 outarray = Array<OneD, int>(nq0 * nq1);
854 }
855
856 for (int i = 0; i < nquad1 * nquad2; i++)
857 {
858 outarray[i] = i * nquad0;
859 }
860 break;
861 case 5:
862 nq0 = nquad0;
863 nq1 = nquad1;
864 // Directions A and B positive
865 if (outarray.size() != nq0 * nq1)
866 {
867 outarray = Array<OneD, int>(nq0 * nq1);
868 }
869
870 for (int i = 0; i < nquad0 * nquad1; i++)
871 {
872 outarray[i] = nquad0 * nquad1 * (nquad2 - 1) + i;
873 }
874
875 break;
876 default:
877 ASSERTL0(false, "face value (> 5) is out of range");
878 break;
879 }
880}
881
883{
884 int i;
885 const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
886 GetGeom()->GetMetricInfo();
887 SpatialDomains::GeomType type = geomFactors->GetGtype();
888
890 for (i = 0; i < ptsKeys.size(); ++i)
891 {
892 // Need at least 2 points for computing normals
893 if (ptsKeys[i].GetNumPoints() == 1)
894 {
895 LibUtilities::PointsKey pKey(2, ptsKeys[i].GetPointsType());
896 ptsKeys[i] = pKey;
897 }
898 }
899
901 geomFactors->GetDerivFactors(ptsKeys);
902 const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
903
904 LibUtilities::BasisKey tobasis0 = GetTraceBasisKey(face, 0);
905 LibUtilities::BasisKey tobasis1 = GetTraceBasisKey(face, 1);
906
907 // Number of quadrature points in face expansion.
908 int nq_face = tobasis0.GetNumPoints() * tobasis1.GetNumPoints();
909
910 int vCoordDim = GetCoordim();
911
914 for (i = 0; i < vCoordDim; ++i)
915 {
916 normal[i] = Array<OneD, NekDouble>(nq_face);
917 }
918
919 size_t nqb = nq_face;
920 size_t nbnd = face;
923
924 // Regular geometry case
925 if ((type == SpatialDomains::eRegular) ||
927 {
928 NekDouble fac;
929 // Set up normals
930 switch (face)
931 {
932 case 0:
933 for (i = 0; i < vCoordDim; ++i)
934 {
935 normal[i][0] = -df[3 * i + 2][0];
936 }
937 break;
938 case 1:
939 for (i = 0; i < vCoordDim; ++i)
940 {
941 normal[i][0] = -df[3 * i + 1][0];
942 }
943 break;
944 case 2:
945 for (i = 0; i < vCoordDim; ++i)
946 {
947 normal[i][0] = df[3 * i][0];
948 }
949 break;
950 case 3:
951 for (i = 0; i < vCoordDim; ++i)
952 {
953 normal[i][0] = df[3 * i + 1][0];
954 }
955 break;
956 case 4:
957 for (i = 0; i < vCoordDim; ++i)
958 {
959 normal[i][0] = -df[3 * i][0];
960 }
961 break;
962 case 5:
963 for (i = 0; i < vCoordDim; ++i)
964 {
965 normal[i][0] = df[3 * i + 2][0];
966 }
967 break;
968 default:
969 ASSERTL0(false, "face is out of range (edge < 5)");
970 }
971
972 // normalise
973 fac = 0.0;
974 for (i = 0; i < vCoordDim; ++i)
975 {
976 fac += normal[i][0] * normal[i][0];
977 }
978 fac = 1.0 / sqrt(fac);
979
980 Vmath::Fill(nqb, fac, length, 1);
981 for (i = 0; i < vCoordDim; ++i)
982 {
983 Vmath::Fill(nq_face, fac * normal[i][0], normal[i], 1);
984 }
985 }
986 else // Set up deformed normals
987 {
988 int j, k;
989
990 int nqe0 = ptsKeys[0].GetNumPoints();
991 int nqe1 = ptsKeys[1].GetNumPoints();
992 int nqe2 = ptsKeys[2].GetNumPoints();
993 int nqe01 = nqe0 * nqe1;
994 int nqe02 = nqe0 * nqe2;
995 int nqe12 = nqe1 * nqe2;
996
997 int nqe;
998 if (face == 0 || face == 5)
999 {
1000 nqe = nqe01;
1001 }
1002 else if (face == 1 || face == 3)
1003 {
1004 nqe = nqe02;
1005 }
1006 else
1007 {
1008 nqe = nqe12;
1009 }
1010
1013
1014 Array<OneD, NekDouble> faceJac(nqe);
1015 Array<OneD, NekDouble> normals(vCoordDim * nqe, 0.0);
1016
1017 // Extract Jacobian along face and recover local
1018 // derivates (dx/dr) for polynomial interpolation by
1019 // multiplying m_gmat by jacobian
1020 switch (face)
1021 {
1022 case 0:
1023 for (j = 0; j < nqe; ++j)
1024 {
1025 normals[j] = -df[2][j] * jac[j];
1026 normals[nqe + j] = -df[5][j] * jac[j];
1027 normals[2 * nqe + j] = -df[8][j] * jac[j];
1028 faceJac[j] = jac[j];
1029 }
1030
1031 points0 = ptsKeys[0];
1032 points1 = ptsKeys[1];
1033 break;
1034 case 1:
1035 for (j = 0; j < nqe0; ++j)
1036 {
1037 for (k = 0; k < nqe2; ++k)
1038 {
1039 int idx = j + nqe01 * k;
1040 normals[j + k * nqe0] = -df[1][idx] * jac[idx];
1041 normals[nqe + j + k * nqe0] = -df[4][idx] * jac[idx];
1042 normals[2 * nqe + j + k * nqe0] =
1043 -df[7][idx] * jac[idx];
1044 faceJac[j + k * nqe0] = jac[idx];
1045 }
1046 }
1047 points0 = ptsKeys[0];
1048 points1 = ptsKeys[2];
1049 break;
1050 case 2:
1051 for (j = 0; j < nqe1; ++j)
1052 {
1053 for (k = 0; k < nqe2; ++k)
1054 {
1055 int idx = nqe0 - 1 + nqe0 * j + nqe01 * k;
1056 normals[j + k * nqe1] = df[0][idx] * jac[idx];
1057 normals[nqe + j + k * nqe1] = df[3][idx] * jac[idx];
1058 normals[2 * nqe + j + k * nqe1] = df[6][idx] * jac[idx];
1059 faceJac[j + k * nqe1] = jac[idx];
1060 }
1061 }
1062 points0 = ptsKeys[1];
1063 points1 = ptsKeys[2];
1064 break;
1065 case 3:
1066 for (j = 0; j < nqe0; ++j)
1067 {
1068 for (k = 0; k < nqe2; ++k)
1069 {
1070 int idx = nqe0 * (nqe1 - 1) + j + nqe01 * k;
1071 normals[j + k * nqe0] = df[1][idx] * jac[idx];
1072 normals[nqe + j + k * nqe0] = df[4][idx] * jac[idx];
1073 normals[2 * nqe + j + k * nqe0] = df[7][idx] * jac[idx];
1074 faceJac[j + k * nqe0] = jac[idx];
1075 }
1076 }
1077 points0 = ptsKeys[0];
1078 points1 = ptsKeys[2];
1079 break;
1080 case 4:
1081 for (j = 0; j < nqe1; ++j)
1082 {
1083 for (k = 0; k < nqe2; ++k)
1084 {
1085 int idx = j * nqe0 + nqe01 * k;
1086 normals[j + k * nqe1] = -df[0][idx] * jac[idx];
1087 normals[nqe + j + k * nqe1] = -df[3][idx] * jac[idx];
1088 normals[2 * nqe + j + k * nqe1] =
1089 -df[6][idx] * jac[idx];
1090 faceJac[j + k * nqe1] = jac[idx];
1091 }
1092 }
1093 points0 = ptsKeys[1];
1094 points1 = ptsKeys[2];
1095 break;
1096 case 5:
1097 for (j = 0; j < nqe01; ++j)
1098 {
1099 int idx = j + nqe01 * (nqe2 - 1);
1100 normals[j] = df[2][idx] * jac[idx];
1101 normals[nqe + j] = df[5][idx] * jac[idx];
1102 normals[2 * nqe + j] = df[8][idx] * jac[idx];
1103 faceJac[j] = jac[idx];
1104 }
1105 points0 = ptsKeys[0];
1106 points1 = ptsKeys[1];
1107 break;
1108 default:
1109 ASSERTL0(false, "face is out of range (face < 5)");
1110 }
1111
1112 Array<OneD, NekDouble> work(nq_face, 0.0);
1113 // Interpolate Jacobian and invert
1114 LibUtilities::Interp2D(points0, points1, faceJac,
1115 tobasis0.GetPointsKey(), tobasis1.GetPointsKey(),
1116 work);
1117
1118 Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
1119
1120 // interpolate
1121 for (i = 0; i < GetCoordim(); ++i)
1122 {
1123 LibUtilities::Interp2D(points0, points1, &normals[i * nqe],
1124 tobasis0.GetPointsKey(),
1125 tobasis1.GetPointsKey(), &normal[i][0]);
1126 Vmath::Vmul(nq_face, work, 1, normal[i], 1, normal[i], 1);
1127 }
1128
1129 // normalise normal vectors
1130 Vmath::Zero(nq_face, work, 1);
1131 for (i = 0; i < GetCoordim(); ++i)
1132 {
1133 Vmath::Vvtvp(nq_face, normal[i], 1, normal[i], 1, work, 1, work, 1);
1134 }
1135
1136 Vmath::Vsqrt(nq_face, work, 1, work, 1);
1137 Vmath::Sdiv(nq_face, 1.0, work, 1, work, 1);
1138
1139 Vmath::Vcopy(nqb, work, 1, length, 1);
1140
1141 for (i = 0; i < GetCoordim(); ++i)
1142 {
1143 Vmath::Vmul(nq_face, normal[i], 1, work, 1, normal[i], 1);
1144 }
1145 }
1146}
1147
1148//-----------------------------
1149// Operator creation functions
1150//-----------------------------
1152 Array<OneD, NekDouble> &outarray,
1153 const StdRegions::StdMatrixKey &mkey)
1154{
1155 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1156}
1157
1159 Array<OneD, NekDouble> &outarray,
1160 const StdRegions::StdMatrixKey &mkey)
1161{
1162 HexExp::v_LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1163}
1164
1165void HexExp::v_LaplacianMatrixOp(const int k1, const int k2,
1166 const Array<OneD, const NekDouble> &inarray,
1167 Array<OneD, NekDouble> &outarray,
1168 const StdRegions::StdMatrixKey &mkey)
1169{
1170 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1171}
1172
1174 const Array<OneD, const NekDouble> &inarray,
1175 Array<OneD, NekDouble> &outarray,
1176 const StdRegions::StdMatrixKey &mkey)
1177{
1178 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1179}
1180
1182 const Array<OneD, const NekDouble> &inarray,
1183 Array<OneD, NekDouble> &outarray, const StdRegions::StdMatrixKey &mkey)
1184{
1185 StdExpansion::WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
1186}
1187
1189 const Array<OneD, const NekDouble> &inarray,
1190 Array<OneD, NekDouble> &outarray, const StdRegions::StdMatrixKey &mkey)
1191{
1192 StdExpansion::MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
1193}
1194
1196 Array<OneD, NekDouble> &outarray,
1197 const StdRegions::StdMatrixKey &mkey)
1198{
1199 HexExp::v_HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1200}
1201
1202/**
1203 * This function is used to compute exactly the advective numerical flux
1204 * on the interface of two elements with different expansions, hence an
1205 * appropriate number of Gauss points has to be used. The number of
1206 * Gauss points has to be equal to the number used by the highest
1207 * polynomial degree of the two adjacent elements
1208 *
1209 * @param numMin Is the reduced polynomial order
1210 * @param inarray Input array of coefficients
1211 * @param dumpVar Output array of reduced coefficients.
1212 */
1214 const Array<OneD, const NekDouble> &inarray,
1215 Array<OneD, NekDouble> &outarray)
1216{
1217 int n_coeffs = inarray.size();
1218 int nmodes0 = m_base[0]->GetNumModes();
1219 int nmodes1 = m_base[1]->GetNumModes();
1220 int nmodes2 = m_base[2]->GetNumModes();
1221 int numMax = nmodes0;
1222
1223 Array<OneD, NekDouble> coeff(n_coeffs);
1224 Array<OneD, NekDouble> coeff_tmp1(nmodes0 * nmodes1, 0.0);
1225 Array<OneD, NekDouble> coeff_tmp2(n_coeffs, 0.0);
1226 Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
1227
1228 Vmath::Vcopy(n_coeffs, inarray, 1, coeff_tmp2, 1);
1229
1230 const LibUtilities::PointsKey Pkey0(nmodes0,
1232 const LibUtilities::PointsKey Pkey1(nmodes1,
1234 const LibUtilities::PointsKey Pkey2(nmodes2,
1236
1237 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
1238 LibUtilities::BasisKey b1(m_base[1]->GetBasisType(), nmodes1, Pkey1);
1239 LibUtilities::BasisKey b2(m_base[2]->GetBasisType(), nmodes2, Pkey2);
1240 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
1241 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A, nmodes1, Pkey1);
1242 LibUtilities::BasisKey bortho2(LibUtilities::eOrtho_A, nmodes2, Pkey2);
1243
1244 LibUtilities::InterpCoeff3D(b0, b1, b2, coeff_tmp2, bortho0, bortho1,
1245 bortho2, coeff);
1246
1247 Vmath::Zero(n_coeffs, coeff_tmp2, 1);
1248
1249 int cnt = 0, cnt2 = 0;
1250
1251 for (int u = 0; u < numMin + 1; ++u)
1252 {
1253 for (int i = 0; i < numMin; ++i)
1254 {
1255 Vmath::Vcopy(numMin, tmp = coeff + cnt + cnt2, 1,
1256 tmp2 = coeff_tmp1 + cnt, 1);
1257
1258 cnt = i * numMax;
1259 }
1260
1261 Vmath::Vcopy(nmodes0 * nmodes1, tmp3 = coeff_tmp1, 1,
1262 tmp4 = coeff_tmp2 + cnt2, 1);
1263
1264 cnt2 = u * nmodes0 * nmodes1;
1265 }
1266
1267 LibUtilities::InterpCoeff3D(bortho0, bortho1, bortho2, coeff_tmp2, b0, b1,
1268 b2, outarray);
1269}
1270
1272 const StdRegions::StdMatrixKey &mkey)
1273{
1274 int nq = GetTotPoints();
1275
1276 // Calculate sqrt of the Jacobian
1278 Array<OneD, NekDouble> sqrt_jac(nq);
1279 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1280 {
1281 Vmath::Vsqrt(nq, jac, 1, sqrt_jac, 1);
1282 }
1283 else
1284 {
1285 Vmath::Fill(nq, sqrt(jac[0]), sqrt_jac, 1);
1286 }
1287
1288 // Multiply array by sqrt(Jac)
1289 Vmath::Vmul(nq, sqrt_jac, 1, array, 1, array, 1);
1290
1291 // Apply std region filter
1292 StdHexExp::v_SVVLaplacianFilter(array, mkey);
1293
1294 // Divide by sqrt(Jac)
1295 Vmath::Vdiv(nq, array, 1, sqrt_jac, 1, array, 1);
1296}
1297
1298//-----------------------------
1299// Matrix creation functions
1300//-----------------------------
1302{
1303 DNekMatSharedPtr returnval;
1304
1305 switch (mkey.GetMatrixType())
1306 {
1314 returnval = Expansion3D::v_GenMatrix(mkey);
1315 break;
1316 default:
1317 returnval = StdHexExp::v_GenMatrix(mkey);
1318 }
1319
1320 return returnval;
1321}
1322
1324{
1325 LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
1326 LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
1327 LibUtilities::BasisKey bkey2 = m_base[2]->GetBasisKey();
1328
1331
1332 return tmp->GetStdMatrix(mkey);
1333}
1334
1336{
1337 return m_matrixManager[mkey];
1338}
1339
1341{
1342 m_matrixManager.DeleteObject(mkey);
1343}
1344
1346{
1347 return m_staticCondMatrixManager[mkey];
1348}
1349
1351{
1352 m_staticCondMatrixManager.DeleteObject(mkey);
1353}
1354
1356 const Array<OneD, const NekDouble> &inarray,
1358{
1359 // This implementation is only valid when there are no
1360 // coefficients associated to the Laplacian operator
1361 if (m_metrics.count(eMetricLaplacian00) == 0)
1362 {
1364 }
1365
1366 int nquad0 = m_base[0]->GetNumPoints();
1367 int nquad1 = m_base[1]->GetNumPoints();
1368 int nquad2 = m_base[2]->GetNumPoints();
1369 int nqtot = nquad0 * nquad1 * nquad2;
1370
1371 ASSERTL1(wsp.size() >= 6 * nqtot, "Insufficient workspace size.");
1372
1373 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
1374 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
1375 const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
1376 const Array<OneD, const NekDouble> &dbase0 = m_base[0]->GetDbdata();
1377 const Array<OneD, const NekDouble> &dbase1 = m_base[1]->GetDbdata();
1378 const Array<OneD, const NekDouble> &dbase2 = m_base[2]->GetDbdata();
1379 const Array<OneD, const NekDouble> &metric00 =
1381 const Array<OneD, const NekDouble> &metric01 =
1383 const Array<OneD, const NekDouble> &metric02 =
1385 const Array<OneD, const NekDouble> &metric11 =
1387 const Array<OneD, const NekDouble> &metric12 =
1389 const Array<OneD, const NekDouble> &metric22 =
1391
1392 // Allocate temporary storage
1393 Array<OneD, NekDouble> wsp0(wsp);
1394 Array<OneD, NekDouble> wsp1(wsp + 1 * nqtot);
1395 Array<OneD, NekDouble> wsp2(wsp + 2 * nqtot);
1396 Array<OneD, NekDouble> wsp3(wsp + 3 * nqtot);
1397 Array<OneD, NekDouble> wsp4(wsp + 4 * nqtot);
1398 Array<OneD, NekDouble> wsp5(wsp + 5 * nqtot);
1399
1400 StdExpansion3D::PhysTensorDeriv(inarray, wsp0, wsp1, wsp2);
1401
1402 // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1403 // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1404 // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1405 // especially for this purpose
1406 Vmath::Vvtvvtp(nqtot, &metric00[0], 1, &wsp0[0], 1, &metric01[0], 1,
1407 &wsp1[0], 1, &wsp3[0], 1);
1408 Vmath::Vvtvp(nqtot, &metric02[0], 1, &wsp2[0], 1, &wsp3[0], 1, &wsp3[0], 1);
1409 Vmath::Vvtvvtp(nqtot, &metric01[0], 1, &wsp0[0], 1, &metric11[0], 1,
1410 &wsp1[0], 1, &wsp4[0], 1);
1411 Vmath::Vvtvp(nqtot, &metric12[0], 1, &wsp2[0], 1, &wsp4[0], 1, &wsp4[0], 1);
1412 Vmath::Vvtvvtp(nqtot, &metric02[0], 1, &wsp0[0], 1, &metric12[0], 1,
1413 &wsp1[0], 1, &wsp5[0], 1);
1414 Vmath::Vvtvp(nqtot, &metric22[0], 1, &wsp2[0], 1, &wsp5[0], 1, &wsp5[0], 1);
1415
1416 // outarray = m = (D_xi1 * B)^T * k
1417 // wsp1 = n = (D_xi2 * B)^T * l
1418 IProductWRTBase_SumFacKernel(dbase0, base1, base2, wsp3, outarray, wsp0,
1419 false, true, true);
1420 IProductWRTBase_SumFacKernel(base0, dbase1, base2, wsp4, wsp2, wsp0, true,
1421 false, true);
1422 Vmath::Vadd(m_ncoeffs, wsp2.get(), 1, outarray.get(), 1, outarray.get(), 1);
1423 IProductWRTBase_SumFacKernel(base0, base1, dbase2, wsp5, wsp2, wsp0, true,
1424 true, false);
1425 Vmath::Vadd(m_ncoeffs, wsp2.get(), 1, outarray.get(), 1, outarray.get(), 1);
1426}
1427
1429{
1430 if (m_metrics.count(eMetricQuadrature) == 0)
1431 {
1433 }
1434
1435 const SpatialDomains::GeomType type = m_metricinfo->GetGtype();
1436 const unsigned int nqtot = GetTotPoints();
1437 const unsigned int dim = 3;
1438 const MetricType m[3][3] = {
1442
1443 for (unsigned int i = 0; i < dim; ++i)
1444 {
1445 for (unsigned int j = i; j < dim; ++j)
1446 {
1447 m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1448 const Array<TwoD, const NekDouble> &gmat =
1449 m_metricinfo->GetGmat(GetPointsKeys());
1450 if (type == SpatialDomains::eDeformed)
1451 {
1452 Vmath::Vcopy(nqtot, &gmat[i * dim + j][0], 1,
1453 &m_metrics[m[i][j]][0], 1);
1454 }
1455 else
1456 {
1457 Vmath::Fill(nqtot, gmat[i * dim + j][0], &m_metrics[m[i][j]][0],
1458 1);
1459 }
1461 }
1462 }
1463}
1464
1465/** @brief: This method gets all of the factors which are
1466 required as part of the Gradient Jump Penalty
1467 stabilisation and involves the product of the normal and
1468 geometric factors along the element trace.
1469*/
1471 Array<OneD, Array<OneD, NekDouble>> &d0factors,
1472 Array<OneD, Array<OneD, NekDouble>> &d1factors,
1473 Array<OneD, Array<OneD, NekDouble>> &d2factors)
1474{
1475 int nquad0 = GetNumPoints(0);
1476 int nquad1 = GetNumPoints(1);
1477 int nquad2 = GetNumPoints(2);
1478
1480 m_metricinfo->GetDerivFactors(GetPointsKeys());
1481
1482 if (d0factors.size() != 6)
1483 {
1484 d0factors = Array<OneD, Array<OneD, NekDouble>>(6);
1485 d1factors = Array<OneD, Array<OneD, NekDouble>>(6);
1486 d2factors = Array<OneD, Array<OneD, NekDouble>>(6);
1487 }
1488
1489 if (d0factors[0].size() != nquad0 * nquad1)
1490 {
1491 d0factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1492 d0factors[5] = Array<OneD, NekDouble>(nquad0 * nquad1);
1493 d1factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1494 d1factors[5] = Array<OneD, NekDouble>(nquad0 * nquad1);
1495 d2factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1496 d2factors[5] = Array<OneD, NekDouble>(nquad0 * nquad1);
1497 }
1498
1499 if (d0factors[1].size() != nquad0 * nquad2)
1500 {
1501 d0factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1502 d0factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1503 d1factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1504 d1factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1505 d2factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1506 d2factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1507 }
1508
1509 if (d0factors[2].size() != nquad1 * nquad2)
1510 {
1511 d0factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1512 d0factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1513 d1factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1514 d1factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1515 d2factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1516 d2factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1517 }
1518
1519 // Outwards normals
1521 GetTraceNormal(0);
1523 GetTraceNormal(1);
1525 GetTraceNormal(2);
1527 GetTraceNormal(3);
1529 GetTraceNormal(4);
1531 GetTraceNormal(5);
1532
1533 int ncoords = normal_0.size();
1534
1535 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1536 {
1537 // faces 0 and 5
1538 for (int i = 0; i < nquad0 * nquad1; ++i)
1539 {
1540 d0factors[0][i] = df[0][i] * normal_0[0][i];
1541 d1factors[0][i] = df[1][i] * normal_0[0][i];
1542 d2factors[0][i] = df[2][i] * normal_0[0][i];
1543
1544 d0factors[5][i] =
1545 df[0][nquad0 * nquad1 * (nquad2 - 1) + i] * normal_5[0][i];
1546 d1factors[5][i] =
1547 df[1][nquad0 * nquad1 * (nquad2 - 1) + i] * normal_5[0][i];
1548 d2factors[5][i] =
1549 df[2][nquad0 * nquad1 * (nquad2 - 1) + i] * normal_5[0][i];
1550 }
1551
1552 for (int n = 1; n < ncoords; ++n)
1553 {
1554 for (int i = 0; i < nquad0 * nquad1; ++i)
1555 {
1556 d0factors[0][i] += df[3 * n][i] * normal_0[n][i];
1557 d1factors[0][i] += df[3 * n + 1][i] * normal_0[n][i];
1558 d2factors[0][i] += df[3 * n + 2][i] * normal_0[n][i];
1559
1560 d0factors[5][i] +=
1561 df[3 * n][nquad0 * nquad1 * (nquad2 - 1) + i] *
1562 normal_5[n][i];
1563 d1factors[5][i] +=
1564 df[3 * n + 1][nquad0 * nquad1 * (nquad2 - 1) + i] *
1565 normal_5[n][i];
1566 d2factors[5][i] +=
1567 df[3 * n + 2][nquad0 * nquad1 * (nquad2 - 1) + i] *
1568 normal_5[n][i];
1569 }
1570 }
1571
1572 // faces 1 and 3
1573 for (int j = 0; j < nquad2; ++j)
1574 {
1575 for (int i = 0; i < nquad0; ++i)
1576 {
1577 d0factors[1][j * nquad0 + i] = df[0][j * nquad0 * nquad1 + i] *
1578 normal_1[0][j * nquad0 + i];
1579 d1factors[1][j * nquad0 + i] = df[1][j * nquad0 * nquad1 + i] *
1580 normal_1[0][j * nquad0 + i];
1581 d2factors[1][j * nquad0 + i] = df[2][j * nquad0 * nquad1 + i] *
1582 normal_1[0][j * nquad0 + i];
1583
1584 d0factors[3][j * nquad0 + i] =
1585 df[0][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1586 normal_3[0][j * nquad0 + i];
1587 d1factors[3][j * nquad0 + i] =
1588 df[1][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1589 normal_3[0][j * nquad0 + i];
1590 d2factors[3][j * nquad0 + i] =
1591 df[2][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1592 normal_3[0][j * nquad0 + i];
1593 }
1594 }
1595
1596 for (int n = 1; n < ncoords; ++n)
1597 {
1598 for (int j = 0; j < nquad2; ++j)
1599 {
1600 for (int i = 0; i < nquad0; ++i)
1601 {
1602 d0factors[1][j * nquad0 + i] +=
1603 df[3 * n][j * nquad0 * nquad1 + i] *
1604 normal_1[0][j * nquad0 + i];
1605 d1factors[1][j * nquad0 + i] +=
1606 df[3 * n + 1][j * nquad0 * nquad1 + i] *
1607 normal_1[0][j * nquad0 + i];
1608 d2factors[1][j * nquad0 + i] +=
1609 df[3 * n + 2][j * nquad0 * nquad1 + i] *
1610 normal_1[0][j * nquad0 + i];
1611
1612 d0factors[3][j * nquad0 + i] +=
1613 df[3 * n][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1614 normal_3[0][j * nquad0 + i];
1615 d1factors[3][j * nquad0 + i] +=
1616 df[3 * n + 1][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1617 normal_3[0][j * nquad0 + i];
1618 d2factors[3][j * nquad0 + i] +=
1619 df[3 * n + 2][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1620 normal_3[0][j * nquad0 + i];
1621 }
1622 }
1623 }
1624
1625 // faces 2 and 4
1626 for (int j = 0; j < nquad2; ++j)
1627 {
1628 for (int i = 0; i < nquad1; ++i)
1629 {
1630 d0factors[2][j * nquad1 + i] =
1631 df[0][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1632 normal_2[0][j * nquad1 + i];
1633 d1factors[2][j * nquad1 + i] =
1634 df[1][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1635 normal_2[0][j * nquad1 + i];
1636 d2factors[2][j * nquad1 + i] =
1637 df[2][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1638 normal_2[0][j * nquad1 + i];
1639
1640 d0factors[4][j * nquad1 + i] =
1641 df[0][j * nquad0 * nquad1 + i * nquad0] *
1642 normal_4[0][j * nquad1 + i];
1643 d1factors[4][j * nquad1 + i] =
1644 df[1][j * nquad0 * nquad1 + i * nquad0] *
1645 normal_4[0][j * nquad1 + i];
1646 d2factors[4][j * nquad1 + i] =
1647 df[2][j * nquad0 * nquad1 + i * nquad0] *
1648 normal_4[0][j * nquad1 + i];
1649 }
1650 }
1651
1652 for (int n = 1; n < ncoords; ++n)
1653 {
1654 for (int j = 0; j < nquad2; ++j)
1655 {
1656 for (int i = 0; i < nquad1; ++i)
1657 {
1658 d0factors[2][j * nquad1 + i] +=
1659 df[3 * n][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1660 normal_2[n][j * nquad1 + i];
1661 d1factors[2][j * nquad1 + i] +=
1662 df[3 * n + 1]
1663 [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1664 normal_2[n][j * nquad1 + i];
1665 d2factors[2][j * nquad1 + i] +=
1666 df[3 * n + 2]
1667 [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1668 normal_2[n][j * nquad1 + i];
1669
1670 d0factors[4][j * nquad1 + i] +=
1671 df[3 * n][i * nquad0 + j * nquad0 * nquad1] *
1672 normal_4[n][j * nquad1 + i];
1673 d1factors[4][j * nquad1 + i] +=
1674 df[3 * n + 1][i * nquad0 + j * nquad0 * nquad1] *
1675 normal_4[n][j * nquad1 + i];
1676 d2factors[4][j * nquad1 + i] +=
1677 df[3 * n + 2][i * nquad0 + j * nquad0 * nquad1] *
1678 normal_4[n][j * nquad1 + i];
1679 }
1680 }
1681 }
1682 }
1683 else
1684 {
1685 // Faces 0 and 5
1686 for (int i = 0; i < nquad0 * nquad1; ++i)
1687 {
1688 d0factors[0][i] = df[0][0] * normal_0[0][i];
1689 d0factors[5][i] = df[0][0] * normal_5[0][i];
1690
1691 d1factors[0][i] = df[1][0] * normal_0[0][i];
1692 d1factors[5][i] = df[1][0] * normal_5[0][i];
1693
1694 d2factors[0][i] = df[2][0] * normal_0[0][i];
1695 d2factors[5][i] = df[2][0] * normal_5[0][i];
1696 }
1697
1698 for (int n = 1; n < ncoords; ++n)
1699 {
1700 for (int i = 0; i < nquad0 * nquad1; ++i)
1701 {
1702 d0factors[0][i] += df[3 * n][0] * normal_0[n][i];
1703 d0factors[5][i] += df[3 * n][0] * normal_5[n][i];
1704
1705 d1factors[0][i] += df[3 * n + 1][0] * normal_0[n][i];
1706 d1factors[5][i] += df[3 * n + 1][0] * normal_5[n][i];
1707
1708 d2factors[0][i] += df[3 * n + 2][0] * normal_0[n][i];
1709 d2factors[5][i] += df[3 * n + 2][0] * normal_5[n][i];
1710 }
1711 }
1712
1713 // faces 1 and 3
1714 for (int i = 0; i < nquad0 * nquad2; ++i)
1715 {
1716 d0factors[1][i] = df[0][0] * normal_1[0][i];
1717 d0factors[3][i] = df[0][0] * normal_3[0][i];
1718
1719 d1factors[1][i] = df[1][0] * normal_1[0][i];
1720 d1factors[3][i] = df[1][0] * normal_3[0][i];
1721
1722 d2factors[1][i] = df[2][0] * normal_1[0][i];
1723 d2factors[3][i] = df[2][0] * normal_3[0][i];
1724 }
1725
1726 for (int n = 1; n < ncoords; ++n)
1727 {
1728 for (int i = 0; i < nquad0 * nquad2; ++i)
1729 {
1730 d0factors[1][i] += df[3 * n][0] * normal_1[n][i];
1731 d0factors[3][i] += df[3 * n][0] * normal_3[n][i];
1732
1733 d1factors[1][i] += df[3 * n + 1][0] * normal_1[n][i];
1734 d1factors[3][i] += df[3 * n + 1][0] * normal_3[n][i];
1735
1736 d2factors[1][i] += df[3 * n + 2][0] * normal_1[n][i];
1737 d2factors[3][i] += df[3 * n + 2][0] * normal_3[n][i];
1738 }
1739 }
1740
1741 // faces 2 and 4
1742 for (int i = 0; i < nquad1 * nquad2; ++i)
1743 {
1744 d0factors[2][i] = df[0][0] * normal_2[0][i];
1745 d0factors[4][i] = df[0][0] * normal_4[0][i];
1746
1747 d1factors[2][i] = df[1][0] * normal_2[0][i];
1748 d1factors[4][i] = df[1][0] * normal_4[0][i];
1749
1750 d2factors[2][i] = df[2][0] * normal_2[0][i];
1751 d2factors[4][i] = df[2][0] * normal_4[0][i];
1752 }
1753
1754 for (int n = 1; n < ncoords; ++n)
1755 {
1756 for (int i = 0; i < nquad1 * nquad2; ++i)
1757 {
1758 d0factors[2][i] += df[3 * n][0] * normal_2[n][i];
1759 d0factors[4][i] += df[3 * n][0] * normal_4[n][i];
1760
1761 d1factors[2][i] += df[3 * n + 1][0] * normal_2[n][i];
1762 d1factors[4][i] += df[3 * n + 1][0] * normal_4[n][i];
1763
1764 d2factors[2][i] += df[3 * n + 2][0] * normal_2[n][i];
1765 d2factors[4][i] += df[3 * n + 2][0] * normal_4[n][i];
1766 }
1767 }
1768 }
1769}
1770} // namespace LocalRegions
1771} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
Describes the specification for a Basis.
Definition: Basis.h:47
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:122
PointsKey GetPointsKey() const
Return distribution of points.
Definition: Basis.h:139
Defines a specification for a set of points.
Definition: Points.h:55
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
std::map< int, NormalVector > m_traceNormals
Definition: Expansion.h:278
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
Definition: Expansion.h:288
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:171
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:275
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
Definition: Expansion.cpp:608
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: Expansion.cpp:535
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:276
const NormalVector & GetTraceNormal(const int id)
Definition: Expansion.cpp:255
virtual void v_WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: HexExp.cpp:1181
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrate the physical point list inarray over region.
Definition: HexExp.cpp:103
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: HexExp.cpp:647
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: HexExp.h:246
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: HexExp.cpp:1151
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: HexExp.cpp:1158
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
Definition: HexExp.cpp:1301
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdRegions::StdMatrixKey &mkey) override
Definition: HexExp.cpp:1323
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: HexExp.h:248
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: HexExp.cpp:1195
void v_DropLocStaticCondMatrix(const MatrixKey &mkey) override
Definition: HexExp.cpp:1350
virtual DNekScalMatSharedPtr v_GetLocMatrix(const MatrixKey &mkey) override
Definition: HexExp.cpp:1335
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Calculate the inner product of inarray with respect to the given basis B = base0 * base1 * base2.
Definition: HexExp.cpp:372
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
Definition: HexExp.cpp:299
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculates the inner product .
Definition: HexExp.cpp:429
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey) override
Definition: HexExp.cpp:1271
void v_ComputeTraceNormal(const int face) override
Definition: HexExp.cpp:882
HexExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, const SpatialDomains::HexGeomSharedPtr &geom)
Constructor using BasisKey class for quadrature points and order definition.
Definition: HexExp.cpp:59
virtual void v_ComputeLaplacianMetric() override
Definition: HexExp.cpp:1428
virtual void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType) override
Definition: HexExp.cpp:664
virtual void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp) override
Definition: HexExp.cpp:1355
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
Calculate the derivative of the physical points.
Definition: HexExp.cpp:144
void v_DropLocMatrix(const MatrixKey &mkey) override
Definition: HexExp.cpp:1340
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: HexExp.cpp:402
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
Retrieves the physical coordinates of a given set of reference coordinates.
Definition: HexExp.cpp:630
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
This function evaluates the expansion at a single (arbitrary) point of the domain.
Definition: HexExp.cpp:583
virtual void v_IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: HexExp.cpp:516
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
Definition: HexExp.cpp:474
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &out) override
Physical derivative along a direction vector.
Definition: HexExp.cpp:257
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
Definition: HexExp.cpp:575
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const MatrixKey &mkey) override
Definition: HexExp.cpp:1345
virtual void v_GetTracePhysMap(const int face, Array< OneD, int > &outarray) override
Definition: HexExp.cpp:767
virtual LibUtilities::ShapeType v_DetShapeType() const override
Return the region shape using the enum-list of ShapeType.
Definition: HexExp.cpp:659
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: HexExp.cpp:1213
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: HexExp.cpp:1173
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const override
Definition: HexExp.cpp:603
virtual StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const override
Definition: HexExp.cpp:610
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculate the inner product of inarray with respect to the elements basis.
Definition: HexExp.cpp:334
virtual void v_NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors) override
: This method gets all of the factors which are required as part of the Gradient Jump Penalty stabili...
Definition: HexExp.cpp:1470
virtual void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: HexExp.cpp:1188
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)
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:162
const LibUtilities::PointsKeyVector GetPointsKeys() const
const LibUtilities::BasisKey GetTraceBasisKey(const int i, int k=-1) const
This function returns the basis key belonging to the i-th trace.
Definition: StdExpansion.h:305
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
this function calculates the inner product of a given function f with the different modes of the expa...
Definition: StdExpansion.h:534
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
Definition: StdExpansion.h:430
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:211
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:224
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:729
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:855
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Class representing a hexehedral element in reference space.
Definition: StdHexExp.h:48
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:87
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:137
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis,...
Definition: Interp.cpp:164
void InterpCoeff3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
Definition: Interp.cpp:103
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:236
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< HexGeom > HexGeomSharedPtr
Definition: HexGeom.h:87
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::shared_ptr< StdHexExp > StdHexExpSharedPtr
Definition: StdHexExp.h:255
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:529
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:207
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:569
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:354
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:245
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
Definition: Vmath.cpp:319
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.cpp:280
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:43
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.cpp:687
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294