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
40namespace Nektar
41{
42namespace LocalRegions
43{
45 const LibUtilities::BasisKey &Bb,
46 const LibUtilities::PointsType Ntype,
48 : StdExpansion(LibUtilities::StdTriData::getNumberOfCoefficients(
49 Ba.GetNumModes(), (Bb.GetNumModes())),
50 2, Ba, Bb),
51 StdExpansion2D(LibUtilities::StdTriData::getNumberOfCoefficients(
52 Ba.GetNumModes(), (Bb.GetNumModes())),
53 Ba, Bb),
54 StdNodalTriExp(Ba, Bb, Ntype), Expansion(geom), Expansion2D(geom),
55 m_matrixManager(
56 std::bind(&Expansion2D::CreateMatrix, this, std::placeholders::_1),
57 std::string("NodalTriExpMatrix")),
58 m_staticCondMatrixManager(std::bind(&Expansion::CreateStaticCondMatrix,
59 this, std::placeholders::_1),
60 std::string("NodalTriExpStaticCondMatrix"))
61{
62}
63
65 : StdExpansion(T),
66 StdExpansion2D(T), StdRegions::StdTriExp(T), StdRegions::StdNodalTriExp(
67 T),
68 Expansion(T), Expansion2D(T), m_matrixManager(T.m_matrixManager),
69 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
70{
71}
72
73//----------------------------
74// Integration Methods
75//----------------------------
76
77/** \brief Integrate the physical point list \a inarray over region
78 and return the value
79
80 Inputs:\n
81
82 - \a inarray: definition of function to be returned at quadrature point
83 of expansion.
84
85 Outputs:\n
86
87 - returns \f$\int^1_{-1}\int^1_{-1} u(\xi_1, \xi_2) J[i,j] d
88 \xi_1 d \xi_2 \f$ where \f$inarray[i,j] = u(\xi_{1i},\xi_{2j})
89 \f$ and \f$ J[i,j] \f$ is the Jacobian evaluated at the
90 quadrature point.
91*/
92
94{
95 int nquad0 = m_base[0]->GetNumPoints();
96 int nquad1 = m_base[1]->GetNumPoints();
98 NekDouble ival;
99 Array<OneD, NekDouble> tmp(nquad0 * nquad1);
100
101 // multiply inarray with Jacobian
102 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
103 {
104 Vmath::Vmul(nquad0 * nquad1, jac, 1, inarray, 1, tmp, 1);
105 }
106 else
107 {
108 Vmath::Smul(nquad0 * nquad1, jac[0], inarray, 1, tmp, 1);
109 }
110
111 // call StdQuadExp version;
112 ival = StdNodalTriExp::v_Integral(tmp);
113 return ival;
114}
115
117 const Array<OneD, const NekDouble> &inarray,
118 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
119{
120 int nquad0 = m_base[0]->GetNumPoints();
121 int nquad1 = m_base[1]->GetNumPoints();
122 int order1 = m_base[1]->GetNumModes();
123
124 if (multiplybyweights)
125 {
126 Array<OneD, NekDouble> tmp(nquad0 * nquad1 + nquad0 * order1);
127 Array<OneD, NekDouble> wsp(tmp + nquad0 * nquad1);
128
129 MultiplyByQuadratureMetric(inarray, tmp);
130 StdTriExp::IProductWRTBase_SumFacKernel(
131 m_base[0]->GetBdata(), m_base[1]->GetBdata(), tmp, outarray, wsp);
132 NodalToModalTranspose(outarray, outarray);
133 }
134 else
135 {
136 Array<OneD, NekDouble> wsp(nquad0 * order1);
137
138 StdTriExp::IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
139 m_base[1]->GetBdata(), inarray,
140 outarray, wsp);
141 NodalToModalTranspose(outarray, outarray);
142 }
143}
144
146 const Array<OneD, const NekDouble> &inarray,
147 Array<OneD, NekDouble> &outarray)
148{
149 int nq = GetTotPoints();
151 DNekScalMatSharedPtr iprodmat = m_matrixManager[iprodmatkey];
152
153 Blas::Dgemv('N', m_ncoeffs, nq, iprodmat->Scale(),
154 (iprodmat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
155 inarray.get(), 1, 0.0, outarray.get(), 1);
156}
157
159 const int dir, const Array<OneD, const NekDouble> &inarray,
160 Array<OneD, NekDouble> &outarray)
161{
162 int nquad0 = m_base[0]->GetNumPoints();
163 int nquad1 = m_base[1]->GetNumPoints();
164 int nqtot = nquad0 * nquad1;
165 int wspsize = max(nqtot, m_ncoeffs);
166
167 Array<OneD, NekDouble> tmp0(4 * wspsize);
168 Array<OneD, NekDouble> tmp1(tmp0 + wspsize);
169 Array<OneD, NekDouble> tmp2(tmp0 + 2 * wspsize);
170 Array<OneD, NekDouble> tmp3(tmp0 + 3 * wspsize);
171
173 tmp2D[0] = tmp1;
174 tmp2D[1] = tmp2;
175
176 AlignVectorToCollapsedDir(dir, inarray, tmp2D);
177
178 MultiplyByQuadratureMetric(tmp1, tmp1);
179 MultiplyByQuadratureMetric(tmp2, tmp2);
180
181 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
182 tmp1, tmp3, tmp0);
183 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
184 tmp2, outarray, tmp0);
185 Vmath::Vadd(m_ncoeffs, tmp3, 1, outarray, 1, outarray, 1);
186
187 NodalToModalTranspose(outarray, outarray);
188}
189
191 const int dir, const Array<OneD, const NekDouble> &inarray,
193{
194 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
195 ASSERTL1((dir == 2) ? (m_geom->GetCoordim() == 3) : true,
196 "Invalid direction.");
197
198 int nquad0 = m_base[0]->GetNumPoints();
199 int nquad1 = m_base[1]->GetNumPoints();
200 int nqtot = nquad0 * nquad1;
201 int wspsize = max(nqtot, m_ncoeffs);
202
204 m_metricinfo->GetDerivFactors(GetPointsKeys());
205
206 Array<OneD, NekDouble> tmp0(4 * wspsize);
207 Array<OneD, NekDouble> tmp3(tmp0 + wspsize);
208 Array<OneD, NekDouble> gfac0(tmp0 + 2 * wspsize);
209 Array<OneD, NekDouble> gfac1(tmp0 + 3 * wspsize);
210
211 Array<OneD, NekDouble> tmp1 = outarray[0];
212 Array<OneD, NekDouble> tmp2 = outarray[1];
213
214 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
215 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
216
217 // set up geometric factor: 2/(1-z1)
218 for (int i = 0; i < nquad1; ++i)
219 {
220 gfac0[i] = 2.0 / (1 - z1[i]);
221 }
222 for (int i = 0; i < nquad0; ++i)
223 {
224 gfac1[i] = 0.5 * (1 + z0[i]);
225 }
226
227 for (int i = 0; i < nquad1; ++i)
228 {
229 Vmath::Smul(nquad0, gfac0[i], &inarray[0] + i * nquad0, 1,
230 &tmp0[0] + i * nquad0, 1);
231 }
232
233 for (int i = 0; i < nquad1; ++i)
234 {
235 Vmath::Vmul(nquad0, &gfac1[0], 1, &tmp0[0] + i * nquad0, 1,
236 &tmp1[0] + i * nquad0, 1);
237 }
238
239 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
240 {
241 Vmath::Vmul(nqtot, &df[2 * dir][0], 1, &tmp0[0], 1, &tmp0[0], 1);
242 Vmath::Vmul(nqtot, &df[2 * dir + 1][0], 1, &tmp1[0], 1, &tmp1[0], 1);
243 Vmath::Vmul(nqtot, &df[2 * dir + 1][0], 1, &inarray[0], 1, &tmp2[0], 1);
244 }
245 else
246 {
247 Vmath::Smul(nqtot, df[2 * dir][0], tmp0, 1, tmp0, 1);
248 Vmath::Smul(nqtot, df[2 * dir + 1][0], tmp1, 1, tmp1, 1);
249 Vmath::Smul(nqtot, df[2 * dir + 1][0], inarray, 1, tmp2, 1);
250 }
251 Vmath::Vadd(nqtot, tmp0, 1, tmp1, 1, tmp1, 1);
252}
253
255 const int dir, const Array<OneD, const NekDouble> &inarray,
256 Array<OneD, NekDouble> &outarray)
257{
258 int nq = GetTotPoints();
260
261 switch (dir)
262 {
263 case 0:
264 {
266 }
267 break;
268 case 1:
269 {
271 }
272 break;
273 case 2:
274 {
276 }
277 break;
278 default:
279 {
280 ASSERTL1(false, "input dir is out of range");
281 }
282 break;
283 }
284
285 MatrixKey iprodmatkey(mtype, DetShapeType(), *this);
286 DNekScalMatSharedPtr iprodmat = m_matrixManager[iprodmatkey];
287
288 Blas::Dgemv('N', m_ncoeffs, nq, iprodmat->Scale(),
289 (iprodmat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
290 inarray.get(), 1, 0.0, outarray.get(), 1);
291}
292
293///////////////////////////////
294/// Differentiation Methods
295///////////////////////////////
296
297/**
298 \brief Calculate the deritive of the physical points
299**/
304{
305 int nquad0 = m_base[0]->GetNumPoints();
306 int nquad1 = m_base[1]->GetNumPoints();
307 int nqtot = nquad0 * nquad1;
309 m_metricinfo->GetDerivFactors(GetPointsKeys());
310
311 Array<OneD, NekDouble> diff0(2 * nqtot);
312 Array<OneD, NekDouble> diff1(diff0 + nqtot);
313
314 StdNodalTriExp::v_PhysDeriv(inarray, diff0, diff1);
315
316 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
317 {
318 if (out_d0.size())
319 {
320 Vmath::Vmul(nqtot, df[0], 1, diff0, 1, out_d0, 1);
321 Vmath::Vvtvp(nqtot, df[1], 1, diff1, 1, out_d0, 1, out_d0, 1);
322 }
323
324 if (out_d1.size())
325 {
326 Vmath::Vmul(nqtot, df[2], 1, diff0, 1, out_d1, 1);
327 Vmath::Vvtvp(nqtot, df[3], 1, diff1, 1, out_d1, 1, out_d1, 1);
328 }
329
330 if (out_d2.size())
331 {
332 Vmath::Vmul(nqtot, df[4], 1, diff0, 1, out_d2, 1);
333 Vmath::Vvtvp(nqtot, df[5], 1, diff1, 1, out_d2, 1, out_d2, 1);
334 }
335 }
336 else // regular geometry
337 {
338 if (out_d0.size())
339 {
340 Vmath::Smul(nqtot, df[0][0], diff0, 1, out_d0, 1);
341 Blas::Daxpy(nqtot, df[1][0], diff1, 1, out_d0, 1);
342 }
343
344 if (out_d1.size())
345 {
346 Vmath::Smul(nqtot, df[2][0], diff0, 1, out_d1, 1);
347 Blas::Daxpy(nqtot, df[3][0], diff1, 1, out_d1, 1);
348 }
349
350 if (out_d2.size())
351 {
352 Vmath::Smul(nqtot, df[4][0], diff0, 1, out_d2, 1);
353 Blas::Daxpy(nqtot, df[5][0], diff1, 1, out_d2, 1);
354 }
355 }
356}
357
358/** \brief Forward transform from physical quadrature space
359 stored in \a inarray and evaluate the expansion coefficients and
360 store in \a (this)->m_coeffs
361
362 Inputs:\n
363
364 - \a inarray: array of physical quadrature points to be transformed
365
366 Outputs:\n
367
368 - (this)->_coeffs: updated array of expansion coefficients.
369
370*/
372 Array<OneD, NekDouble> &outarray)
373{
374 IProductWRTBase(inarray, outarray);
375
376 // get Mass matrix inverse
381 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
382
383 // copy inarray in case inarray == outarray
384 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
386
387 out = (*matsys) * in;
388}
389
391 const Array<OneD, const NekDouble> &inarray,
393{
395
396 if (inarray.get() == outarray.get())
397 {
399 Vmath::Vcopy(m_ncoeffs, inarray.get(), 1, tmp.get(), 1);
400
401 Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, mat->Scale(),
402 (mat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
403 tmp.get(), 1, 0.0, outarray.get(), 1);
404 }
405 else
406 {
407 Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, mat->Scale(),
408 (mat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
409 inarray.get(), 1, 0.0, outarray.get(), 1);
410 }
411}
412
414 Array<OneD, NekDouble> &coords_1,
415 Array<OneD, NekDouble> &coords_2)
416{
417 Expansion::v_GetCoords(coords_0, coords_1, coords_2);
418}
419
420// get the coordinates "coords" at the local coordinates "Lcoords"
423{
424 int i;
425
426 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[1] <= 1.0 && Lcoords[1] >= -1.0 &&
427 Lcoords[1] <= 1.0,
428 "Local coordinates are not in region [-1,1]");
429
430 m_geom->FillGeom();
431
432 for (i = 0; i < m_geom->GetCoordim(); ++i)
433 {
434 coords[i] = m_geom->GetCoord(i, Lcoords);
435 }
436}
437
439 const StdRegions::StdMatrixKey &mkey)
440{
441 LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
442 LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
446
447 return tmp->GetStdMatrix(mkey);
448}
449
451 const Array<OneD, const NekDouble> &coord,
452 const Array<OneD, const NekDouble> &physvals)
453
454{
456
457 ASSERTL0(m_geom, "m_geom not defined");
458 m_geom->GetLocCoords(coord, Lcoord);
459
460 return StdExpansion2D::v_PhysEvaluate(Lcoord, physvals);
461}
462
464 const Array<OneD, NekDouble> &coord,
465 const Array<OneD, const NekDouble> &inarray,
466 std::array<NekDouble, 3> &firstOrderDerivs)
467{
468 Array<OneD, NekDouble> Lcoord(2);
469 ASSERTL0(m_geom, "m_geom not defined");
470 m_geom->GetLocCoords(coord, Lcoord);
471 return StdExpansion2D::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
472}
473
475{
476
478 m_base[0]->GetBasisKey(), m_base[1]->GetBasisKey(),
480}
481
483{
485 m_base[0]->GetPointsKey());
487 m_base[1]->GetPointsKey());
488
490 bkey0, bkey1, m_nodalPointsKey.GetPointsType());
491}
492
494{
495 DNekMatSharedPtr returnval;
496
497 switch (mkey.GetMatrixType())
498 {
505 returnval = Expansion2D::v_GenMatrix(mkey);
506 break;
507 default:
508 returnval = StdNodalTriExp::v_GenMatrix(mkey);
509 break;
510 }
511 return returnval;
512}
513
515{
516 int i;
517 const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
518 GetGeom()->GetMetricInfo();
519 const SpatialDomains::GeomType type = geomFactors->GetGtype();
520
523 geomFactors->GetDerivFactors(ptsKeys);
524 const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
525 int nqe = m_base[0]->GetNumPoints();
526 int dim = GetCoordim();
527
530 for (i = 0; i < dim; ++i)
531 {
532 normal[i] = Array<OneD, NekDouble>(nqe);
533 }
534
535 size_t nqb = nqe;
536 size_t nbnd = edge;
539
540 // Regular geometry case
541 if ((type == SpatialDomains::eRegular) ||
543 {
544 NekDouble fac;
545 // Set up normals
546 switch (edge)
547 {
548 case 0:
549 for (i = 0; i < GetCoordim(); ++i)
550 {
551 Vmath::Fill(nqe, -df[2 * i + 1][0], normal[i], 1);
552 }
553 break;
554 case 1:
555 for (i = 0; i < GetCoordim(); ++i)
556 {
557 Vmath::Fill(nqe, df[2 * i + 1][0] + df[2 * i][0], normal[i],
558 1);
559 }
560 break;
561 case 2:
562 for (i = 0; i < GetCoordim(); ++i)
563 {
564 Vmath::Fill(nqe, -df[2 * i][0], normal[i], 1);
565 }
566 break;
567 default:
568 ASSERTL0(false, "Edge is out of range (edge < 3)");
569 }
570
571 // normalise
572 fac = 0.0;
573 for (i = 0; i < GetCoordim(); ++i)
574 {
575 fac += normal[i][0] * normal[i][0];
576 }
577 fac = 1.0 / sqrt(fac);
578
579 Vmath::Fill(nqb, fac, length, 1);
580
581 for (i = 0; i < GetCoordim(); ++i)
582 {
583 Vmath::Smul(nqe, fac, normal[i], 1, normal[i], 1);
584 }
585 }
586 else // Set up deformed normals
587 {
588 int j;
589
590 int nquad0 = ptsKeys[0].GetNumPoints();
591 int nquad1 = ptsKeys[1].GetNumPoints();
592
594
595 Array<OneD, NekDouble> normals(GetCoordim() * max(nquad0, nquad1), 0.0);
596 Array<OneD, NekDouble> edgejac(GetCoordim() * max(nquad0, nquad1), 0.0);
597
598 // Extract Jacobian along edges and recover local
599 // derivates (dx/dr) for polynomial interpolation by
600 // multiplying m_gmat by jacobian
601 switch (edge)
602 {
603 case 0:
604 for (j = 0; j < nquad0; ++j)
605 {
606 edgejac[j] = jac[j];
607 for (i = 0; i < GetCoordim(); ++i)
608 {
609 normals[i * nquad0 + j] =
610 -df[2 * i + 1][j] * edgejac[j];
611 }
612 }
613 from_key = ptsKeys[0];
614 break;
615 case 1:
616 for (j = 0; j < nquad1; ++j)
617 {
618 edgejac[j] = jac[nquad0 * j + nquad0 - 1];
619 for (i = 0; i < GetCoordim(); ++i)
620 {
621 normals[i * nquad1 + j] =
622 (df[2 * i][nquad0 * j + nquad0 - 1] +
623 df[2 * i + 1][nquad0 * j + nquad0 - 1]) *
624 edgejac[j];
625 }
626 }
627 from_key = ptsKeys[1];
628 break;
629 case 2:
630 for (j = 0; j < nquad1; ++j)
631 {
632 edgejac[j] = jac[nquad0 * j];
633 for (i = 0; i < GetCoordim(); ++i)
634 {
635 normals[i * nquad1 + j] =
636 -df[2 * i][nquad0 * j] * edgejac[j];
637 }
638 }
639 from_key = ptsKeys[1];
640 break;
641 default:
642 ASSERTL0(false, "edge is out of range (edge < 3)");
643 }
644
645 int nq = from_key.GetNumPoints();
646 Array<OneD, NekDouble> work(nqe, 0.0);
647
648 // interpolate Jacobian and invert
649 LibUtilities::Interp1D(from_key, jac, m_base[0]->GetPointsKey(), work);
650 Vmath::Sdiv(nq, 1.0, &work[0], 1, &work[0], 1);
651
652 // interpolate
653 for (i = 0; i < GetCoordim(); ++i)
654 {
655 LibUtilities::Interp1D(from_key, &normals[i * nq],
656 m_base[0]->GetPointsKey(), &normal[i][0]);
657 Vmath::Vmul(nqe, work, 1, normal[i], 1, normal[i], 1);
658 }
659
660 // normalise normal vectors
661 Vmath::Zero(nqe, work, 1);
662 for (i = 0; i < GetCoordim(); ++i)
663 {
664 Vmath::Vvtvp(nqe, normal[i], 1, normal[i], 1, work, 1, work, 1);
665 }
666
667 Vmath::Vsqrt(nqe, work, 1, work, 1);
668 Vmath::Sdiv(nqe, 1.0, work, 1, work, 1);
669
670 Vmath::Vcopy(nqb, work, 1, length, 1);
671
672 for (i = 0; i < GetCoordim(); ++i)
673 {
674 Vmath::Vmul(nqe, normal[i], 1, work, 1, normal[i], 1);
675 }
676
677 // Reverse direction so that points are in
678 // anticlockwise direction if edge >=2
679 if (edge >= 2)
680 {
681 for (i = 0; i < GetCoordim(); ++i)
682 {
683 Vmath::Reverse(nqe, normal[i], 1, normal[i], 1);
684 }
685 }
686 }
687}
688
689} // namespace LocalRegions
690} // 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
Defines a specification for a set of points.
Definition: Points.h:55
PointsType GetPointsType() const
Definition: Points.h:95
size_t GetNumPoints() const
Definition: Points.h:90
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 AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Expansion.h:153
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
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:88
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: NodalTriExp.h:198
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)
virtual 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:86
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)
virtual 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:93
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)
virtual 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:44
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: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
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
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
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:213
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
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:236
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< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::shared_ptr< StdNodalTriExp > StdNodalTriExpSharedPtr
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:409
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:353
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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.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 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 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