Nektar++
NodalTriExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NodalTriExp.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: NodalTriExp routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
37
38using namespace std;
39
41{
43 const LibUtilities::BasisKey &Bb,
44 const LibUtilities::PointsType Ntype,
46 : StdExpansion(LibUtilities::StdTriData::getNumberOfCoefficients(
47 Ba.GetNumModes(), (Bb.GetNumModes())),
48 2, Ba, Bb),
49 StdExpansion2D(LibUtilities::StdTriData::getNumberOfCoefficients(
50 Ba.GetNumModes(), (Bb.GetNumModes())),
51 Ba, Bb),
52 StdNodalTriExp(Ba, Bb, Ntype), Expansion(geom), Expansion2D(geom),
53 m_matrixManager(
54 std::bind(&Expansion2D::CreateMatrix, this, std::placeholders::_1),
55 std::string("NodalTriExpMatrix")),
56 m_staticCondMatrixManager(std::bind(&Expansion::CreateStaticCondMatrix,
57 this, std::placeholders::_1),
58 std::string("NodalTriExpStaticCondMatrix"))
59{
60}
61
63 : StdExpansion(T),
64 StdExpansion2D(T), StdRegions::StdTriExp(T), StdRegions::StdNodalTriExp(
65 T),
66 Expansion(T), Expansion2D(T), m_matrixManager(T.m_matrixManager),
67 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
68{
69}
70
71//----------------------------
72// Integration Methods
73//----------------------------
74
75/** \brief Integrate the physical point list \a inarray over region
76 and return the value
77
78 Inputs:\n
79
80 - \a inarray: definition of function to be returned at quadrature point
81 of expansion.
82
83 Outputs:\n
84
85 - returns \f$\int^1_{-1}\int^1_{-1} u(\xi_1, \xi_2) J[i,j] d
86 \xi_1 d \xi_2 \f$ where \f$inarray[i,j] = u(\xi_{1i},\xi_{2j})
87 \f$ and \f$ J[i,j] \f$ is the Jacobian evaluated at the
88 quadrature point.
89*/
90
92{
93 int nquad0 = m_base[0]->GetNumPoints();
94 int nquad1 = m_base[1]->GetNumPoints();
96 NekDouble ival;
97 Array<OneD, NekDouble> tmp(nquad0 * nquad1);
98
99 // multiply inarray with Jacobian
100 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
101 {
102 Vmath::Vmul(nquad0 * nquad1, jac, 1, inarray, 1, tmp, 1);
103 }
104 else
105 {
106 Vmath::Smul(nquad0 * nquad1, jac[0], inarray, 1, tmp, 1);
107 }
108
109 // call StdQuadExp version;
110 ival = StdNodalTriExp::v_Integral(tmp);
111 return ival;
112}
113
115 const Array<OneD, const NekDouble> &inarray,
116 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
117{
118 int nquad0 = m_base[0]->GetNumPoints();
119 int nquad1 = m_base[1]->GetNumPoints();
120 int order1 = m_base[1]->GetNumModes();
121
122 if (multiplybyweights)
123 {
124 Array<OneD, NekDouble> tmp(nquad0 * nquad1 + nquad0 * order1);
125 Array<OneD, NekDouble> wsp(tmp + nquad0 * nquad1);
126
127 MultiplyByQuadratureMetric(inarray, tmp);
128 StdTriExp::IProductWRTBase_SumFacKernel(
129 m_base[0]->GetBdata(), m_base[1]->GetBdata(), tmp, outarray, wsp);
130 NodalToModalTranspose(outarray, outarray);
131 }
132 else
133 {
134 Array<OneD, NekDouble> wsp(nquad0 * order1);
135
136 StdTriExp::IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
137 m_base[1]->GetBdata(), inarray,
138 outarray, wsp);
139 NodalToModalTranspose(outarray, outarray);
140 }
141}
142
144 const Array<OneD, const NekDouble> &inarray,
145 Array<OneD, NekDouble> &outarray)
146{
147 int nq = GetTotPoints();
149 DNekScalMatSharedPtr iprodmat = m_matrixManager[iprodmatkey];
150
151 Blas::Dgemv('N', m_ncoeffs, nq, iprodmat->Scale(),
152 (iprodmat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
153 inarray.get(), 1, 0.0, outarray.get(), 1);
154}
155
157 const int dir, const Array<OneD, const NekDouble> &inarray,
158 Array<OneD, NekDouble> &outarray)
159{
160 int nquad0 = m_base[0]->GetNumPoints();
161 int nquad1 = m_base[1]->GetNumPoints();
162 int nqtot = nquad0 * nquad1;
163 int wspsize = max(nqtot, m_ncoeffs);
164
165 Array<OneD, NekDouble> tmp0(4 * wspsize);
166 Array<OneD, NekDouble> tmp1(tmp0 + wspsize);
167 Array<OneD, NekDouble> tmp2(tmp0 + 2 * wspsize);
168 Array<OneD, NekDouble> tmp3(tmp0 + 3 * wspsize);
169
171 tmp2D[0] = tmp1;
172 tmp2D[1] = tmp2;
173
174 AlignVectorToCollapsedDir(dir, inarray, tmp2D);
175
176 MultiplyByQuadratureMetric(tmp1, tmp1);
177 MultiplyByQuadratureMetric(tmp2, tmp2);
178
179 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
180 tmp1, tmp3, tmp0);
181 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
182 tmp2, outarray, tmp0);
183 Vmath::Vadd(m_ncoeffs, tmp3, 1, outarray, 1, outarray, 1);
184
185 NodalToModalTranspose(outarray, outarray);
186}
187
189 const int dir, const Array<OneD, const NekDouble> &inarray,
191{
192 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
193 ASSERTL1((dir == 2) ? (m_geom->GetCoordim() == 3) : true,
194 "Invalid direction.");
195
196 int nquad0 = m_base[0]->GetNumPoints();
197 int nquad1 = m_base[1]->GetNumPoints();
198 int nqtot = nquad0 * nquad1;
199 int wspsize = max(nqtot, m_ncoeffs);
200
202 m_metricinfo->GetDerivFactors(GetPointsKeys());
203
204 Array<OneD, NekDouble> tmp0(4 * wspsize);
205 Array<OneD, NekDouble> tmp3(tmp0 + wspsize);
206 Array<OneD, NekDouble> gfac0(tmp0 + 2 * wspsize);
207 Array<OneD, NekDouble> gfac1(tmp0 + 3 * wspsize);
208
209 Array<OneD, NekDouble> tmp1 = outarray[0];
210 Array<OneD, NekDouble> tmp2 = outarray[1];
211
212 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
213 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
214
215 // set up geometric factor: 2/(1-z1)
216 for (int i = 0; i < nquad1; ++i)
217 {
218 gfac0[i] = 2.0 / (1 - z1[i]);
219 }
220 for (int i = 0; i < nquad0; ++i)
221 {
222 gfac1[i] = 0.5 * (1 + z0[i]);
223 }
224
225 for (int i = 0; i < nquad1; ++i)
226 {
227 Vmath::Smul(nquad0, gfac0[i], &inarray[0] + i * nquad0, 1,
228 &tmp0[0] + i * nquad0, 1);
229 }
230
231 for (int i = 0; i < nquad1; ++i)
232 {
233 Vmath::Vmul(nquad0, &gfac1[0], 1, &tmp0[0] + i * nquad0, 1,
234 &tmp1[0] + i * nquad0, 1);
235 }
236
237 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
238 {
239 Vmath::Vmul(nqtot, &df[2 * dir][0], 1, &tmp0[0], 1, &tmp0[0], 1);
240 Vmath::Vmul(nqtot, &df[2 * dir + 1][0], 1, &tmp1[0], 1, &tmp1[0], 1);
241 Vmath::Vmul(nqtot, &df[2 * dir + 1][0], 1, &inarray[0], 1, &tmp2[0], 1);
242 }
243 else
244 {
245 Vmath::Smul(nqtot, df[2 * dir][0], tmp0, 1, tmp0, 1);
246 Vmath::Smul(nqtot, df[2 * dir + 1][0], tmp1, 1, tmp1, 1);
247 Vmath::Smul(nqtot, df[2 * dir + 1][0], inarray, 1, tmp2, 1);
248 }
249 Vmath::Vadd(nqtot, tmp0, 1, tmp1, 1, tmp1, 1);
250}
251
253 const int dir, const Array<OneD, const NekDouble> &inarray,
254 Array<OneD, NekDouble> &outarray)
255{
256 int nq = GetTotPoints();
258
259 switch (dir)
260 {
261 case 0:
262 {
264 }
265 break;
266 case 1:
267 {
269 }
270 break;
271 case 2:
272 {
274 }
275 break;
276 default:
277 {
278 ASSERTL1(false, "input dir is out of range");
279 }
280 break;
281 }
282
283 MatrixKey iprodmatkey(mtype, DetShapeType(), *this);
284 DNekScalMatSharedPtr iprodmat = m_matrixManager[iprodmatkey];
285
286 Blas::Dgemv('N', m_ncoeffs, nq, iprodmat->Scale(),
287 (iprodmat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
288 inarray.get(), 1, 0.0, outarray.get(), 1);
289}
290
291///////////////////////////////
292/// Differentiation Methods
293///////////////////////////////
294
295/**
296 \brief Calculate the deritive of the physical points
297**/
302{
303 int nquad0 = m_base[0]->GetNumPoints();
304 int nquad1 = m_base[1]->GetNumPoints();
305 int nqtot = nquad0 * nquad1;
307 m_metricinfo->GetDerivFactors(GetPointsKeys());
308
309 Array<OneD, NekDouble> diff0(2 * nqtot);
310 Array<OneD, NekDouble> diff1(diff0 + nqtot);
311
312 StdNodalTriExp::v_PhysDeriv(inarray, diff0, diff1);
313
314 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
315 {
316 if (out_d0.size())
317 {
318 Vmath::Vmul(nqtot, df[0], 1, diff0, 1, out_d0, 1);
319 Vmath::Vvtvp(nqtot, df[1], 1, diff1, 1, out_d0, 1, out_d0, 1);
320 }
321
322 if (out_d1.size())
323 {
324 Vmath::Vmul(nqtot, df[2], 1, diff0, 1, out_d1, 1);
325 Vmath::Vvtvp(nqtot, df[3], 1, diff1, 1, out_d1, 1, out_d1, 1);
326 }
327
328 if (out_d2.size())
329 {
330 Vmath::Vmul(nqtot, df[4], 1, diff0, 1, out_d2, 1);
331 Vmath::Vvtvp(nqtot, df[5], 1, diff1, 1, out_d2, 1, out_d2, 1);
332 }
333 }
334 else // regular geometry
335 {
336 if (out_d0.size())
337 {
338 Vmath::Smul(nqtot, df[0][0], diff0, 1, out_d0, 1);
339 Blas::Daxpy(nqtot, df[1][0], diff1, 1, out_d0, 1);
340 }
341
342 if (out_d1.size())
343 {
344 Vmath::Smul(nqtot, df[2][0], diff0, 1, out_d1, 1);
345 Blas::Daxpy(nqtot, df[3][0], diff1, 1, out_d1, 1);
346 }
347
348 if (out_d2.size())
349 {
350 Vmath::Smul(nqtot, df[4][0], diff0, 1, out_d2, 1);
351 Blas::Daxpy(nqtot, df[5][0], diff1, 1, out_d2, 1);
352 }
353 }
354}
355
356/** \brief Forward transform from physical quadrature space
357 stored in \a inarray and evaluate the expansion coefficients and
358 store in \a (this)->m_coeffs
359
360 Inputs:\n
361
362 - \a inarray: array of physical quadrature points to be transformed
363
364 Outputs:\n
365
366 - (this)->_coeffs: updated array of expansion coefficients.
367
368*/
370 Array<OneD, NekDouble> &outarray)
371{
372 IProductWRTBase(inarray, outarray);
373
374 // get Mass matrix inverse
379 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
380
381 // copy inarray in case inarray == outarray
382 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
384
385 out = (*matsys) * in;
386}
387
389 const Array<OneD, const NekDouble> &inarray,
391{
393
394 if (inarray.get() == outarray.get())
395 {
397 Vmath::Vcopy(m_ncoeffs, inarray.get(), 1, tmp.get(), 1);
398
399 Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, mat->Scale(),
400 (mat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
401 tmp.get(), 1, 0.0, outarray.get(), 1);
402 }
403 else
404 {
405 Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, mat->Scale(),
406 (mat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
407 inarray.get(), 1, 0.0, outarray.get(), 1);
408 }
409}
410
412 Array<OneD, NekDouble> &coords_1,
413 Array<OneD, NekDouble> &coords_2)
414{
415 Expansion::v_GetCoords(coords_0, coords_1, coords_2);
416}
417
418// get the coordinates "coords" at the local coordinates "Lcoords"
421{
422 int i;
423
424 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[1] <= 1.0 && Lcoords[1] >= -1.0 &&
425 Lcoords[1] <= 1.0,
426 "Local coordinates are not in region [-1,1]");
427
428 m_geom->FillGeom();
429
430 for (i = 0; i < m_geom->GetCoordim(); ++i)
431 {
432 coords[i] = m_geom->GetCoord(i, Lcoords);
433 }
434}
435
437 const StdRegions::StdMatrixKey &mkey)
438{
439 LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
440 LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
444
445 return tmp->GetStdMatrix(mkey);
446}
447
449 const Array<OneD, const NekDouble> &coord,
450 const Array<OneD, const NekDouble> &physvals)
451
452{
454
455 ASSERTL0(m_geom, "m_geom not defined");
456 m_geom->GetLocCoords(coord, Lcoord);
457
458 return StdExpansion2D::v_PhysEvaluate(Lcoord, physvals);
459}
460
462 const Array<OneD, NekDouble> &coord,
463 const Array<OneD, const NekDouble> &inarray,
464 std::array<NekDouble, 3> &firstOrderDerivs)
465{
466 Array<OneD, NekDouble> Lcoord(2);
467 ASSERTL0(m_geom, "m_geom not defined");
468 m_geom->GetLocCoords(coord, Lcoord);
469 return StdExpansion2D::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
470}
471
473{
474
476 m_base[0]->GetBasisKey(), m_base[1]->GetBasisKey(),
478}
479
481{
483 m_base[0]->GetPointsKey());
485 m_base[1]->GetPointsKey());
486
488 bkey0, bkey1, m_nodalPointsKey.GetPointsType());
489}
490
492{
493 DNekMatSharedPtr returnval;
494
495 switch (mkey.GetMatrixType())
496 {
503 returnval = Expansion2D::v_GenMatrix(mkey);
504 break;
505 default:
506 returnval = StdNodalTriExp::v_GenMatrix(mkey);
507 break;
508 }
509 return returnval;
510}
511
513{
514 int i;
515 const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
516 GetGeom()->GetMetricInfo();
517 const SpatialDomains::GeomType type = geomFactors->GetGtype();
518
521 geomFactors->GetDerivFactors(ptsKeys);
522 const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
523 int nqe = m_base[0]->GetNumPoints();
524 int dim = GetCoordim();
525
528 for (i = 0; i < dim; ++i)
529 {
530 normal[i] = Array<OneD, NekDouble>(nqe);
531 }
532
533 size_t nqb = nqe;
534 size_t nbnd = edge;
537
538 // Regular geometry case
539 if ((type == SpatialDomains::eRegular) ||
541 {
542 NekDouble fac;
543 // Set up normals
544 switch (edge)
545 {
546 case 0:
547 for (i = 0; i < GetCoordim(); ++i)
548 {
549 Vmath::Fill(nqe, -df[2 * i + 1][0], normal[i], 1);
550 }
551 break;
552 case 1:
553 for (i = 0; i < GetCoordim(); ++i)
554 {
555 Vmath::Fill(nqe, df[2 * i + 1][0] + df[2 * i][0], normal[i],
556 1);
557 }
558 break;
559 case 2:
560 for (i = 0; i < GetCoordim(); ++i)
561 {
562 Vmath::Fill(nqe, -df[2 * i][0], normal[i], 1);
563 }
564 break;
565 default:
566 ASSERTL0(false, "Edge is out of range (edge < 3)");
567 }
568
569 // normalise
570 fac = 0.0;
571 for (i = 0; i < GetCoordim(); ++i)
572 {
573 fac += normal[i][0] * normal[i][0];
574 }
575 fac = 1.0 / sqrt(fac);
576
577 Vmath::Fill(nqb, fac, length, 1);
578
579 for (i = 0; i < GetCoordim(); ++i)
580 {
581 Vmath::Smul(nqe, fac, normal[i], 1, normal[i], 1);
582 }
583 }
584 else // Set up deformed normals
585 {
586 int j;
587
588 int nquad0 = ptsKeys[0].GetNumPoints();
589 int nquad1 = ptsKeys[1].GetNumPoints();
590
592
593 Array<OneD, NekDouble> normals(GetCoordim() * max(nquad0, nquad1), 0.0);
594 Array<OneD, NekDouble> edgejac(GetCoordim() * max(nquad0, nquad1), 0.0);
595
596 // Extract Jacobian along edges and recover local
597 // derivates (dx/dr) for polynomial interpolation by
598 // multiplying m_gmat by jacobian
599 switch (edge)
600 {
601 case 0:
602 for (j = 0; j < nquad0; ++j)
603 {
604 edgejac[j] = jac[j];
605 for (i = 0; i < GetCoordim(); ++i)
606 {
607 normals[i * nquad0 + j] =
608 -df[2 * i + 1][j] * edgejac[j];
609 }
610 }
611 from_key = ptsKeys[0];
612 break;
613 case 1:
614 for (j = 0; j < nquad1; ++j)
615 {
616 edgejac[j] = jac[nquad0 * j + nquad0 - 1];
617 for (i = 0; i < GetCoordim(); ++i)
618 {
619 normals[i * nquad1 + j] =
620 (df[2 * i][nquad0 * j + nquad0 - 1] +
621 df[2 * i + 1][nquad0 * j + nquad0 - 1]) *
622 edgejac[j];
623 }
624 }
625 from_key = ptsKeys[1];
626 break;
627 case 2:
628 for (j = 0; j < nquad1; ++j)
629 {
630 edgejac[j] = jac[nquad0 * j];
631 for (i = 0; i < GetCoordim(); ++i)
632 {
633 normals[i * nquad1 + j] =
634 -df[2 * i][nquad0 * j] * edgejac[j];
635 }
636 }
637 from_key = ptsKeys[1];
638 break;
639 default:
640 ASSERTL0(false, "edge is out of range (edge < 3)");
641 }
642
643 int nq = from_key.GetNumPoints();
644 Array<OneD, NekDouble> work(nqe, 0.0);
645
646 // interpolate Jacobian and invert
647 LibUtilities::Interp1D(from_key, jac, m_base[0]->GetPointsKey(), work);
648 Vmath::Sdiv(nq, 1.0, &work[0], 1, &work[0], 1);
649
650 // interpolate
651 for (i = 0; i < GetCoordim(); ++i)
652 {
653 LibUtilities::Interp1D(from_key, &normals[i * nq],
654 m_base[0]->GetPointsKey(), &normal[i][0]);
655 Vmath::Vmul(nqe, work, 1, normal[i], 1, normal[i], 1);
656 }
657
658 // normalise normal vectors
659 Vmath::Zero(nqe, work, 1);
660 for (i = 0; i < GetCoordim(); ++i)
661 {
662 Vmath::Vvtvp(nqe, normal[i], 1, normal[i], 1, work, 1, work, 1);
663 }
664
665 Vmath::Vsqrt(nqe, work, 1, work, 1);
666 Vmath::Sdiv(nqe, 1.0, work, 1, work, 1);
667
668 Vmath::Vcopy(nqb, work, 1, length, 1);
669
670 for (i = 0; i < GetCoordim(); ++i)
671 {
672 Vmath::Vmul(nqe, normal[i], 1, work, 1, normal[i], 1);
673 }
674
675 // Reverse direction so that points are in
676 // anticlockwise direction if edge >=2
677 if (edge >= 2)
678 {
679 for (i = 0; i < GetCoordim(); ++i)
680 {
681 Vmath::Reverse(nqe, normal[i], 1, normal[i], 1);
682 }
683 }
684 }
685}
686
687} // namespace Nektar::LocalRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Describes the specification for a Basis.
Definition: Basis.h:45
Defines a specification for a set of points.
Definition: Points.h:50
PointsType GetPointsType() const
Definition: Points.h:90
size_t GetNumPoints() const
Definition: Points.h:85
DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
std::map< int, NormalVector > m_traceNormals
Definition: Expansion.h:276
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
Definition: Expansion.h:286
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:167
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:273
void AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Expansion.h:151
void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: Expansion.cpp:530
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:274
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:84
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: NodalTriExp.h:192
void GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords)
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Differentiation Methods.
DNekMatSharedPtr CreateStdMatrix(const StdRegions::StdMatrixKey &mkey)
DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
NekDouble v_PhysEvaluate(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) final
void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
void GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
void GeneralMatrixOp_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return ...
Definition: NodalTriExp.h:82
void IProductWRTDerivBase_MatOp(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals)
StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const override
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
Integrate the physical point list inarray over region.
Definition: NodalTriExp.cpp:91
void IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
void IProductWRTBase_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const override
NodalTriExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::PointsType Ntype, const SpatialDomains::TriGeomSharedPtr &geom)
Constructor using BasisKey class for quadrature points and order definition.
Definition: NodalTriExp.cpp:42
void v_ComputeTraceNormal(const int edge) override
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 GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
const LibUtilities::PointsKeyVector GetPointsKeys() const
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:723
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:83
void NodalToModalTranspose(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
LibUtilities::PointsKey m_nodalPointsKey
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:211
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:135
void 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:47
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:60
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:56
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::shared_ptr< StdNodalTriExp > StdNodalTriExpSharedPtr
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:403
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:347
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.hpp:340
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
Definition: Vmath.hpp:154
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Reverse(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:844
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294