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