Nektar++
QuadExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: QuadExp.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: Expansion for quadrilateral elements.
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <boost/core/ignore_unused.hpp>
36
43#include <LocalRegions/SegExp.h>
44
45using namespace std;
46
47namespace Nektar
48{
49namespace LocalRegions
50{
52 const LibUtilities::BasisKey &Bb,
54 : StdExpansion(Ba.GetNumModes() * Bb.GetNumModes(), 2, Ba, Bb),
55 StdExpansion2D(Ba.GetNumModes() * Bb.GetNumModes(), Ba, Bb),
56 StdQuadExp(Ba, Bb), Expansion(geom), Expansion2D(geom),
57 m_matrixManager(
58 std::bind(&Expansion2D::CreateMatrix, this, std::placeholders::_1),
59 std::string("QuadExpMatrix")),
60 m_staticCondMatrixManager(std::bind(&Expansion::CreateStaticCondMatrix,
61 this, std::placeholders::_1),
62 std::string("QuadExpStaticCondMatrix"))
63{
64}
65
67 : StdExpansion(T), StdExpansion2D(T), StdQuadExp(T), Expansion(T),
68 Expansion2D(T), m_matrixManager(T.m_matrixManager),
69 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
70{
71}
72
74{
75 int nquad0 = m_base[0]->GetNumPoints();
76 int nquad1 = m_base[1]->GetNumPoints();
78 NekDouble ival;
79 Array<OneD, NekDouble> tmp(nquad0 * nquad1);
80
81 // multiply inarray with Jacobian
82 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
83 {
84 Vmath::Vmul(nquad0 * nquad1, jac, 1, inarray, 1, tmp, 1);
85 }
86 else
87 {
88 Vmath::Smul(nquad0 * nquad1, jac[0], inarray, 1, tmp, 1);
89 }
90
91 // call StdQuadExp version;
92 ival = StdQuadExp::v_Integral(tmp);
93 return ival;
94}
95
100{
101 int nquad0 = m_base[0]->GetNumPoints();
102 int nquad1 = m_base[1]->GetNumPoints();
103 int nqtot = nquad0 * nquad1;
105 m_metricinfo->GetDerivFactors(GetPointsKeys());
106 Array<OneD, NekDouble> diff0(2 * nqtot);
107 Array<OneD, NekDouble> diff1(diff0 + nqtot);
108
109 StdQuadExp::v_PhysDeriv(inarray, diff0, diff1);
110
111 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
112 {
113 if (out_d0.size())
114 {
115 Vmath::Vmul(nqtot, df[0], 1, diff0, 1, out_d0, 1);
116 Vmath::Vvtvp(nqtot, df[1], 1, diff1, 1, out_d0, 1, out_d0, 1);
117 }
118
119 if (out_d1.size())
120 {
121 Vmath::Vmul(nqtot, df[2], 1, diff0, 1, out_d1, 1);
122 Vmath::Vvtvp(nqtot, df[3], 1, diff1, 1, out_d1, 1, out_d1, 1);
123 }
124
125 if (out_d2.size())
126 {
127 Vmath::Vmul(nqtot, df[4], 1, diff0, 1, out_d2, 1);
128 Vmath::Vvtvp(nqtot, df[5], 1, diff1, 1, out_d2, 1, out_d2, 1);
129 }
130 }
131 else // regular geometry
132 {
133 if (out_d0.size())
134 {
135 Vmath::Smul(nqtot, df[0][0], diff0, 1, out_d0, 1);
136 Blas::Daxpy(nqtot, df[1][0], diff1, 1, out_d0, 1);
137 }
138
139 if (out_d1.size())
140 {
141 Vmath::Smul(nqtot, df[2][0], diff0, 1, out_d1, 1);
142 Blas::Daxpy(nqtot, df[3][0], diff1, 1, out_d1, 1);
143 }
144
145 if (out_d2.size())
146 {
147 Vmath::Smul(nqtot, df[4][0], diff0, 1, out_d2, 1);
148 Blas::Daxpy(nqtot, df[5][0], diff1, 1, out_d2, 1);
149 }
150 }
151}
152
153void QuadExp::v_PhysDeriv(const int dir,
154 const Array<OneD, const NekDouble> &inarray,
155 Array<OneD, NekDouble> &outarray)
156{
157 switch (dir)
158 {
159 case 0:
160 {
161 v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
163 }
164 break;
165 case 1:
166 {
167 v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
169 }
170 break;
171 case 2:
172 {
174 outarray);
175 }
176 break;
177 default:
178 {
179 ASSERTL1(false, "input dir is out of range");
180 }
181 break;
182 }
183}
184
186 const Array<OneD, const NekDouble> &inarray,
188{
189 int nquad0 = m_base[0]->GetNumPoints();
190 int nquad1 = m_base[1]->GetNumPoints();
191 int nqtot = nquad0 * nquad1;
192
194 m_metricinfo->GetDerivFactors(GetPointsKeys());
195
196 Array<OneD, NekDouble> diff0(2 * nqtot);
197 Array<OneD, NekDouble> diff1(diff0 + nqtot);
198
199 StdQuadExp::v_PhysDeriv(inarray, diff0, diff1);
200
201 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
202 {
204
205 // d/dx_v^s = v_x*ds/dx + v_y*ds/dy + v_z*dx/dz
206 for (int i = 0; i < 2; ++i)
207 {
208 tangmat[i] = Array<OneD, NekDouble>(nqtot, 0.0);
209 for (int k = 0; k < (m_geom->GetCoordim()); ++k)
210 {
211 Vmath::Vvtvp(nqtot, &df[2 * k + i][0], 1, &direction[k * nqtot],
212 1, &tangmat[i][0], 1, &tangmat[i][0], 1);
213 }
214 }
215
216 /// D_v = d/dx_v^s + d/dx_v^r
217 if (out.size())
218 {
219 Vmath::Vmul(nqtot, &tangmat[0][0], 1, &diff0[0], 1, &out[0], 1);
220 Vmath::Vvtvp(nqtot, &tangmat[1][0], 1, &diff1[0], 1, &out[0], 1,
221 &out[0], 1);
222 }
223 }
224 else
225 {
227 "Wrong route");
228 }
229}
230
232 Array<OneD, NekDouble> &outarray)
233{
234 if ((m_base[0]->Collocation()) && (m_base[1]->Collocation()))
235 {
236 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
237 }
238 else
239 {
240 IProductWRTBase(inarray, outarray);
241
242 // get Mass matrix inverse
244 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
245
246 // copy inarray in case inarray == outarray
247 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
249
250 out = (*matsys) * in;
251 }
252}
253
255 const Array<OneD, const NekDouble> &inarray,
256 Array<OneD, NekDouble> &outarray)
257{
258 if ((m_base[0]->Collocation()) && (m_base[1]->Collocation()))
259 {
260 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
261 }
262 else
263 {
264 int i, j;
265 int npoints[2] = {m_base[0]->GetNumPoints(), m_base[1]->GetNumPoints()};
266 int nmodes[2] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes()};
267
268 fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
269
270 if (nmodes[0] == 1 && nmodes[1] == 1)
271 {
272 outarray[0] = inarray[0];
273 return;
274 }
275
276 Array<OneD, NekDouble> physEdge[4];
277 Array<OneD, NekDouble> coeffEdge[4];
278 StdRegions::Orientation orient[4];
279 for (i = 0; i < 4; i++)
280 {
281 physEdge[i] = Array<OneD, NekDouble>(npoints[i % 2]);
282 coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i % 2]);
283 orient[i] = GetTraceOrient(i);
284 }
285
286 for (i = 0; i < npoints[0]; i++)
287 {
288 physEdge[0][i] = inarray[i];
289 physEdge[2][i] = inarray[npoints[0] * (npoints[1] - 1) + i];
290 }
291
292 for (i = 0; i < npoints[1]; i++)
293 {
294 physEdge[1][i] = inarray[npoints[0] - 1 + i * npoints[0]];
295 physEdge[3][i] = inarray[i * npoints[0]];
296 }
297
298 for (i = 0; i < 4; i++)
299 {
300 if (orient[i] == StdRegions::eBackwards)
301 {
302 reverse((physEdge[i]).get(),
303 (physEdge[i]).get() + npoints[i % 2]);
304 }
305 }
306
307 SegExpSharedPtr segexp[4];
308 for (i = 0; i < 4; i++)
309 {
311 m_base[i % 2]->GetBasisKey(), GetGeom2D()->GetEdge(i));
312 }
313
315 Array<OneD, int> signArray;
317
318 for (i = 0; i < 4; i++)
319 {
320 segexp[i % 2]->FwdTransBndConstrained(physEdge[i], coeffEdge[i]);
321
322 GetTraceToElementMap(i, mapArray, signArray, orient[i]);
323 for (j = 0; j < nmodes[i % 2]; j++)
324 {
325 sign = (NekDouble)signArray[j];
326 outarray[mapArray[j]] = sign * coeffEdge[i][j];
327 }
328 }
329
330 int nBoundaryDofs = NumBndryCoeffs();
331 int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
332
333 if (nInteriorDofs > 0)
334 {
337
339 DetShapeType(), *this);
340 MassMatrixOp(outarray, tmp0, stdmasskey);
341 IProductWRTBase(inarray, tmp1);
342
343 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
344
345 // get Mass matrix inverse (only of interior DOF)
346 // use block (1,1) of the static condensed system
347 // note: this block alreay contains the inverse matrix
348 MatrixKey masskey(StdRegions::eMass, DetShapeType(), *this);
349 DNekScalMatSharedPtr matsys =
350 (m_staticCondMatrixManager[masskey])->GetBlock(1, 1);
351
352 Array<OneD, NekDouble> rhs(nInteriorDofs);
353 Array<OneD, NekDouble> result(nInteriorDofs);
354
355 GetInteriorMap(mapArray);
356
357 for (i = 0; i < nInteriorDofs; i++)
358 {
359 rhs[i] = tmp1[mapArray[i]];
360 }
361
362 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, matsys->Scale(),
363 &((matsys->GetOwnedMatrix())->GetPtr())[0],
364 nInteriorDofs, rhs.get(), 1, 0.0, result.get(), 1);
365
366 for (i = 0; i < nInteriorDofs; i++)
367 {
368 outarray[mapArray[i]] = result[i];
369 }
370 }
371 }
372}
373
375 Array<OneD, NekDouble> &outarray)
376{
377 if (m_base[0]->Collocation() && m_base[1]->Collocation())
378 {
379 MultiplyByQuadratureMetric(inarray, outarray);
380 }
381 else
382 {
383 IProductWRTBase_SumFac(inarray, outarray);
384 }
385}
386
388 const int dir, const Array<OneD, const NekDouble> &inarray,
389 Array<OneD, NekDouble> &outarray)
390{
391 IProductWRTDerivBase_SumFac(dir, inarray, outarray);
392}
393
395 const Array<OneD, const NekDouble> &inarray,
396 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
397{
398 int nquad0 = m_base[0]->GetNumPoints();
399 int nquad1 = m_base[1]->GetNumPoints();
400 int order0 = m_base[0]->GetNumModes();
401
402 if (multiplybyweights)
403 {
404 Array<OneD, NekDouble> tmp(nquad0 * nquad1 + nquad1 * order0);
405 Array<OneD, NekDouble> wsp(tmp + nquad0 * nquad1);
406
407 MultiplyByQuadratureMetric(inarray, tmp);
408 StdQuadExp::IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
409 m_base[1]->GetBdata(), tmp,
410 outarray, wsp, true, true);
411 }
412 else
413 {
414 Array<OneD, NekDouble> wsp(nquad1 * order0);
415
416 StdQuadExp::IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
417 m_base[1]->GetBdata(), inarray,
418 outarray, wsp, true, true);
419 }
420}
421
423 const int dir, const Array<OneD, const NekDouble> &inarray,
424 Array<OneD, NekDouble> &outarray)
425{
426 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
427 ASSERTL1((dir == 2) ? (m_geom->GetCoordim() == 3) : true,
428 "Invalid direction.");
429
430 int nquad0 = m_base[0]->GetNumPoints();
431 int nquad1 = m_base[1]->GetNumPoints();
432 int nqtot = nquad0 * nquad1;
433 int nmodes0 = m_base[0]->GetNumModes();
434
435 Array<OneD, NekDouble> tmp1(2 * nqtot + m_ncoeffs + nmodes0 * nquad1);
436 Array<OneD, NekDouble> tmp2(tmp1 + nqtot);
437 Array<OneD, NekDouble> tmp3(tmp1 + 2 * nqtot);
438 Array<OneD, NekDouble> tmp4(tmp1 + 2 * nqtot + m_ncoeffs);
439
441 tmp2D[0] = tmp1;
442 tmp2D[1] = tmp2;
443
444 QuadExp::v_AlignVectorToCollapsedDir(dir, inarray, tmp2D);
445
446 MultiplyByQuadratureMetric(tmp1, tmp1);
447 MultiplyByQuadratureMetric(tmp2, tmp2);
448
449 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
450 tmp1, tmp3, tmp4, false, true);
451 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
452 tmp2, outarray, tmp4, true, false);
453 Vmath::Vadd(m_ncoeffs, tmp3, 1, outarray, 1, outarray, 1);
454}
455
457 const int dir, const Array<OneD, const NekDouble> &inarray,
459{
460 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
461 ASSERTL1((dir == 2) ? (m_geom->GetCoordim() == 3) : true,
462 "Invalid direction.");
463
464 int nquad0 = m_base[0]->GetNumPoints();
465 int nquad1 = m_base[1]->GetNumPoints();
466 int nqtot = nquad0 * nquad1;
467 int nmodes0 = m_base[0]->GetNumModes();
468
470 m_metricinfo->GetDerivFactors(GetPointsKeys());
471
472 Array<OneD, NekDouble> tmp1 = outarray[0];
473 Array<OneD, NekDouble> tmp2 = outarray[1];
475 Array<OneD, NekDouble> tmp4(nmodes0 * nquad1);
476
477 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
478 {
479 Vmath::Vmul(nqtot, &df[2 * dir][0], 1, inarray.get(), 1, tmp1.get(), 1);
480 Vmath::Vmul(nqtot, &df[2 * dir + 1][0], 1, inarray.get(), 1, tmp2.get(),
481 1);
482 }
483 else
484 {
485 Vmath::Smul(nqtot, df[2 * dir][0], inarray.get(), 1, tmp1.get(), 1);
486 Vmath::Smul(nqtot, df[2 * dir + 1][0], inarray.get(), 1, tmp2.get(), 1);
487 }
488}
489
494{
495 int nq = m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints();
497
499 GetLeftAdjacentElementExp()->GetTraceNormal(
501
502 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
503 {
504 Vmath::Vvtvvtp(nq, &normals[0][0], 1, &Fx[0], 1, &normals[1][0], 1,
505 &Fy[0], 1, &Fn[0], 1);
506 Vmath::Vvtvp(nq, &normals[2][0], 1, &Fz[0], 1, &Fn[0], 1, &Fn[0], 1);
507 }
508 else
509 {
510 Vmath::Svtsvtp(nq, normals[0][0], &Fx[0], 1, normals[1][0], &Fy[0], 1,
511 &Fn[0], 1);
512 Vmath::Svtvp(nq, normals[2][0], &Fz[0], 1, &Fn[0], 1, &Fn[0], 1);
513 }
514
515 IProductWRTBase(Fn, outarray);
516}
517
519 const Array<OneD, const Array<OneD, NekDouble>> &Fvec,
520 Array<OneD, NekDouble> &outarray)
521{
522 NormVectorIProductWRTBase(Fvec[0], Fvec[1], Fvec[2], outarray);
523}
524
526{
528 m_base[0]->GetBasisKey(), m_base[1]->GetBasisKey());
529}
530
532{
534 m_base[0]->GetPointsKey());
536 m_base[1]->GetPointsKey());
537
539 bkey1);
540}
541
543 Array<OneD, NekDouble> &coords_1,
544 Array<OneD, NekDouble> &coords_2)
545{
546 Expansion::v_GetCoords(coords_0, coords_1, coords_2);
547}
548
551{
552 int i;
553
554 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[1] <= 1.0 && Lcoords[1] >= -1.0 &&
555 Lcoords[1] <= 1.0,
556 "Local coordinates are not in region [-1,1]");
557
558 m_geom->FillGeom();
559 for (i = 0; i < m_geom->GetCoordim(); ++i)
560 {
561 coords[i] = m_geom->GetCoord(i, Lcoords);
562 }
563}
564
565/**
566 * Given the local cartesian coordinate \a Lcoord evaluate the
567 * value of physvals at this point by calling through to the
568 * StdExpansion method
569 */
571 const Array<OneD, const NekDouble> &Lcoord,
572 const Array<OneD, const NekDouble> &physvals)
573{
574 // Evaluate point in local (eta) coordinates.
575 return StdExpansion2D::v_PhysEvaluate(Lcoord, physvals);
576}
577
579 const Array<OneD, const NekDouble> &physvals)
580{
582
583 ASSERTL0(m_geom, "m_geom not defined");
584 m_geom->GetLocCoords(coord, Lcoord);
585
586 return StdExpansion2D::v_PhysEvaluate(Lcoord, physvals);
587}
588
590 const Array<OneD, const NekDouble> &inarray,
591 std::array<NekDouble, 3> &firstOrderDerivs)
592{
593 Array<OneD, NekDouble> Lcoord(2);
594 ASSERTL0(m_geom, "m_geom not defined");
595 m_geom->GetLocCoords(coord, Lcoord);
596 return StdQuadExp::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
597}
598
600 const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp,
601 const Array<OneD, const NekDouble> &inarray,
603{
604 int nquad0 = m_base[0]->GetNumPoints();
605 int nquad1 = m_base[1]->GetNumPoints();
606
607 // Implementation for all the basis except Gauss points
610 {
611 switch (edge)
612 {
613 case 0:
614 Vmath::Vcopy(nquad0, &(inarray[0]), 1, &(outarray[0]), 1);
615 break;
616 case 1:
617 Vmath::Vcopy(nquad1, &(inarray[0]) + (nquad0 - 1), nquad0,
618 &(outarray[0]), 1);
619 break;
620 case 2:
621 Vmath::Vcopy(nquad0, &(inarray[0]) + nquad0 * (nquad1 - 1), 1,
622 &(outarray[0]), 1);
623 break;
624 case 3:
625 Vmath::Vcopy(nquad1, &(inarray[0]), nquad0, &(outarray[0]), 1);
626 break;
627 default:
628 ASSERTL0(false, "edge value (< 3) is out of range");
629 break;
630 }
631 }
632 else
633 {
634 QuadExp::GetEdgeInterpVals(edge, inarray, outarray);
635 }
636
637 // Interpolate if required
638 if (m_base[edge % 2]->GetPointsKey() !=
639 EdgeExp->GetBasis(0)->GetPointsKey())
640 {
641 Array<OneD, NekDouble> outtmp(max(nquad0, nquad1));
642
643 outtmp = outarray;
644
645 LibUtilities::Interp1D(m_base[edge % 2]->GetPointsKey(), outtmp,
646 EdgeExp->GetBasis(0)->GetPointsKey(), outarray);
647 }
648
649 if (orient == StdRegions::eNoOrientation)
650 {
651 orient = GetTraceOrient(edge);
652 }
653 // Reverse data if necessary
654 if (orient == StdRegions::eBackwards)
655 {
656 Vmath::Reverse(EdgeExp->GetNumPoints(0), &outarray[0], 1, &outarray[0],
657 1);
658 }
659}
660
661void QuadExp::GetEdgeInterpVals(const int edge,
662 const Array<OneD, const NekDouble> &inarray,
663 Array<OneD, NekDouble> &outarray)
664{
665 int i;
666 int nq0 = m_base[0]->GetNumPoints();
667 int nq1 = m_base[1]->GetNumPoints();
668
671
673 *this, factors);
674
675 DNekScalMatSharedPtr mat_gauss = m_matrixManager[key];
676
677 switch (edge)
678 {
679 case 0:
680 {
681 for (i = 0; i < nq0; i++)
682 {
683 outarray[i] =
684 Blas::Ddot(nq1, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
685 1, &inarray[i], nq0);
686 }
687 break;
688 }
689 case 1:
690 {
691 for (i = 0; i < nq1; i++)
692 {
693 outarray[i] =
694 Blas::Ddot(nq0, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
695 1, &inarray[i * nq0], 1);
696 }
697 break;
698 }
699 case 2:
700 {
701 for (i = 0; i < nq0; i++)
702 {
703 outarray[i] =
704 Blas::Ddot(nq1, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
705 1, &inarray[i], nq0);
706 }
707 break;
708 }
709 case 3:
710 {
711 for (i = 0; i < nq1; i++)
712 {
713 outarray[i] =
714 Blas::Ddot(nq0, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
715 1, &inarray[i * nq0], 1);
716 }
717 break;
718 }
719 default:
720 ASSERTL0(false, "edge value (< 3) is out of range");
721 break;
722 }
723}
724
725void QuadExp::v_GetTracePhysMap(const int edge, Array<OneD, int> &outarray)
726{
727 int nquad0 = m_base[0]->GetNumPoints();
728 int nquad1 = m_base[1]->GetNumPoints();
729
730 // Get points in Cartesian orientation
731 switch (edge)
732 {
733 case 0:
734 outarray = Array<OneD, int>(nquad0);
735 for (int i = 0; i < nquad0; ++i)
736 {
737 outarray[i] = i;
738 }
739 break;
740 case 1:
741 outarray = Array<OneD, int>(nquad1);
742 for (int i = 0; i < nquad1; ++i)
743 {
744 outarray[i] = (nquad0 - 1) + i * nquad0;
745 }
746 break;
747 case 2:
748 outarray = Array<OneD, int>(nquad0);
749 for (int i = 0; i < nquad0; ++i)
750 {
751 outarray[i] = i + nquad0 * (nquad1 - 1);
752 }
753 break;
754 case 3:
755 outarray = Array<OneD, int>(nquad1);
756 for (int i = 0; i < nquad1; ++i)
757 {
758 outarray[i] = i * nquad0;
759 }
760 break;
761 default:
762 ASSERTL0(false, "edge value (< 3) is out of range");
763 break;
764 }
765}
766
767void QuadExp::v_GetTraceQFactors(const int edge,
768 Array<OneD, NekDouble> &outarray)
769{
770 int i;
771 int nquad0 = m_base[0]->GetNumPoints();
772 int nquad1 = m_base[1]->GetNumPoints();
773
775 const Array<OneD, const NekDouble> &jac = m_metricinfo->GetJac(ptsKeys);
777 m_metricinfo->GetDerivFactors(ptsKeys);
778
779 Array<OneD, NekDouble> j(max(nquad0, nquad1), 0.0);
780 Array<OneD, NekDouble> g0(max(nquad0, nquad1), 0.0);
781 Array<OneD, NekDouble> g1(max(nquad0, nquad1), 0.0);
782 Array<OneD, NekDouble> g2(max(nquad0, nquad1), 0.0);
783 Array<OneD, NekDouble> g3(max(nquad0, nquad1), 0.0);
784
785 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
786 {
787 // Implementation for all the basis except Gauss points
790 {
791 switch (edge)
792 {
793 case 0:
794 Vmath::Vcopy(nquad0, &(df[1][0]), 1, &(g1[0]), 1);
795 Vmath::Vcopy(nquad0, &(df[3][0]), 1, &(g3[0]), 1);
796 Vmath::Vcopy(nquad0, &(jac[0]), 1, &(j[0]), 1);
797
798 for (i = 0; i < nquad0; ++i)
799 {
800 outarray[i] =
801 j[i] * sqrt(g1[i] * g1[i] + g3[i] * g3[i]);
802 }
803 break;
804 case 1:
805 Vmath::Vcopy(nquad1, &(df[0][0]) + (nquad0 - 1), nquad0,
806 &(g0[0]), 1);
807
808 Vmath::Vcopy(nquad1, &(df[2][0]) + (nquad0 - 1), nquad0,
809 &(g2[0]), 1);
810
811 Vmath::Vcopy(nquad1, &(jac[0]) + (nquad0 - 1), nquad0,
812 &(j[0]), 1);
813
814 for (i = 0; i < nquad1; ++i)
815 {
816 outarray[i] =
817 j[i] * sqrt(g0[i] * g0[i] + g2[i] * g2[i]);
818 }
819 break;
820 case 2:
821
822 Vmath::Vcopy(nquad0, &(df[1][0]) + (nquad0 * (nquad1 - 1)),
823 1, &(g1[0]), 1);
824
825 Vmath::Vcopy(nquad0, &(df[3][0]) + (nquad0 * (nquad1 - 1)),
826 1, &(g3[0]), 1);
827
828 Vmath::Vcopy(nquad0, &(jac[0]) + (nquad0 * (nquad1 - 1)), 1,
829 &(j[0]), 1);
830
831 for (i = 0; i < nquad0; ++i)
832 {
833 outarray[i] =
834 j[i] * sqrt(g1[i] * g1[i] + g3[i] * g3[i]);
835 }
836 break;
837 case 3:
838
839 Vmath::Vcopy(nquad1, &(df[0][0]), nquad0, &(g0[0]), 1);
840 Vmath::Vcopy(nquad1, &(df[2][0]), nquad0, &(g2[0]), 1);
841 Vmath::Vcopy(nquad1, &(jac[0]), nquad0, &(j[0]), 1);
842
843 for (i = 0; i < nquad1; ++i)
844 {
845 outarray[i] =
846 j[i] * sqrt(g0[i] * g0[i] + g2[i] * g2[i]);
847 }
848 break;
849 default:
850 ASSERTL0(false, "edge value (< 3) is out of range");
851 break;
852 }
853 }
854 else
855 {
856 int nqtot = nquad0 * nquad1;
857 Array<OneD, NekDouble> tmp_gmat0(nqtot, 0.0);
858 Array<OneD, NekDouble> tmp_gmat1(nqtot, 0.0);
859 Array<OneD, NekDouble> tmp_gmat2(nqtot, 0.0);
860 Array<OneD, NekDouble> tmp_gmat3(nqtot, 0.0);
861 Array<OneD, NekDouble> g0_edge(max(nquad0, nquad1), 0.0);
862 Array<OneD, NekDouble> g1_edge(max(nquad0, nquad1), 0.0);
863 Array<OneD, NekDouble> g2_edge(max(nquad0, nquad1), 0.0);
864 Array<OneD, NekDouble> g3_edge(max(nquad0, nquad1), 0.0);
865 Array<OneD, NekDouble> jac_edge(max(nquad0, nquad1), 0.0);
866
867 switch (edge)
868 {
869 case 0:
870 Vmath::Vmul(nqtot, &(df[1][0]), 1, &jac[0], 1,
871 &(tmp_gmat1[0]), 1);
872 Vmath::Vmul(nqtot, &(df[3][0]), 1, &jac[0], 1,
873 &(tmp_gmat3[0]), 1);
874 QuadExp::GetEdgeInterpVals(edge, tmp_gmat1, g1_edge);
875 QuadExp::GetEdgeInterpVals(edge, tmp_gmat3, g3_edge);
876
877 for (i = 0; i < nquad0; ++i)
878 {
879 outarray[i] = sqrt(g1_edge[i] * g1_edge[i] +
880 g3_edge[i] * g3_edge[i]);
881 }
882 break;
883
884 case 1:
885 Vmath::Vmul(nqtot, &(df[0][0]), 1, &jac[0], 1,
886 &(tmp_gmat0[0]), 1);
887 Vmath::Vmul(nqtot, &(df[2][0]), 1, &jac[0], 1,
888 &(tmp_gmat2[0]), 1);
889 QuadExp::GetEdgeInterpVals(edge, tmp_gmat0, g0_edge);
890 QuadExp::GetEdgeInterpVals(edge, tmp_gmat2, g2_edge);
891
892 for (i = 0; i < nquad1; ++i)
893 {
894 outarray[i] = sqrt(g0_edge[i] * g0_edge[i] +
895 g2_edge[i] * g2_edge[i]);
896 }
897
898 break;
899 case 2:
900
901 Vmath::Vmul(nqtot, &(df[1][0]), 1, &jac[0], 1,
902 &(tmp_gmat1[0]), 1);
903 Vmath::Vmul(nqtot, &(df[3][0]), 1, &jac[0], 1,
904 &(tmp_gmat3[0]), 1);
905 QuadExp::GetEdgeInterpVals(edge, tmp_gmat1, g1_edge);
906 QuadExp::GetEdgeInterpVals(edge, tmp_gmat3, g3_edge);
907
908 for (i = 0; i < nquad0; ++i)
909 {
910 outarray[i] = sqrt(g1_edge[i] * g1_edge[i] +
911 g3_edge[i] * g3_edge[i]);
912 }
913
914 Vmath::Reverse(nquad0, &outarray[0], 1, &outarray[0], 1);
915
916 break;
917 case 3:
918 Vmath::Vmul(nqtot, &(df[0][0]), 1, &jac[0], 1,
919 &(tmp_gmat0[0]), 1);
920 Vmath::Vmul(nqtot, &(df[2][0]), 1, &jac[0], 1,
921 &(tmp_gmat2[0]), 1);
922 QuadExp::GetEdgeInterpVals(edge, tmp_gmat0, g0_edge);
923 QuadExp::GetEdgeInterpVals(edge, tmp_gmat2, g2_edge);
924
925 for (i = 0; i < nquad1; ++i)
926 {
927 outarray[i] = sqrt(g0_edge[i] * g0_edge[i] +
928 g2_edge[i] * g2_edge[i]);
929 }
930
931 Vmath::Reverse(nquad1, &outarray[0], 1, &outarray[0], 1);
932
933 break;
934 default:
935 ASSERTL0(false, "edge value (< 3) is out of range");
936 break;
937 }
938 }
939 }
940 else
941 {
942
943 switch (edge)
944 {
945 case 0:
946
947 for (i = 0; i < nquad0; ++i)
948 {
949 outarray[i] = jac[0] * sqrt(df[1][0] * df[1][0] +
950 df[3][0] * df[3][0]);
951 }
952 break;
953 case 1:
954 for (i = 0; i < nquad1; ++i)
955 {
956 outarray[i] = jac[0] * sqrt(df[0][0] * df[0][0] +
957 df[2][0] * df[2][0]);
958 }
959 break;
960 case 2:
961 for (i = 0; i < nquad0; ++i)
962 {
963 outarray[i] = jac[0] * sqrt(df[1][0] * df[1][0] +
964 df[3][0] * df[3][0]);
965 }
966 break;
967 case 3:
968 for (i = 0; i < nquad1; ++i)
969 {
970 outarray[i] = jac[0] * sqrt(df[0][0] * df[0][0] +
971 df[2][0] * df[2][0]);
972 }
973 break;
974 default:
975 ASSERTL0(false, "edge value (< 3) is out of range");
976 break;
977 }
978 }
979}
980
982{
983 int i;
984 const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
985 GetGeom()->GetMetricInfo();
986 SpatialDomains::GeomType type = geomFactors->GetGtype();
987
989 for (i = 0; i < ptsKeys.size(); ++i)
990 {
991 // Need at least 2 points for computing normals
992 if (ptsKeys[i].GetNumPoints() == 1)
993 {
994 LibUtilities::PointsKey pKey(2, ptsKeys[i].GetPointsType());
995 ptsKeys[i] = pKey;
996 }
997 }
998
1000 geomFactors->GetDerivFactors(ptsKeys);
1001 const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
1002 int nqe;
1003 if (edge == 0 || edge == 2)
1004 {
1005 nqe = m_base[0]->GetNumPoints();
1006 }
1007 else
1008 {
1009 nqe = m_base[1]->GetNumPoints();
1010 }
1011 int vCoordDim = GetCoordim();
1012
1015 for (i = 0; i < vCoordDim; ++i)
1016 {
1017 normal[i] = Array<OneD, NekDouble>(nqe);
1018 }
1019
1020 size_t nqb = nqe;
1021 size_t nbnd = edge;
1024
1025 // Regular geometry case
1026 if ((type == SpatialDomains::eRegular) ||
1028 {
1029 NekDouble fac;
1030 // Set up normals
1031 switch (edge)
1032 {
1033 case 0:
1034 for (i = 0; i < vCoordDim; ++i)
1035 {
1036 Vmath::Fill(nqe, -df[2 * i + 1][0], normal[i], 1);
1037 }
1038 break;
1039 case 1:
1040 for (i = 0; i < vCoordDim; ++i)
1041 {
1042 Vmath::Fill(nqe, df[2 * i][0], normal[i], 1);
1043 }
1044 break;
1045 case 2:
1046 for (i = 0; i < vCoordDim; ++i)
1047 {
1048 Vmath::Fill(nqe, df[2 * i + 1][0], normal[i], 1);
1049 }
1050 break;
1051 case 3:
1052 for (i = 0; i < vCoordDim; ++i)
1053 {
1054 Vmath::Fill(nqe, -df[2 * i][0], normal[i], 1);
1055 }
1056 break;
1057 default:
1058 ASSERTL0(false, "edge is out of range (edge < 4)");
1059 }
1060
1061 // normalise
1062 fac = 0.0;
1063 for (i = 0; i < vCoordDim; ++i)
1064 {
1065 fac += normal[i][0] * normal[i][0];
1066 }
1067 fac = 1.0 / sqrt(fac);
1068
1069 Vmath::Fill(nqb, fac, length, 1);
1070
1071 for (i = 0; i < vCoordDim; ++i)
1072 {
1073 Vmath::Smul(nqe, fac, normal[i], 1, normal[i], 1);
1074 }
1075 }
1076 else // Set up deformed normals
1077 {
1078 int j;
1079
1080 int nquad0 = ptsKeys[0].GetNumPoints();
1081 int nquad1 = ptsKeys[1].GetNumPoints();
1082
1083 LibUtilities::PointsKey from_key;
1084
1085 Array<OneD, NekDouble> normals(vCoordDim * max(nquad0, nquad1), 0.0);
1086 Array<OneD, NekDouble> edgejac(vCoordDim * max(nquad0, nquad1), 0.0);
1087
1088 // Extract Jacobian along edges and recover local
1089 // derivates (dx/dr) for polynomial interpolation by
1090 // multiplying m_gmat by jacobian
1091
1092 // Implementation for all the basis except Gauss points
1095 {
1096 switch (edge)
1097 {
1098 case 0:
1099 for (j = 0; j < nquad0; ++j)
1100 {
1101 edgejac[j] = jac[j];
1102 for (i = 0; i < vCoordDim; ++i)
1103 {
1104 normals[i * nquad0 + j] =
1105 -df[2 * i + 1][j] * edgejac[j];
1106 }
1107 }
1108 from_key = ptsKeys[0];
1109 break;
1110 case 1:
1111 for (j = 0; j < nquad1; ++j)
1112 {
1113 edgejac[j] = jac[nquad0 * j + nquad0 - 1];
1114 for (i = 0; i < vCoordDim; ++i)
1115 {
1116 normals[i * nquad1 + j] =
1117 df[2 * i][nquad0 * j + nquad0 - 1] * edgejac[j];
1118 }
1119 }
1120 from_key = ptsKeys[1];
1121 break;
1122 case 2:
1123 for (j = 0; j < nquad0; ++j)
1124 {
1125 edgejac[j] = jac[nquad0 * (nquad1 - 1) + j];
1126 for (i = 0; i < vCoordDim; ++i)
1127 {
1128 normals[i * nquad0 + j] =
1129 (df[2 * i + 1][nquad0 * (nquad1 - 1) + j]) *
1130 edgejac[j];
1131 }
1132 }
1133 from_key = ptsKeys[0];
1134 break;
1135 case 3:
1136 for (j = 0; j < nquad1; ++j)
1137 {
1138 edgejac[j] = jac[nquad0 * j];
1139 for (i = 0; i < vCoordDim; ++i)
1140 {
1141 normals[i * nquad1 + j] =
1142 -df[2 * i][nquad0 * j] * edgejac[j];
1143 }
1144 }
1145 from_key = ptsKeys[1];
1146 break;
1147 default:
1148 ASSERTL0(false, "edge is out of range (edge < 3)");
1149 }
1150 }
1151 else
1152 {
1153 int nqtot = nquad0 * nquad1;
1154 Array<OneD, NekDouble> tmp_gmat(nqtot, 0.0);
1155 Array<OneD, NekDouble> tmp_gmat_edge(nqe, 0.0);
1156
1157 switch (edge)
1158 {
1159 case 0:
1160 for (j = 0; j < nquad0; ++j)
1161 {
1162 for (i = 0; i < vCoordDim; ++i)
1163 {
1164 Vmath::Vmul(nqtot, &(df[2 * i + 1][0]), 1, &jac[0],
1165 1, &(tmp_gmat[0]), 1);
1166 QuadExp::GetEdgeInterpVals(edge, tmp_gmat,
1167 tmp_gmat_edge);
1168 normals[i * nquad0 + j] = -tmp_gmat_edge[j];
1169 }
1170 }
1171 from_key = ptsKeys[0];
1172 break;
1173 case 1:
1174 for (j = 0; j < nquad1; ++j)
1175 {
1176 for (i = 0; i < vCoordDim; ++i)
1177 {
1178 Vmath::Vmul(nqtot, &(df[2 * i][0]), 1, &jac[0], 1,
1179 &(tmp_gmat[0]), 1);
1180 QuadExp::GetEdgeInterpVals(edge, tmp_gmat,
1181 tmp_gmat_edge);
1182 normals[i * nquad1 + j] = tmp_gmat_edge[j];
1183 }
1184 }
1185 from_key = ptsKeys[1];
1186 break;
1187 case 2:
1188 for (j = 0; j < nquad0; ++j)
1189 {
1190 for (i = 0; i < vCoordDim; ++i)
1191 {
1192 Vmath::Vmul(nqtot, &(df[2 * i + 1][0]), 1, &jac[0],
1193 1, &(tmp_gmat[0]), 1);
1194 QuadExp::GetEdgeInterpVals(edge, tmp_gmat,
1195 tmp_gmat_edge);
1196 normals[i * nquad0 + j] = tmp_gmat_edge[j];
1197 }
1198 }
1199 from_key = ptsKeys[0];
1200 break;
1201 case 3:
1202 for (j = 0; j < nquad1; ++j)
1203 {
1204 for (i = 0; i < vCoordDim; ++i)
1205 {
1206 Vmath::Vmul(nqtot, &(df[2 * i][0]), 1, &jac[0], 1,
1207 &(tmp_gmat[0]), 1);
1208 QuadExp::GetEdgeInterpVals(edge, tmp_gmat,
1209 tmp_gmat_edge);
1210 normals[i * nquad1 + j] = -tmp_gmat_edge[j];
1211 }
1212 }
1213 from_key = ptsKeys[1];
1214 break;
1215 default:
1216 ASSERTL0(false, "edge is out of range (edge < 3)");
1217 }
1218 }
1219
1220 int nq = from_key.GetNumPoints();
1221 Array<OneD, NekDouble> work(nqe, 0.0);
1222
1223 // interpolate Jacobian and invert
1224 LibUtilities::Interp1D(from_key, jac, m_base[0]->GetPointsKey(), work);
1225 Vmath::Sdiv(nqe, 1.0, &work[0], 1, &work[0], 1);
1226
1227 // interpolate
1228 for (i = 0; i < GetCoordim(); ++i)
1229 {
1230 LibUtilities::Interp1D(from_key, &normals[i * nq],
1231 m_base[0]->GetPointsKey(), &normal[i][0]);
1232 Vmath::Vmul(nqe, work, 1, normal[i], 1, normal[i], 1);
1233 }
1234
1235 // normalise normal vectors
1236 Vmath::Zero(nqe, work, 1);
1237 for (i = 0; i < GetCoordim(); ++i)
1238 {
1239 Vmath::Vvtvp(nqe, normal[i], 1, normal[i], 1, work, 1, work, 1);
1240 }
1241
1242 Vmath::Vsqrt(nqe, work, 1, work, 1);
1243 Vmath::Sdiv(nqe, 1.0, work, 1, work, 1);
1244
1245 Vmath::Vcopy(nqb, work, 1, length, 1);
1246
1247 for (i = 0; i < GetCoordim(); ++i)
1248 {
1249 Vmath::Vmul(nqe, normal[i], 1, work, 1, normal[i], 1);
1250 }
1251 }
1252 if (GetGeom()->GetEorient(edge) == StdRegions::eBackwards)
1253 {
1254 for (i = 0; i < vCoordDim; ++i)
1255 {
1256 if (geomFactors->GetGtype() == SpatialDomains::eDeformed)
1257 {
1258 Vmath::Reverse(nqe, normal[i], 1, normal[i], 1);
1259 }
1260 }
1261 }
1262}
1263
1265 const NekDouble *data, const std::vector<unsigned int> &nummodes,
1266 int mode_offset, NekDouble *coeffs,
1267 std::vector<LibUtilities::BasisType> &fromType)
1268{
1269 int data_order0 = nummodes[mode_offset];
1270 int fillorder0 = std::min(m_base[0]->GetNumModes(), data_order0);
1271
1272 int data_order1 = nummodes[mode_offset + 1];
1273 int order1 = m_base[1]->GetNumModes();
1274 int fillorder1 = min(order1, data_order1);
1275
1276 // Check if same basis
1277 if (fromType[0] != m_base[0]->GetBasisType() ||
1278 fromType[1] != m_base[1]->GetBasisType())
1279 {
1280 // Construct a quad with the appropriate basis type at our
1281 // quadrature points, and one more to do a forwards
1282 // transform. We can then copy the output to coeffs.
1283 StdRegions::StdQuadExp tmpQuad(
1284 LibUtilities::BasisKey(fromType[0], data_order0,
1285 m_base[0]->GetPointsKey()),
1286 LibUtilities::BasisKey(fromType[1], data_order1,
1287 m_base[1]->GetPointsKey()));
1288 StdRegions::StdQuadExp tmpQuad2(m_base[0]->GetBasisKey(),
1289 m_base[1]->GetBasisKey());
1290
1291 Array<OneD, const NekDouble> tmpData(tmpQuad.GetNcoeffs(), data);
1292 Array<OneD, NekDouble> tmpBwd(tmpQuad2.GetTotPoints());
1293 Array<OneD, NekDouble> tmpOut(tmpQuad2.GetNcoeffs());
1294
1295 tmpQuad.BwdTrans(tmpData, tmpBwd);
1296 tmpQuad2.FwdTrans(tmpBwd, tmpOut);
1297 Vmath::Vcopy(tmpOut.size(), &tmpOut[0], 1, coeffs, 1);
1298
1299 return;
1300 }
1301
1302 switch (m_base[0]->GetBasisType())
1303 {
1305 {
1306 int i;
1307 int cnt = 0;
1308 int cnt1 = 0;
1309
1311 "Extraction routine not set up for this basis");
1312
1313 Vmath::Zero(m_ncoeffs, coeffs, 1);
1314 for (i = 0; i < fillorder0; ++i)
1315 {
1316 Vmath::Vcopy(fillorder1, data + cnt, 1, coeffs + cnt1, 1);
1317 cnt += data_order1;
1318 cnt1 += order1;
1319 }
1320 }
1321 break;
1323 {
1324 LibUtilities::PointsKey p0(nummodes[0],
1326 LibUtilities::PointsKey p1(nummodes[1],
1328 LibUtilities::PointsKey t0(m_base[0]->GetNumModes(),
1330 LibUtilities::PointsKey t1(m_base[1]->GetNumModes(),
1332 LibUtilities::Interp2D(p0, p1, data, t0, t1, coeffs);
1333 }
1334 break;
1336 {
1337 // Assume that input is also Gll_Lagrange but no way to check;
1338 LibUtilities::PointsKey p0(nummodes[0],
1340 LibUtilities::PointsKey p1(nummodes[1],
1342 LibUtilities::PointsKey t0(m_base[0]->GetNumModes(),
1344 LibUtilities::PointsKey t1(m_base[1]->GetNumModes(),
1346 LibUtilities::Interp2D(p0, p1, data, t0, t1, coeffs);
1347 }
1348 break;
1349 default:
1350 ASSERTL0(false, "basis is either not set up or not hierarchicial");
1351 }
1352}
1353
1355{
1356 return m_geom->GetEorient(edge);
1357}
1358
1360{
1361 DNekMatSharedPtr returnval;
1362 switch (mkey.GetMatrixType())
1363 {
1371 returnval = Expansion2D::v_GenMatrix(mkey);
1372 break;
1373 default:
1374 returnval = StdQuadExp::v_GenMatrix(mkey);
1375 }
1376 return returnval;
1377}
1378
1380 const StdRegions::StdMatrixKey &mkey)
1381{
1382 LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
1383 LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
1386 return tmp->GetStdMatrix(mkey);
1387}
1388
1390{
1391 return m_matrixManager[mkey];
1392}
1393
1395{
1396 m_matrixManager.DeleteObject(mkey);
1397}
1398
1400{
1401 return m_staticCondMatrixManager[mkey];
1402}
1403
1405{
1406 m_staticCondMatrixManager.DeleteObject(mkey);
1407}
1408
1410 Array<OneD, NekDouble> &outarray,
1411 const StdRegions::StdMatrixKey &mkey)
1412{
1413 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1414}
1415
1417 Array<OneD, NekDouble> &outarray,
1418 const StdRegions::StdMatrixKey &mkey)
1419{
1420 QuadExp::LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1421}
1422
1423void QuadExp::v_LaplacianMatrixOp(const int k1, const int k2,
1424 const Array<OneD, const NekDouble> &inarray,
1425 Array<OneD, NekDouble> &outarray,
1426 const StdRegions::StdMatrixKey &mkey)
1427{
1428 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1429}
1430
1432 const Array<OneD, const NekDouble> &inarray,
1433 Array<OneD, NekDouble> &outarray,
1434 const StdRegions::StdMatrixKey &mkey)
1435{
1436 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1437}
1438
1440 const Array<OneD, const NekDouble> &inarray,
1441 Array<OneD, NekDouble> &outarray, const StdRegions::StdMatrixKey &mkey)
1442{
1443 StdExpansion::WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
1444}
1445
1447 const Array<OneD, const NekDouble> &inarray,
1448 Array<OneD, NekDouble> &outarray, const StdRegions::StdMatrixKey &mkey)
1449{
1450 StdExpansion::MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
1451}
1452
1454 Array<OneD, NekDouble> &outarray,
1455 const StdRegions::StdMatrixKey &mkey)
1456{
1457 QuadExp::HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1458}
1459
1461 const Array<OneD, const NekDouble> &inarray,
1462 Array<OneD, NekDouble> &outarray)
1463{
1464 int n_coeffs = inarray.size();
1465
1466 Array<OneD, NekDouble> coeff(n_coeffs);
1467 Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
1468 Array<OneD, NekDouble> tmp, tmp2;
1469
1470 int nmodes0 = m_base[0]->GetNumModes();
1471 int nmodes1 = m_base[1]->GetNumModes();
1472 int numMax = nmodes0;
1473
1474 Vmath::Vcopy(n_coeffs, inarray, 1, coeff_tmp, 1);
1475
1476 const LibUtilities::PointsKey Pkey0(nmodes0,
1478 const LibUtilities::PointsKey Pkey1(nmodes1,
1480 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
1481 LibUtilities::BasisKey b1(m_base[1]->GetBasisType(), nmodes1, Pkey1);
1482 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
1483 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A, nmodes1, Pkey1);
1484
1485 LibUtilities::InterpCoeff2D(b0, b1, coeff_tmp, bortho0, bortho1, coeff);
1486
1487 Vmath::Zero(n_coeffs, coeff_tmp, 1);
1488
1489 int cnt = 0;
1490 for (int i = 0; i < numMin + 1; ++i)
1491 {
1492 Vmath::Vcopy(numMin, tmp = coeff + cnt, 1, tmp2 = coeff_tmp + cnt, 1);
1493
1494 cnt = i * numMax;
1495 }
1496
1497 LibUtilities::InterpCoeff2D(bortho0, bortho1, coeff_tmp, b0, b1, outarray);
1498}
1499
1501 const Array<OneD, const NekDouble> &inarray,
1503{
1504 if (m_metrics.count(eMetricLaplacian00) == 0)
1505 {
1507 }
1508
1509 int nquad0 = m_base[0]->GetNumPoints();
1510 int nquad1 = m_base[1]->GetNumPoints();
1511 int nqtot = nquad0 * nquad1;
1512 int nmodes0 = m_base[0]->GetNumModes();
1513 int nmodes1 = m_base[1]->GetNumModes();
1514 int wspsize =
1515 max(max(max(nqtot, m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
1516
1517 ASSERTL1(wsp.size() >= 3 * wspsize, "Workspace is of insufficient size.");
1518
1519 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
1520 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
1521 const Array<OneD, const NekDouble> &dbase0 = m_base[0]->GetDbdata();
1522 const Array<OneD, const NekDouble> &dbase1 = m_base[1]->GetDbdata();
1523 const Array<OneD, const NekDouble> &metric00 =
1525 const Array<OneD, const NekDouble> &metric01 =
1527 const Array<OneD, const NekDouble> &metric11 =
1529
1530 // Allocate temporary storage
1531 Array<OneD, NekDouble> wsp0(wsp);
1532 Array<OneD, NekDouble> wsp1(wsp + wspsize);
1533 Array<OneD, NekDouble> wsp2(wsp + 2 * wspsize);
1534
1535 StdExpansion2D::PhysTensorDeriv(inarray, wsp1, wsp2);
1536
1537 // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1538 // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1539 // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1540 // especially for this purpose
1541 Vmath::Vvtvvtp(nqtot, &metric00[0], 1, &wsp1[0], 1, &metric01[0], 1,
1542 &wsp2[0], 1, &wsp0[0], 1);
1543 Vmath::Vvtvvtp(nqtot, &metric01[0], 1, &wsp1[0], 1, &metric11[0], 1,
1544 &wsp2[0], 1, &wsp2[0], 1);
1545
1546 // outarray = m = (D_xi1 * B)^T * k
1547 // wsp1 = n = (D_xi2 * B)^T * l
1548 IProductWRTBase_SumFacKernel(dbase0, base1, wsp0, outarray, wsp1, false,
1549 true);
1550 IProductWRTBase_SumFacKernel(base0, dbase1, wsp2, wsp1, wsp0, true, false);
1551
1552 // outarray = outarray + wsp1
1553 // = L * u_hat
1554 Vmath::Vadd(m_ncoeffs, wsp1.get(), 1, outarray.get(), 1, outarray.get(), 1);
1555}
1556
1558{
1559 if (m_metrics.count(eMetricQuadrature) == 0)
1560 {
1562 }
1563
1564 const SpatialDomains::GeomType type = m_metricinfo->GetGtype();
1565 const unsigned int nqtot = GetTotPoints();
1566 const unsigned int dim = 2;
1567 const MetricType m[3][3] = {
1571
1572 const Array<TwoD, const NekDouble> gmat =
1573 m_metricinfo->GetGmat(GetPointsKeys());
1574 for (unsigned int i = 0; i < dim; ++i)
1575 {
1576 for (unsigned int j = i; j < dim; ++j)
1577 {
1578 m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1579 if (type == SpatialDomains::eDeformed)
1580 {
1581 Vmath::Vcopy(nqtot, &gmat[i * dim + j][0], 1,
1582 &m_metrics[m[i][j]][0], 1);
1583 }
1584 else
1585 {
1586 Vmath::Fill(nqtot, gmat[i * dim + j][0], &m_metrics[m[i][j]][0],
1587 1);
1588 }
1590 }
1591 }
1592}
1593
1595 const StdRegions::StdMatrixKey &mkey)
1596{
1597 int nq = GetTotPoints();
1598
1599 // Calculate sqrt of the Jacobian
1601 Array<OneD, NekDouble> sqrt_jac(nq);
1602 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1603 {
1604 Vmath::Vsqrt(nq, jac, 1, sqrt_jac, 1);
1605 }
1606 else
1607 {
1608 Vmath::Fill(nq, sqrt(jac[0]), sqrt_jac, 1);
1609 }
1610
1611 // Multiply array by sqrt(Jac)
1612 Vmath::Vmul(nq, sqrt_jac, 1, array, 1, array, 1);
1613
1614 // Apply std region filter
1615 StdQuadExp::v_SVVLaplacianFilter(array, mkey);
1616
1617 // Divide by sqrt(Jac)
1618 Vmath::Vdiv(nq, array, 1, sqrt_jac, 1, array, 1);
1619}
1620
1621/** @brief: This method gets all of the factors which are
1622 required as part of the Gradient Jump Penalty (GJP)
1623 stabilisation and involves the product of the normal and
1624 geometric factors along the element trace.
1625*/
1627 Array<OneD, Array<OneD, NekDouble>> &d0factors,
1628 Array<OneD, Array<OneD, NekDouble>> &d1factors,
1629 Array<OneD, Array<OneD, NekDouble>> &d2factors)
1630{
1631 boost::ignore_unused(d2factors); // for 3D shapes
1632 int nquad0 = GetNumPoints(0);
1633 int nquad1 = GetNumPoints(1);
1634
1636 m_metricinfo->GetDerivFactors(GetPointsKeys());
1637
1638 if (d0factors.size() != 4)
1639 {
1640 d0factors = Array<OneD, Array<OneD, NekDouble>>(4);
1641 d1factors = Array<OneD, Array<OneD, NekDouble>>(4);
1642 }
1643
1644 if (d0factors[0].size() != nquad0)
1645 {
1646 d0factors[0] = Array<OneD, NekDouble>(nquad0);
1647 d0factors[2] = Array<OneD, NekDouble>(nquad0);
1648 d1factors[0] = Array<OneD, NekDouble>(nquad0);
1649 d1factors[2] = Array<OneD, NekDouble>(nquad0);
1650 }
1651
1652 if (d0factors[1].size() != nquad1)
1653 {
1654 d0factors[1] = Array<OneD, NekDouble>(nquad1);
1655 d0factors[3] = Array<OneD, NekDouble>(nquad1);
1656 d1factors[1] = Array<OneD, NekDouble>(nquad1);
1657 d1factors[3] = Array<OneD, NekDouble>(nquad1);
1658 }
1659
1660 // Outwards normals
1662 GetTraceNormal(0);
1664 GetTraceNormal(1);
1666 GetTraceNormal(2);
1668 GetTraceNormal(3);
1669
1670 int ncoords = normal_0.size();
1671
1672 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1673 {
1674 // needs checking for 3D coords
1675
1676 // factors 0 and 2
1677 for (int i = 0; i < nquad0; ++i)
1678 {
1679 d0factors[0][i] = df[0][i] * normal_0[0][i];
1680 d0factors[2][i] = df[0][nquad0 * (nquad1 - 1) + i] * normal_2[0][i];
1681
1682 d1factors[0][i] = df[1][i] * normal_0[0][i];
1683 d1factors[2][i] = df[1][nquad0 * (nquad1 - 1) + i] * normal_2[0][i];
1684 }
1685
1686 for (int n = 1; n < ncoords; ++n)
1687 {
1688 // d xi_1/dy n_y
1689 // needs checking for 3D coords
1690 for (int i = 0; i < nquad0; ++i)
1691 {
1692 d0factors[0][i] += df[2 * n][i] * normal_0[n][i];
1693 d0factors[2][i] +=
1694 df[2 * n][nquad0 * (nquad1 - 1) + i] * normal_2[n][i];
1695
1696 d1factors[0][i] += df[2 * n + 1][i] * normal_0[n][i];
1697 d1factors[2][i] +=
1698 df[2 * n + 1][nquad0 * (nquad1 - 1) + i] * normal_2[n][i];
1699 }
1700 }
1701
1702 // faces 1 and 3
1703 for (int i = 0; i < nquad1; ++i)
1704 {
1705 d0factors[1][i] = df[0][(i + 1) * nquad0 - 1] * normal_1[0][i];
1706 d0factors[3][i] = df[0][i * nquad0] * normal_3[0][i];
1707
1708 d1factors[1][i] = df[1][(i + 1) * nquad0 - 1] * normal_1[0][i];
1709 d1factors[3][i] = df[1][i * nquad0] * normal_3[0][i];
1710 }
1711
1712 for (int n = 1; n < ncoords; ++n)
1713 {
1714 for (int i = 0; i < nquad1; ++i)
1715 {
1716 d0factors[1][i] +=
1717 df[2 * n][(i + 1) * nquad0 - 1] * normal_1[n][i];
1718 d0factors[3][i] += df[2 * n][i * nquad0] * normal_3[n][i];
1719
1720 d1factors[1][i] +=
1721 df[2 * n + 1][(i + 1) * nquad0 - 1] * normal_1[n][i];
1722 d1factors[3][i] += df[2 * n + 1][i * nquad0] * normal_3[n][i];
1723 }
1724 }
1725 }
1726 else
1727 {
1728 // d xi_2/dx n_x
1729 for (int i = 0; i < nquad0; ++i)
1730 {
1731 d1factors[0][i] = df[1][0] * normal_0[0][i];
1732 d1factors[2][i] = df[1][0] * normal_2[0][i];
1733 }
1734
1735 // d xi_1/dx n_x
1736 for (int i = 0; i < nquad1; ++i)
1737 {
1738 d0factors[1][i] = df[0][0] * normal_1[0][i];
1739 d0factors[3][i] = df[0][0] * normal_3[0][i];
1740 }
1741
1742 for (int n = 1; n < ncoords; ++n)
1743 {
1744 // d xi_2/dy n_y
1745 // needs checking for 3D coords
1746 for (int i = 0; i < nquad0; ++i)
1747 {
1748 d1factors[0][i] += df[2 * n + 1][0] * normal_0[n][i];
1749 d1factors[2][i] += df[2 * n + 1][0] * normal_2[n][i];
1750 }
1751
1752 // d xi_1/dy n_y
1753 // needs checking for 3D coords
1754 for (int i = 0; i < nquad1; ++i)
1755 {
1756 d0factors[1][i] += df[2 * n][0] * normal_1[n][i];
1757 d0factors[3][i] += df[2 * n][0] * normal_3[n][i];
1758 }
1759 }
1760
1761 // d1factors
1762 // d xi_1/dx n_x
1763 for (int i = 0; i < nquad0; ++i)
1764 {
1765 d0factors[0][i] = df[0][0] * normal_0[0][i];
1766 d0factors[2][i] = df[0][0] * normal_2[0][i];
1767 }
1768
1769 // d xi_2/dx n_x
1770 for (int i = 0; i < nquad1; ++i)
1771 {
1772 d1factors[1][i] = df[1][0] * normal_1[0][i];
1773 d1factors[3][i] = df[1][0] * normal_3[0][i];
1774 }
1775
1776 for (int n = 1; n < ncoords; ++n)
1777 {
1778 // d xi_1/dy n_y
1779 // needs checking for 3D coords
1780 for (int i = 0; i < nquad0; ++i)
1781 {
1782 d0factors[0][i] += df[2 * n][0] * normal_0[n][i];
1783 d0factors[2][i] += df[2 * n][0] * normal_2[n][i];
1784 }
1785
1786 // d xi_2/dy n_y
1787 // needs checking for 3D coords
1788 for (int i = 0; i < nquad1; ++i)
1789 {
1790 d1factors[1][i] += df[2 * n + 1][0] * normal_1[n][i];
1791 d1factors[3][i] += df[2 * n + 1][0] * normal_3[n][i];
1792 }
1793 }
1794 }
1795}
1796} // namespace LocalRegions
1797} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:49
Describes the specification for a Basis.
Definition: Basis.h:47
Defines a specification for a set of points.
Definition: Points.h:55
size_t GetNumPoints() const
Definition: Points.h:90
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
SpatialDomains::Geometry2DSharedPtr GetGeom2D() const
Definition: Expansion2D.h:172
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
ExpansionSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion.h:443
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: Expansion.cpp:535
int GetLeftAdjacentElementTrace() const
Definition: Expansion.h:456
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:276
StdRegions::Orientation GetTraceOrient(int trace)
Definition: Expansion.h:170
const NormalVector & GetTraceNormal(const int id)
Definition: Expansion.cpp:255
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: QuadExp.cpp:1409
virtual StdRegions::Orientation v_GetTraceOrient(int edge) override
Definition: QuadExp.cpp:1354
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const MatrixKey &mkey) override
Definition: QuadExp.cpp:1399
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Transform a given function from physical quadrature space to coefficient space.
Definition: QuadExp.cpp:231
virtual DNekScalMatSharedPtr v_GetLocMatrix(const MatrixKey &mkey) override
Definition: QuadExp.cpp:1389
virtual StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const override
Definition: QuadExp.cpp:531
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: QuadExp.cpp:1453
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
Definition: QuadExp.cpp:549
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculate the inner product of inarray with respect to the basis B=base0*base1 and put into outarray.
Definition: QuadExp.cpp:374
virtual void v_ComputeLaplacianMetric() override
Definition: QuadExp.cpp:1557
virtual void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp) override
Definition: QuadExp.cpp:1500
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: QuadExp.cpp:394
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
Definition: QuadExp.cpp:456
void GetEdgeInterpVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: QuadExp.cpp:661
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: QuadExp.cpp:422
virtual void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray) override
Definition: QuadExp.cpp:490
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: QuadExp.cpp:1416
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals) override
This function evaluates the expansion at a single (arbitrary) point of the domain.
Definition: QuadExp.cpp:578
virtual void v_GetTracePhysVals(const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient) override
Definition: QuadExp.cpp:599
virtual void v_WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: QuadExp.cpp:1439
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=NullNekDouble1DArray) override
Calculate the derivative of the physical points.
Definition: QuadExp.cpp:96
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrates the specified function over the domain.
Definition: QuadExp.cpp:73
void v_DropLocMatrix(const MatrixKey &mkey) override
Definition: QuadExp.cpp:1394
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: QuadExp.h:247
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: QuadExp.cpp:1264
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: QuadExp.cpp:542
virtual void v_ComputeTraceNormal(const int edge) override
Definition: QuadExp.cpp:981
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: QuadExp.cpp:185
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: QuadExp.h:249
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const override
Definition: QuadExp.cpp:525
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
Definition: QuadExp.cpp:1359
void v_DropLocStaticCondMatrix(const MatrixKey &mkey) override
Definition: QuadExp.cpp:1404
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: QuadExp.cpp:387
virtual void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: QuadExp.cpp:1446
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdRegions::StdMatrixKey &mkey) override
Definition: QuadExp.cpp:1379
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey) override
Definition: QuadExp.cpp:1594
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 (GJP) s...
Definition: QuadExp.cpp:1626
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
Definition: QuadExp.cpp:570
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: QuadExp.cpp:254
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: QuadExp.cpp:1460
virtual void v_GetTracePhysMap(const int edge, Array< OneD, int > &outarray) override
Definition: QuadExp.cpp:725
virtual void v_GetTraceQFactors(const int edge, Array< OneD, NekDouble > &outarray) override
Definition: QuadExp.cpp:767
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: QuadExp.cpp:1431
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 > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
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
void HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
const LibUtilities::PointsKeyVector GetPointsKeys() const
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:758
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:619
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
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
Definition: StdExpansion.h:690
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:680
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 IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:729
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:87
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
Definition: Blas.hpp:213
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:165
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 Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:49
void InterpCoeff2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:69
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
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:48
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:59
@ 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< SegExp > SegExpSharedPtr
Definition: SegExp.h:251
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< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:65
std::shared_ptr< StdQuadExp > StdQuadExpSharedPtr
Definition: StdQuadExp.h:243
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
StdRegions::ConstFactorMap factors
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 Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
svtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:746
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:207
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:617
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:569
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:354
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:245
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
Definition: Vmath.cpp:319
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.cpp:280
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:43
void Reverse(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1222
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
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:414
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294