Nektar++
TriExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: TriExp.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// The above copyright notice and this permission notice shall be included
20// in all copies or substantial portions of the Software.
21//
22// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
23// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
24// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
25// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
26// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
27// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
28// DEALINGS IN THE SOFTWARE.
29//
30// Description: Expasion for triangular elements.
31//
32///////////////////////////////////////////////////////////////////////////////
33
34#include <boost/core/ignore_unused.hpp>
35
39#include <LocalRegions/SegExp.h>
40#include <LocalRegions/TriExp.h>
42
43using namespace std;
44
45namespace Nektar
46{
47namespace LocalRegions
48{
50 const LibUtilities::BasisKey &Bb,
52 : StdExpansion(LibUtilities::StdTriData::getNumberOfCoefficients(
53 Ba.GetNumModes(), (Bb.GetNumModes())),
54 2, Ba, Bb),
55 StdExpansion2D(LibUtilities::StdTriData::getNumberOfCoefficients(
56 Ba.GetNumModes(), (Bb.GetNumModes())),
57 Ba, Bb),
58 StdTriExp(Ba, Bb), Expansion(geom), Expansion2D(geom),
59 m_matrixManager(
60 std::bind(&Expansion2D::CreateMatrix, this, std::placeholders::_1),
61 std::string("TriExpMatrix")),
62 m_staticCondMatrixManager(std::bind(&Expansion::CreateStaticCondMatrix,
63 this, std::placeholders::_1),
64 std::string("TriExpStaticCondMatrix"))
65{
66}
67
69 : StdExpansion(T), StdExpansion2D(T), StdTriExp(T), Expansion(T),
70 Expansion2D(T), m_matrixManager(T.m_matrixManager),
71 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
72{
73}
74
76{
77 int nquad0 = m_base[0]->GetNumPoints();
78 int nquad1 = m_base[1]->GetNumPoints();
80 NekDouble ival;
81 Array<OneD, NekDouble> tmp(nquad0 * nquad1);
82
83 // multiply inarray with Jacobian
84 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
85 {
86 Vmath::Vmul(nquad0 * nquad1, jac, 1, inarray, 1, tmp, 1);
87 }
88 else
89 {
90 Vmath::Smul(nquad0 * nquad1, jac[0], inarray, 1, tmp, 1);
91 }
92
93 // call StdQuadExp version;
94 ival = StdTriExp::v_Integral(tmp);
95 return ival;
96}
97
102{
103 int nquad0 = m_base[0]->GetNumPoints();
104 int nquad1 = m_base[1]->GetNumPoints();
105 int nqtot = nquad0 * nquad1;
107 m_metricinfo->GetDerivFactors(GetPointsKeys());
108
109 Array<OneD, NekDouble> diff0(2 * nqtot);
110 Array<OneD, NekDouble> diff1(diff0 + nqtot);
111
112 StdTriExp::v_PhysDeriv(inarray, diff0, diff1);
113
114 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
115 {
116 if (out_d0.size())
117 {
118 Vmath::Vmul(nqtot, df[0], 1, diff0, 1, out_d0, 1);
119 Vmath::Vvtvp(nqtot, df[1], 1, diff1, 1, out_d0, 1, out_d0, 1);
120 }
121
122 if (out_d1.size())
123 {
124 Vmath::Vmul(nqtot, df[2], 1, diff0, 1, out_d1, 1);
125 Vmath::Vvtvp(nqtot, df[3], 1, diff1, 1, out_d1, 1, out_d1, 1);
126 }
127
128 if (out_d2.size())
129 {
130 Vmath::Vmul(nqtot, df[4], 1, diff0, 1, out_d2, 1);
131 Vmath::Vvtvp(nqtot, df[5], 1, diff1, 1, out_d2, 1, out_d2, 1);
132 }
133 }
134 else // regular geometry
135 {
136 if (out_d0.size())
137 {
138 Vmath::Smul(nqtot, df[0][0], diff0, 1, out_d0, 1);
139 Blas::Daxpy(nqtot, df[1][0], diff1, 1, out_d0, 1);
140 }
141
142 if (out_d1.size())
143 {
144 Vmath::Smul(nqtot, df[2][0], diff0, 1, out_d1, 1);
145 Blas::Daxpy(nqtot, df[3][0], diff1, 1, out_d1, 1);
146 }
147
148 if (out_d2.size())
149 {
150 Vmath::Smul(nqtot, df[4][0], diff0, 1, out_d2, 1);
151 Blas::Daxpy(nqtot, df[5][0], diff1, 1, out_d2, 1);
152 }
153 }
154}
155
156void TriExp::v_PhysDeriv(const int dir,
157 const Array<OneD, const NekDouble> &inarray,
158 Array<OneD, NekDouble> &outarray)
159{
160 switch (dir)
161 {
162 case 0:
163 {
164 PhysDeriv(inarray, outarray, NullNekDouble1DArray,
166 }
167 break;
168 case 1:
169 {
170 PhysDeriv(inarray, NullNekDouble1DArray, outarray,
172 }
173 break;
174 case 2:
175 {
177 outarray);
178 }
179 break;
180 default:
181 {
182 ASSERTL1(false, "input dir is out of range");
183 }
184 break;
185 }
186}
187
189 const Array<OneD, const NekDouble> &inarray,
191{
192 if (!out.size())
193 {
194 return;
195 }
196
197 int nquad0 = m_base[0]->GetNumPoints();
198 int nquad1 = m_base[1]->GetNumPoints();
199 int nqtot = nquad0 * nquad1;
200
202 m_metricinfo->GetDerivFactors(GetPointsKeys());
203
204 Array<OneD, NekDouble> diff0(2 * nqtot);
205 Array<OneD, NekDouble> diff1(diff0 + nqtot);
206
207 // diff0 = du/d_xi, diff1 = du/d_eta
208 StdTriExp::v_PhysDeriv(inarray, diff0, diff1);
209
210 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
211 {
213
214 // D^v_xi = v_x*d_xi/dx + v_y*d_xi/dy + v_z*d_xi/dz
215 // D^v_eta = v_x*d_eta/dx + v_y*d_eta/dy + v_z*d_eta/dz
216 for (int i = 0; i < 2; ++i)
217 {
218 tangmat[i] = Array<OneD, NekDouble>(nqtot, 0.0);
219 for (int k = 0; k < (m_geom->GetCoordim()); ++k)
220 {
221 Vmath::Vvtvp(nqtot, &df[2 * k + i][0], 1, &direction[k * nqtot],
222 1, &tangmat[i][0], 1, &tangmat[i][0], 1);
223 }
224 }
225
226 /// D_v = D^v_xi * du/d_xi + D^v_eta * du/d_eta
227 Vmath::Vmul(nqtot, &tangmat[0][0], 1, &diff0[0], 1, &out[0], 1);
228 Vmath::Vvtvp(nqtot, &tangmat[1][0], 1, &diff1[0], 1, &out[0], 1,
229 &out[0], 1);
230 }
231 else
232 {
234
235 for (int i = 0; i < 2; ++i)
236 {
237 tangmat[i] = Array<OneD, NekDouble>(nqtot, 0.0);
238 for (int k = 0; k < (m_geom->GetCoordim()); ++k)
239 {
240 Vmath::Svtvp(nqtot, df[2 * k + i][0], &direction[k * nqtot], 1,
241 &tangmat[i][0], 1, &tangmat[i][0], 1);
242 }
243 }
244
245 /// D_v = D^v_xi * du/d_xi + D^v_eta * du/d_eta
246 Vmath::Vmul(nqtot, &tangmat[0][0], 1, &diff0[0], 1, &out[0], 1);
247
248 Vmath::Vvtvp(nqtot, &tangmat[1][0], 1, &diff1[0], 1, &out[0], 1,
249 &out[0], 1);
250 }
251}
252
254 Array<OneD, NekDouble> &outarray)
255{
256 IProductWRTBase(inarray, outarray);
257
258 // get Mass matrix inverse
260 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
261
262 // copy inarray in case inarray == outarray
263 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
265
266 out = (*matsys) * in;
267}
268
270 const Array<OneD, const NekDouble> &inarray,
271 Array<OneD, NekDouble> &outarray)
272{
273 int i, j;
274 int npoints[2] = {m_base[0]->GetNumPoints(), m_base[1]->GetNumPoints()};
275 int nmodes[2] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes()};
276
277 fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
278
279 if (nmodes[0] == 1 && nmodes[1] == 1)
280 {
281 outarray[0] = inarray[0];
282 return;
283 }
284
285 Array<OneD, NekDouble> physEdge[3];
286 Array<OneD, NekDouble> coeffEdge[3];
287 for (i = 0; i < 3; i++)
288 {
289 // define physEdge and add 1 so can interpolate grl10 points if
290 // necessary
291 physEdge[i] = Array<OneD, NekDouble>(max(npoints[i != 0], npoints[0]));
292 coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i != 0]);
293 }
294
295 for (i = 0; i < npoints[0]; i++)
296 {
297 physEdge[0][i] = inarray[i];
298 }
299
300 // extract data in cartesian directions
301 for (i = 0; i < npoints[1]; i++)
302 {
303 physEdge[1][i] = inarray[npoints[0] - 1 + i * npoints[0]];
304 physEdge[2][i] = inarray[i * npoints[0]];
305 }
306
307 SegExpSharedPtr segexp[3];
309 m_base[0]->GetBasisKey(), GetGeom2D()->GetEdge(0));
310
312 {
313 for (i = 1; i < 3; i++)
314 {
316 m_base[i != 0]->GetBasisKey(), GetGeom2D()->GetEdge(i));
317 }
318 }
319 else // interploate using edge 0 GLL distribution
320 {
321 for (i = 1; i < 3; i++)
322 {
324 m_base[0]->GetBasisKey(), GetGeom2D()->GetEdge(i));
325
326 LibUtilities::Interp1D(m_base[1]->GetPointsKey(), physEdge[i],
327 m_base[0]->GetPointsKey(), physEdge[i]);
328 }
329 npoints[1] = npoints[0];
330 }
331
333 Array<OneD, int> signArray;
335 // define an orientation to get EdgeToElmtMapping from Cartesian data
336 StdRegions::Orientation orient[3] = {
338
339 for (i = 0; i < 3; i++)
340 {
341 segexp[i]->FwdTransBndConstrained(physEdge[i], coeffEdge[i]);
342
343 // this orient goes with the one above and so could
344 // probably set both to eForwards
345 GetTraceToElementMap(i, mapArray, signArray, orient[i]);
346 for (j = 0; j < nmodes[i != 0]; j++)
347 {
348 sign = (NekDouble)signArray[j];
349 outarray[mapArray[j]] = sign * coeffEdge[i][j];
350 }
351 }
352
353 int nBoundaryDofs = NumBndryCoeffs();
354 int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
355
356 if (nInteriorDofs > 0)
357 {
360
362 *this);
363 MassMatrixOp(outarray, tmp0, stdmasskey);
364 IProductWRTBase(inarray, tmp1);
365
366 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
367
368 // get Mass matrix inverse (only of interior DOF)
369 // use block (1,1) of the static condensed system
370 // note: this block alreay contains the inverse matrix
371 MatrixKey masskey(StdRegions::eMass, DetShapeType(), *this);
372 DNekScalMatSharedPtr matsys =
373 (m_staticCondMatrixManager[masskey])->GetBlock(1, 1);
374
375 Array<OneD, NekDouble> rhs(nInteriorDofs);
376 Array<OneD, NekDouble> result(nInteriorDofs);
377
378 GetInteriorMap(mapArray);
379
380 for (i = 0; i < nInteriorDofs; i++)
381 {
382 rhs[i] = tmp1[mapArray[i]];
383 }
384
385 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, matsys->Scale(),
386 &((matsys->GetOwnedMatrix())->GetPtr())[0], nInteriorDofs,
387 rhs.get(), 1, 0.0, result.get(), 1);
388
389 for (i = 0; i < nInteriorDofs; i++)
390 {
391 outarray[mapArray[i]] = result[i];
392 }
393 }
394}
395
397 Array<OneD, NekDouble> &outarray)
398{
399 IProductWRTBase_SumFac(inarray, outarray);
400}
401
403 const Array<OneD, const NekDouble> &inarray,
404 Array<OneD, NekDouble> &outarray)
405{
406 IProductWRTDerivBase_SumFac(dir, inarray, outarray);
407}
408
410 const Array<OneD, const NekDouble> &inarray,
411 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
412{
413 int nquad0 = m_base[0]->GetNumPoints();
414 int nquad1 = m_base[1]->GetNumPoints();
415 int order0 = m_base[0]->GetNumModes();
416
417 if (multiplybyweights)
418 {
419 Array<OneD, NekDouble> tmp(nquad0 * nquad1 + nquad1 * order0);
420 Array<OneD, NekDouble> wsp(tmp + nquad0 * nquad1);
421
422 MultiplyByQuadratureMetric(inarray, tmp);
424 m_base[1]->GetBdata(), tmp, outarray, wsp);
425 }
426 else
427 {
428 Array<OneD, NekDouble> wsp(+nquad1 * order0);
429
431 m_base[1]->GetBdata(), inarray, outarray,
432 wsp);
433 }
434}
435
437 const int dir, const Array<OneD, const NekDouble> &inarray,
438 Array<OneD, NekDouble> &outarray)
439{
440 int nquad0 = m_base[0]->GetNumPoints();
441 int nquad1 = m_base[1]->GetNumPoints();
442 int nqtot = nquad0 * nquad1;
443 int nmodes0 = m_base[0]->GetNumModes();
444 int wspsize = max(max(nqtot, m_ncoeffs), nquad1 * nmodes0);
445
446 Array<OneD, NekDouble> tmp0(4 * wspsize);
447 Array<OneD, NekDouble> tmp1(tmp0 + wspsize);
448 Array<OneD, NekDouble> tmp2(tmp0 + 2 * wspsize);
449 Array<OneD, NekDouble> tmp3(tmp0 + 3 * wspsize);
450
452 tmp2D[0] = tmp1;
453 tmp2D[1] = tmp2;
454
455 TriExp::v_AlignVectorToCollapsedDir(dir, inarray, tmp2D);
456
457 MultiplyByQuadratureMetric(tmp1, tmp1);
458 MultiplyByQuadratureMetric(tmp2, tmp2);
459
460 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
461 tmp1, tmp3, tmp0);
462 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
463 tmp2, outarray, tmp0);
464 Vmath::Vadd(m_ncoeffs, tmp3, 1, outarray, 1, outarray, 1);
465}
466
468 const int dir, const Array<OneD, const NekDouble> &inarray,
470{
471 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
472 ASSERTL1((dir == 2) ? (m_geom->GetCoordim() == 3) : true,
473 "Invalid direction.");
474
475 int nquad0 = m_base[0]->GetNumPoints();
476 int nquad1 = m_base[1]->GetNumPoints();
477 int nqtot = nquad0 * nquad1;
478 int nmodes0 = m_base[0]->GetNumModes();
479 int wspsize = max(max(nqtot, m_ncoeffs), nquad1 * nmodes0);
480
482 m_metricinfo->GetDerivFactors(GetPointsKeys());
483
484 Array<OneD, NekDouble> tmp0(wspsize);
485 Array<OneD, NekDouble> tmp3(wspsize);
486 Array<OneD, NekDouble> gfac0(wspsize);
487 Array<OneD, NekDouble> gfac1(wspsize);
488
489 Array<OneD, NekDouble> tmp1 = outarray[0];
490 Array<OneD, NekDouble> tmp2 = outarray[1];
491
492 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
493 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
494
495 // set up geometric factor: 2/(1-z1)
496 for (int i = 0; i < nquad1; ++i)
497 {
498 gfac0[i] = 2.0 / (1 - z1[i]);
499 }
500 for (int i = 0; i < nquad0; ++i)
501 {
502 gfac1[i] = 0.5 * (1 + z0[i]);
503 }
504
505 for (int i = 0; i < nquad1; ++i)
506 {
507 Vmath::Smul(nquad0, gfac0[i], &inarray[0] + i * nquad0, 1,
508 &tmp0[0] + i * nquad0, 1);
509 }
510
511 for (int i = 0; i < nquad1; ++i)
512 {
513 Vmath::Vmul(nquad0, &gfac1[0], 1, &tmp0[0] + i * nquad0, 1,
514 &tmp1[0] + i * nquad0, 1);
515 }
516
517 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
518 {
519 Vmath::Vmul(nqtot, &df[2 * dir][0], 1, &tmp0[0], 1, &tmp0[0], 1);
520 Vmath::Vmul(nqtot, &df[2 * dir + 1][0], 1, &tmp1[0], 1, &tmp1[0], 1);
521 Vmath::Vmul(nqtot, &df[2 * dir + 1][0], 1, &inarray[0], 1, &tmp2[0], 1);
522 }
523 else
524 {
525 Vmath::Smul(nqtot, df[2 * dir][0], tmp0, 1, tmp0, 1);
526 Vmath::Smul(nqtot, df[2 * dir + 1][0], tmp1, 1, tmp1, 1);
527 Vmath::Smul(nqtot, df[2 * dir + 1][0], inarray, 1, tmp2, 1);
528 }
529 Vmath::Vadd(nqtot, tmp0, 1, tmp1, 1, tmp1, 1);
530}
531
533 const Array<OneD, const NekDouble> &direction,
534 const Array<OneD, const NekDouble> &inarray,
535 Array<OneD, NekDouble> &outarray)
536{
537 IProductWRTDirectionalDerivBase_SumFac(direction, inarray, outarray);
538}
539
540/**
541 * @brief Directinoal Derivative in the modal space in the dir
542 * direction of varcoeffs.
543 */
545 const Array<OneD, const NekDouble> &direction,
546 const Array<OneD, const NekDouble> &inarray,
547 Array<OneD, NekDouble> &outarray)
548{
549 int i;
550 int shapedim = 2;
551 int nquad0 = m_base[0]->GetNumPoints();
552 int nquad1 = m_base[1]->GetNumPoints();
553 int nqtot = nquad0 * nquad1;
554 int nmodes0 = m_base[0]->GetNumModes();
555 int wspsize = max(max(nqtot, m_ncoeffs), nquad1 * nmodes0);
556
558 m_metricinfo->GetDerivFactors(GetPointsKeys());
559
560 Array<OneD, NekDouble> tmp0(6 * wspsize);
561 Array<OneD, NekDouble> tmp1(tmp0 + wspsize);
562 Array<OneD, NekDouble> tmp2(tmp0 + 2 * wspsize);
563 Array<OneD, NekDouble> tmp3(tmp0 + 3 * wspsize);
564 Array<OneD, NekDouble> gfac0(tmp0 + 4 * wspsize);
565 Array<OneD, NekDouble> gfac1(tmp0 + 5 * wspsize);
566
567 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
568 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
569
570 // set up geometric factor: 2/(1-z1)
571 for (i = 0; i < nquad1; ++i)
572 {
573 gfac0[i] = 2.0 / (1 - z1[i]);
574 }
575 for (i = 0; i < nquad0; ++i)
576 {
577 gfac1[i] = 0.5 * (1 + z0[i]);
578 }
579 for (i = 0; i < nquad1; ++i)
580 {
581 Vmath::Smul(nquad0, gfac0[i], &inarray[0] + i * nquad0, 1,
582 &tmp0[0] + i * nquad0, 1);
583 }
584 for (i = 0; i < nquad1; ++i)
585 {
586 Vmath::Vmul(nquad0, &gfac1[0], 1, &tmp0[0] + i * nquad0, 1,
587 &tmp1[0] + i * nquad0, 1);
588 }
589
590 // Compute gmat \cdot e^j
591 Array<OneD, Array<OneD, NekDouble>> dfdir(shapedim);
592 Expansion::ComputeGmatcdotMF(df, direction, dfdir);
593
594 Vmath::Vmul(nqtot, &dfdir[0][0], 1, &tmp0[0], 1, &tmp0[0], 1);
595 Vmath::Vmul(nqtot, &dfdir[1][0], 1, &tmp1[0], 1, &tmp1[0], 1);
596 Vmath::Vmul(nqtot, &dfdir[1][0], 1, &inarray[0], 1, &tmp2[0], 1);
597
598 Vmath::Vadd(nqtot, &tmp0[0], 1, &tmp1[0], 1, &tmp1[0], 1);
599
600 MultiplyByQuadratureMetric(tmp1, tmp1);
601 MultiplyByQuadratureMetric(tmp2, tmp2);
602
603 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
604 tmp1, tmp3, tmp0);
605 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
606 tmp2, outarray, tmp0);
607 Vmath::Vadd(m_ncoeffs, tmp3, 1, outarray, 1, outarray, 1);
608}
609
613 Array<OneD, NekDouble> &outarray)
614{
615 int nq = m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints();
617
619 GetLeftAdjacentElementExp()->GetTraceNormal(
621
622 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
623 {
624 Vmath::Vvtvvtp(nq, &normals[0][0], 1, &Fx[0], 1, &normals[1][0], 1,
625 &Fy[0], 1, &Fn[0], 1);
626 Vmath::Vvtvp(nq, &normals[2][0], 1, &Fz[0], 1, &Fn[0], 1, &Fn[0], 1);
627 }
628 else
629 {
630 Vmath::Svtsvtp(nq, normals[0][0], &Fx[0], 1, normals[1][0], &Fy[0], 1,
631 &Fn[0], 1);
632 Vmath::Svtvp(nq, normals[2][0], &Fz[0], 1, &Fn[0], 1, &Fn[0], 1);
633 }
634
635 IProductWRTBase(Fn, outarray);
636}
637
639 const Array<OneD, const Array<OneD, NekDouble>> &Fvec,
640 Array<OneD, NekDouble> &outarray)
641{
642 NormVectorIProductWRTBase(Fvec[0], Fvec[1], Fvec[2], outarray);
643}
644
646{
647
649 m_base[0]->GetBasisKey(), m_base[1]->GetBasisKey());
650}
651
653{
655 m_base[0]->GetPointsKey());
657 m_base[1]->GetPointsKey());
658
660 bkey1);
661}
662
665{
666 int i;
667
668 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[1] <= 1.0 && Lcoords[1] >= -1.0 &&
669 Lcoords[1] <= 1.0,
670 "Local coordinates are not in region [-1,1]");
671
672 m_geom->FillGeom();
673
674 for (i = 0; i < m_geom->GetCoordim(); ++i)
675 {
676 coords[i] = m_geom->GetCoord(i, Lcoords);
677 }
678}
679
681 Array<OneD, NekDouble> &coords_1,
682 Array<OneD, NekDouble> &coords_2)
683{
684 Expansion::v_GetCoords(coords_0, coords_1, coords_2);
685}
686
687/**
688 * Given the local cartesian coordinate \a Lcoord evaluate the
689 * value of physvals at this point by calling through to the
690 * StdExpansion method
691 */
693 const Array<OneD, const NekDouble> &Lcoord,
694 const Array<OneD, const NekDouble> &physvals)
695{
696 // Evaluate point in local (eta) coordinates.
697 return StdExpansion2D::v_PhysEvaluate(Lcoord, physvals);
698}
699
701 const Array<OneD, const NekDouble> &physvals)
702{
704
705 ASSERTL0(m_geom, "m_geom not defined");
706 m_geom->GetLocCoords(coord, Lcoord);
707
708 return StdExpansion2D::v_PhysEvaluate(Lcoord, physvals);
709}
710
712 const Array<OneD, const NekDouble> &inarray,
713 std::array<NekDouble, 3> &firstOrderDerivs)
714{
715 Array<OneD, NekDouble> Lcoord(2);
716 ASSERTL0(m_geom, "m_geom not defined");
717 m_geom->GetLocCoords(coord, Lcoord);
718 return StdTriExp::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
719}
720
722 const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp,
723 const Array<OneD, const NekDouble> &inarray,
725{
726 int nquad0 = m_base[0]->GetNumPoints();
727 int nquad1 = m_base[1]->GetNumPoints();
728 int nt = 0;
729 // Extract in Cartesian direction because we have to deal with
730 // e.g. Gauss-Radau points.
731 switch (edge)
732 {
733 case 0:
734 Vmath::Vcopy(nquad0, &(inarray[0]), 1, &(outarray[0]), 1);
735 nt = nquad0;
736 break;
737 case 1:
738 Vmath::Vcopy(nquad1, &(inarray[0]) + (nquad0 - 1), nquad0,
739 &(outarray[0]), 1);
740 nt = nquad1;
741 break;
742 case 2:
743 Vmath::Vcopy(nquad1, &(inarray[0]), nquad0, &(outarray[0]), 1);
744 nt = nquad1;
745 break;
746 default:
747 ASSERTL0(false, "edge value (< 3) is out of range");
748 break;
749 }
750
751 ASSERTL1(EdgeExp->GetBasis(0)->GetPointsType() ==
753 "Edge expansion should be GLL");
754
755 // Interpolate if required
756 if (m_base[edge ? 1 : 0]->GetPointsKey() !=
757 EdgeExp->GetBasis(0)->GetPointsKey())
758 {
759 Array<OneD, NekDouble> outtmp(max(nquad0, nquad1));
760
761 Vmath::Vcopy(nt, outarray, 1, outtmp, 1);
762
763 LibUtilities::Interp1D(m_base[edge ? 1 : 0]->GetPointsKey(), outtmp,
764 EdgeExp->GetBasis(0)->GetPointsKey(), outarray);
765 }
766
767 if (orient == StdRegions::eNoOrientation)
768 {
769 orient = GetTraceOrient(edge);
770 }
771
772 // Reverse data if necessary
773 if (orient == StdRegions::eBackwards)
774 {
775 Vmath::Reverse(EdgeExp->GetNumPoints(0), &outarray[0], 1, &outarray[0],
776 1);
777 }
778}
779
780void TriExp::v_GetTraceQFactors(const int edge,
781 Array<OneD, NekDouble> &outarray)
782{
783 boost::ignore_unused(edge, outarray);
784 ASSERTL0(false, "Routine not implemented for triangular elements");
785}
786
787void TriExp::v_GetTracePhysMap(const int edge, Array<OneD, int> &outarray)
788{
789 int nquad0 = m_base[0]->GetNumPoints();
790 int nquad1 = m_base[1]->GetNumPoints();
791
792 // Get points in Cartesian orientation
793 switch (edge)
794 {
795 case 0:
796 outarray = Array<OneD, int>(nquad0);
797 for (int i = 0; i < nquad0; ++i)
798 {
799 outarray[i] = i;
800 }
801 break;
802 case 1:
803 outarray = Array<OneD, int>(nquad1);
804 for (int i = 0; i < nquad1; ++i)
805 {
806 outarray[i] = (nquad0 - 1) + i * nquad0;
807 }
808 break;
809 case 2:
810 outarray = Array<OneD, int>(nquad1);
811 for (int i = 0; i < nquad1; ++i)
812 {
813 outarray[i] = i * nquad0;
814 }
815 break;
816 default:
817 ASSERTL0(false, "edge value (< 3) is out of range");
818 break;
819 }
820}
821
823{
824 int i;
825 const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
826 GetGeom()->GetMetricInfo();
827
829 for (i = 0; i < ptsKeys.size(); ++i)
830 {
831 // Need at least 2 points for computing normals
832 if (ptsKeys[i].GetNumPoints() == 1)
833 {
834 LibUtilities::PointsKey pKey(2, ptsKeys[i].GetPointsType());
835 ptsKeys[i] = pKey;
836 }
837 }
838
839 const SpatialDomains::GeomType type = geomFactors->GetGtype();
841 geomFactors->GetDerivFactors(ptsKeys);
842 const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
843
844 // The points of normals should follow trace basis, not local basis.
846
847 int nqe = tobasis.GetNumPoints();
848 int dim = GetCoordim();
849
852 for (i = 0; i < dim; ++i)
853 {
854 normal[i] = Array<OneD, NekDouble>(nqe);
855 }
856
857 size_t nqb = nqe;
858 size_t nbnd = edge;
861
862 // Regular geometry case
863 if ((type == SpatialDomains::eRegular) ||
865 {
866 NekDouble fac;
867 // Set up normals
868 switch (edge)
869 {
870 case 0:
871 for (i = 0; i < GetCoordim(); ++i)
872 {
873 Vmath::Fill(nqe, -df[2 * i + 1][0], normal[i], 1);
874 }
875 break;
876 case 1:
877 for (i = 0; i < GetCoordim(); ++i)
878 {
879 Vmath::Fill(nqe, df[2 * i + 1][0] + df[2 * i][0], normal[i],
880 1);
881 }
882 break;
883 case 2:
884 for (i = 0; i < GetCoordim(); ++i)
885 {
886 Vmath::Fill(nqe, -df[2 * i][0], normal[i], 1);
887 }
888 break;
889 default:
890 ASSERTL0(false, "Edge is out of range (edge < 3)");
891 }
892
893 // normalise
894 fac = 0.0;
895 for (i = 0; i < GetCoordim(); ++i)
896 {
897 fac += normal[i][0] * normal[i][0];
898 }
899 fac = 1.0 / sqrt(fac);
900
901 Vmath::Fill(nqb, fac, length, 1);
902
903 for (i = 0; i < GetCoordim(); ++i)
904 {
905 Vmath::Smul(nqe, fac, normal[i], 1, normal[i], 1);
906 }
907 }
908 else // Set up deformed normals
909 {
910 int j;
911
912 int nquad0 = ptsKeys[0].GetNumPoints();
913 int nquad1 = ptsKeys[1].GetNumPoints();
914
916
917 Array<OneD, NekDouble> normals(GetCoordim() * max(nquad0, nquad1), 0.0);
918 Array<OneD, NekDouble> edgejac(GetCoordim() * max(nquad0, nquad1), 0.0);
919
920 // Extract Jacobian along edges and recover local
921 // derivates (dx/dr) for polynomial interpolation by
922 // multiplying m_gmat by jacobian
923 switch (edge)
924 {
925 case 0:
926 for (j = 0; j < nquad0; ++j)
927 {
928 edgejac[j] = jac[j];
929 for (i = 0; i < GetCoordim(); ++i)
930 {
931 normals[i * nquad0 + j] =
932 -df[2 * i + 1][j] * edgejac[j];
933 }
934 }
935 from_key = ptsKeys[0];
936 break;
937 case 1:
938 for (j = 0; j < nquad1; ++j)
939 {
940 edgejac[j] = jac[nquad0 * j + nquad0 - 1];
941 for (i = 0; i < GetCoordim(); ++i)
942 {
943 normals[i * nquad1 + j] =
944 (df[2 * i][nquad0 * j + nquad0 - 1] +
945 df[2 * i + 1][nquad0 * j + nquad0 - 1]) *
946 edgejac[j];
947 }
948 }
949 from_key = ptsKeys[1];
950 break;
951 case 2:
952 for (j = 0; j < nquad1; ++j)
953 {
954 edgejac[j] = jac[nquad0 * j];
955 for (i = 0; i < GetCoordim(); ++i)
956 {
957 normals[i * nquad1 + j] =
958 -df[2 * i][nquad0 * j] * edgejac[j];
959 }
960 }
961 from_key = ptsKeys[1];
962 break;
963 default:
964 ASSERTL0(false, "edge is out of range (edge < 3)");
965 }
966
967 int nq = from_key.GetNumPoints();
968 Array<OneD, NekDouble> work(nqe, 0.0);
969
970 // interpolate Jacobian and invert
971 LibUtilities::Interp1D(from_key, jac, tobasis.GetPointsKey(), work);
972 Vmath::Sdiv(nqe, 1.0, &work[0], 1, &work[0], 1);
973
974 // interpolate
975 for (i = 0; i < GetCoordim(); ++i)
976 {
977 LibUtilities::Interp1D(from_key, &normals[i * nq],
978 tobasis.GetPointsKey(), &normal[i][0]);
979 Vmath::Vmul(nqe, work, 1, normal[i], 1, normal[i], 1);
980 }
981
982 // normalise normal vectors
983 Vmath::Zero(nqe, work, 1);
984 for (i = 0; i < GetCoordim(); ++i)
985 {
986 Vmath::Vvtvp(nqe, normal[i], 1, normal[i], 1, work, 1, work, 1);
987 }
988
989 Vmath::Vsqrt(nqe, work, 1, work, 1);
990 Vmath::Sdiv(nqe, 1.0, work, 1, work, 1);
991
992 Vmath::Vcopy(nqb, work, 1, length, 1);
993
994 for (i = 0; i < GetCoordim(); ++i)
995 {
996 Vmath::Vmul(nqe, normal[i], 1, work, 1, normal[i], 1);
997 }
998 }
999
1000 if (GetGeom()->GetEorient(edge) == StdRegions::eBackwards)
1001 {
1002 for (i = 0; i < GetCoordim(); ++i)
1003 {
1004 if (geomFactors->GetGtype() == SpatialDomains::eDeformed)
1005 {
1006 Vmath::Reverse(nqe, normal[i], 1, normal[i], 1);
1007 }
1008 }
1009 }
1010}
1011
1013 const NekDouble *data, const std::vector<unsigned int> &nummodes,
1014 const int mode_offset, NekDouble *coeffs,
1015 std::vector<LibUtilities::BasisType> &fromType)
1016{
1017 boost::ignore_unused(fromType);
1018
1019 int data_order0 = nummodes[mode_offset];
1020 int fillorder0 = min(m_base[0]->GetNumModes(), data_order0);
1021 int data_order1 = nummodes[mode_offset + 1];
1022 int order1 = m_base[1]->GetNumModes();
1023 int fillorder1 = min(order1, data_order1);
1024
1025 switch (m_base[0]->GetBasisType())
1026 {
1029 {
1030 int i;
1031 int cnt = 0;
1032 int cnt1 = 0;
1033
1036 "Extraction routine not set up for this basis");
1037
1038 Vmath::Zero(m_ncoeffs, coeffs, 1);
1039 for (i = 0; i < fillorder0; ++i)
1040 {
1041 Vmath::Vcopy(fillorder1 - i, &data[cnt], 1, &coeffs[cnt1], 1);
1042 cnt += data_order1 - i;
1043 cnt1 += order1 - i;
1044 }
1045 }
1046 break;
1047 default:
1048 ASSERTL0(false, "basis is either not set up or not hierarchicial");
1049 }
1050}
1051
1053{
1054 return GetGeom2D()->GetEorient(edge);
1055}
1056
1058{
1059 DNekMatSharedPtr returnval;
1060 switch (mkey.GetMatrixType())
1061 {
1069 returnval = Expansion2D::v_GenMatrix(mkey);
1070 break;
1071 default:
1072 returnval = StdTriExp::v_GenMatrix(mkey);
1073 break;
1074 }
1075
1076 return returnval;
1077}
1078
1080{
1081 LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
1082 LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
1085
1086 return tmp->GetStdMatrix(mkey);
1087}
1088
1090{
1091 return m_matrixManager[mkey];
1092}
1093
1095{
1096 m_matrixManager.DeleteObject(mkey);
1097}
1098
1100{
1101 return m_staticCondMatrixManager[mkey];
1102}
1103
1105{
1106 m_staticCondMatrixManager.DeleteObject(mkey);
1107}
1108
1110 Array<OneD, NekDouble> &outarray,
1111 const StdRegions::StdMatrixKey &mkey)
1112{
1113 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1114}
1115
1117 Array<OneD, NekDouble> &outarray,
1118 const StdRegions::StdMatrixKey &mkey)
1119{
1120 TriExp::LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1121}
1122
1123void TriExp::v_LaplacianMatrixOp(const int k1, const int k2,
1124 const Array<OneD, const NekDouble> &inarray,
1125 Array<OneD, NekDouble> &outarray,
1126 const StdRegions::StdMatrixKey &mkey)
1127{
1128 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1129}
1130
1132 const Array<OneD, const NekDouble> &inarray,
1133 Array<OneD, NekDouble> &outarray,
1134 const StdRegions::StdMatrixKey &mkey)
1135{
1136 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1137}
1138
1140 const Array<OneD, const NekDouble> &inarray,
1141 Array<OneD, NekDouble> &outarray, const StdRegions::StdMatrixKey &mkey)
1142{
1143 StdExpansion::WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
1144}
1145
1147 const Array<OneD, const NekDouble> &inarray,
1148 Array<OneD, NekDouble> &outarray, const StdRegions::StdMatrixKey &mkey)
1149{
1150 StdExpansion::MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
1151}
1152
1154 Array<OneD, NekDouble> &outarray,
1155 const StdRegions::StdMatrixKey &mkey)
1156{
1157 TriExp::HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1158}
1159
1161 const Array<OneD, const NekDouble> &inarray,
1163{
1164 if (m_metrics.count(eMetricLaplacian00) == 0)
1165 {
1167 }
1168
1169 int nquad0 = m_base[0]->GetNumPoints();
1170 int nquad1 = m_base[1]->GetNumPoints();
1171 int nqtot = nquad0 * nquad1;
1172 int nmodes0 = m_base[0]->GetNumModes();
1173 int nmodes1 = m_base[1]->GetNumModes();
1174 int wspsize =
1175 max(max(max(nqtot, m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
1176
1177 ASSERTL1(wsp.size() >= 3 * wspsize, "Workspace is of insufficient size.");
1178
1179 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
1180 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
1181 const Array<OneD, const NekDouble> &dbase0 = m_base[0]->GetDbdata();
1182 const Array<OneD, const NekDouble> &dbase1 = m_base[1]->GetDbdata();
1183 const Array<OneD, const NekDouble> &metric00 =
1185 const Array<OneD, const NekDouble> &metric01 =
1187 const Array<OneD, const NekDouble> &metric11 =
1189
1190 // Allocate temporary storage
1191 Array<OneD, NekDouble> wsp0(wsp);
1192 Array<OneD, NekDouble> wsp1(wsp + wspsize);
1193 Array<OneD, NekDouble> wsp2(wsp + 2 * wspsize);
1194
1195 StdExpansion2D::PhysTensorDeriv(inarray, wsp1, wsp2);
1196
1197 // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1198 // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1199 // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1200 // especially for this purpose
1201 Vmath::Vvtvvtp(nqtot, &metric00[0], 1, &wsp1[0], 1, &metric01[0], 1,
1202 &wsp2[0], 1, &wsp0[0], 1);
1203 Vmath::Vvtvvtp(nqtot, &metric01[0], 1, &wsp1[0], 1, &metric11[0], 1,
1204 &wsp2[0], 1, &wsp2[0], 1);
1205
1206 // outarray = m = (D_xi1 * B)^T * k
1207 // wsp1 = n = (D_xi2 * B)^T * l
1208 IProductWRTBase_SumFacKernel(dbase0, base1, wsp0, outarray, wsp1);
1209 IProductWRTBase_SumFacKernel(base0, dbase1, wsp2, wsp1, wsp0);
1210
1211 // outarray = outarray + wsp1
1212 // = L * u_hat
1213 Vmath::Vadd(m_ncoeffs, wsp1.get(), 1, outarray.get(), 1, outarray.get(), 1);
1214}
1215
1217{
1218 if (m_metrics.count(eMetricQuadrature) == 0)
1219 {
1221 }
1222
1223 unsigned int i, j;
1224 const SpatialDomains::GeomType type = m_metricinfo->GetGtype();
1225 const unsigned int nqtot = GetTotPoints();
1226 const unsigned int dim = 2;
1227 const MetricType m[3][3] = {
1231
1232 Array<OneD, NekDouble> dEta_dXi[2] = {Array<OneD, NekDouble>(nqtot, 1.0),
1233 Array<OneD, NekDouble>(nqtot, 1.0)};
1234
1235 for (i = 0; i < dim; ++i)
1236 {
1237 for (j = i; j < dim; ++j)
1238 {
1239 m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1240 }
1241 }
1242
1243 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
1244 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
1245 const unsigned int nquad0 = m_base[0]->GetNumPoints();
1246 const unsigned int nquad1 = m_base[1]->GetNumPoints();
1248 m_metricinfo->GetDerivFactors(GetPointsKeys());
1249
1250 for (i = 0; i < nquad1; i++)
1251 {
1252 Blas::Dscal(nquad0, 2.0 / (1 - z1[i]), &dEta_dXi[0][0] + i * nquad0, 1);
1253 Blas::Dscal(nquad0, 2.0 / (1 - z1[i]), &dEta_dXi[1][0] + i * nquad0, 1);
1254 }
1255 for (i = 0; i < nquad0; i++)
1256 {
1257 Blas::Dscal(nquad1, 0.5 * (1 + z0[i]), &dEta_dXi[1][0] + i, nquad0);
1258 }
1259
1260 Array<OneD, NekDouble> tmp(nqtot);
1261 if ((type == SpatialDomains::eRegular ||
1263 {
1264 Vmath::Smul(nqtot, df[0][0], &dEta_dXi[0][0], 1, &tmp[0], 1);
1265 Vmath::Svtvp(nqtot, df[1][0], &dEta_dXi[1][0], 1, &tmp[0], 1, &tmp[0],
1266 1);
1267
1268 Vmath::Vmul(nqtot, &tmp[0], 1, &tmp[0], 1,
1270 Vmath::Smul(nqtot, df[1][0], &tmp[0], 1,
1272
1273 Vmath::Smul(nqtot, df[2][0], &dEta_dXi[0][0], 1, &tmp[0], 1);
1274 Vmath::Svtvp(nqtot, df[3][0], &dEta_dXi[1][0], 1, &tmp[0], 1, &tmp[0],
1275 1);
1276
1277 Vmath::Vvtvp(nqtot, &tmp[0], 1, &tmp[0], 1,
1280 Vmath::Svtvp(nqtot, df[3][0], &tmp[0], 1,
1283
1284 if (GetCoordim() == 3)
1285 {
1286 Vmath::Smul(nqtot, df[4][0], &dEta_dXi[0][0], 1, &tmp[0], 1);
1287 Vmath::Svtvp(nqtot, df[5][0], &dEta_dXi[1][0], 1, &tmp[0], 1,
1288 &tmp[0], 1);
1289
1290 Vmath::Vvtvp(nqtot, &tmp[0], 1, &tmp[0], 1,
1293 Vmath::Svtvp(nqtot, df[5][0], &tmp[0], 1,
1296 }
1297
1298 NekDouble g2 = df[1][0] * df[1][0] + df[3][0] * df[3][0];
1299 if (GetCoordim() == 3)
1300 {
1301 g2 += df[5][0] * df[5][0];
1302 }
1303 Vmath::Fill(nqtot, g2, &m_metrics[eMetricLaplacian11][0], 1);
1304 }
1305 else
1306 {
1307
1308 Vmath::Vmul(nqtot, &df[0][0], 1, &dEta_dXi[0][0], 1, &tmp[0], 1);
1309 Vmath::Vvtvp(nqtot, &df[1][0], 1, &dEta_dXi[1][0], 1, &tmp[0], 1,
1310 &tmp[0], 1);
1311
1312 Vmath::Vmul(nqtot, &tmp[0], 1, &tmp[0], 1,
1314 Vmath::Vmul(nqtot, &df[1][0], 1, &tmp[0], 1,
1316 Vmath::Vmul(nqtot, &df[1][0], 1, &df[1][0], 1,
1318
1319 Vmath::Vmul(nqtot, &df[2][0], 1, &dEta_dXi[0][0], 1, &tmp[0], 1);
1320 Vmath::Vvtvp(nqtot, &df[3][0], 1, &dEta_dXi[1][0], 1, &tmp[0], 1,
1321 &tmp[0], 1);
1322
1323 Vmath::Vvtvp(nqtot, &tmp[0], 1, &tmp[0], 1,
1326 Vmath::Vvtvp(nqtot, &df[3][0], 1, &tmp[0], 1,
1329 Vmath::Vvtvp(nqtot, &df[3][0], 1, &df[3][0], 1,
1332
1333 if (GetCoordim() == 3)
1334 {
1335 Vmath::Vmul(nqtot, &df[4][0], 1, &dEta_dXi[0][0], 1, &tmp[0], 1);
1336 Vmath::Vvtvp(nqtot, &df[5][0], 1, &dEta_dXi[1][0], 1, &tmp[0], 1,
1337 &tmp[0], 1);
1338
1339 Vmath::Vvtvp(nqtot, &tmp[0], 1, &tmp[0], 1,
1342 Vmath::Vvtvp(nqtot, &df[5][0], 1, &tmp[0], 1,
1345 Vmath::Vvtvp(nqtot, &df[5][0], 1, &df[5][0], 1,
1348 }
1349 }
1350
1351 for (unsigned int i = 0; i < dim; ++i)
1352 {
1353 for (unsigned int j = i; j < dim; ++j)
1354 {
1356 }
1357 }
1358}
1359
1360/**
1361 * Function is used to compute exactly the advective numerical flux on
1362 * theinterface of two elements with different expansions, hence an
1363 * appropriate number of Gauss points has to be used. The number of
1364 * Gauss points has to be equal to the number used by the highest
1365 * polynomial degree of the two adjacent elements. Furthermore, this
1366 * function is used to compute the sensor value in each element.
1367 *
1368 * @param numMin Is the reduced polynomial order
1369 * @param inarray Input array of coefficients
1370 * @param dumpVar Output array of reduced coefficients.
1371 */
1373 const Array<OneD, const NekDouble> &inarray,
1374 Array<OneD, NekDouble> &outarray)
1375{
1376 int n_coeffs = inarray.size();
1377 int nquad0 = m_base[0]->GetNumPoints();
1378 int nquad1 = m_base[1]->GetNumPoints();
1379 int nqtot = nquad0 * nquad1;
1380 int nmodes0 = m_base[0]->GetNumModes();
1381 int nmodes1 = m_base[1]->GetNumModes();
1382 int numMin2 = nmodes0, i;
1383
1384 Array<OneD, NekDouble> coeff(n_coeffs, 0.0);
1385 Array<OneD, NekDouble> phys_tmp(nqtot, 0.0);
1386 Array<OneD, NekDouble> tmp, tmp2;
1387
1388 const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
1389 const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
1390
1391 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
1392 LibUtilities::BasisKey b1(m_base[1]->GetBasisType(), nmodes1, Pkey1);
1393 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
1394 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_B, nmodes1, Pkey1);
1395
1396 // Check if it is also possible to use the same InterCoeff routine
1397 // which is also used for Quadrilateral and Hexagonal shaped
1398 // elements
1399
1400 // For now, set up the used basis on the standard element to
1401 // calculate the phys values, set up the orthogonal basis to do a
1402 // forward transform, to obtain the coefficients in orthogonal
1403 // coefficient space
1404 StdRegions::StdTriExpSharedPtr m_OrthoTriExp;
1406
1409 bortho0, bortho1);
1410
1411 m_TriExp->BwdTrans(inarray, phys_tmp);
1412 m_OrthoTriExp->FwdTrans(phys_tmp, coeff);
1413
1414 for (i = 0; i < n_coeffs; i++)
1415 {
1416 if (i == numMin)
1417 {
1418 coeff[i] = 0.0;
1419 numMin += numMin2 - 1;
1420 numMin2 -= 1.0;
1421 }
1422 }
1423
1424 m_OrthoTriExp->BwdTrans(coeff, phys_tmp);
1425 m_TriExp->FwdTrans(phys_tmp, outarray);
1426}
1427
1429 const StdRegions::StdMatrixKey &mkey)
1430{
1431 int nq = GetTotPoints();
1432
1433 // Calculate sqrt of the Jacobian
1435 Array<OneD, NekDouble> sqrt_jac(nq);
1436 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1437 {
1438 Vmath::Vsqrt(nq, jac, 1, sqrt_jac, 1);
1439 }
1440 else
1441 {
1442 Vmath::Fill(nq, sqrt(jac[0]), sqrt_jac, 1);
1443 }
1444
1445 // Multiply array by sqrt(Jac)
1446 Vmath::Vmul(nq, sqrt_jac, 1, array, 1, array, 1);
1447
1448 // Apply std region filter
1449 StdTriExp::v_SVVLaplacianFilter(array, mkey);
1450
1451 // Divide by sqrt(Jac)
1452 Vmath::Vdiv(nq, array, 1, sqrt_jac, 1, array, 1);
1453}
1454
1455/** @brief: This method gets all of the factors which are
1456 required as part of the Gradient Jump Penalty
1457 stabilisation and involves the product of the normal and
1458 geometric factors along the element trace.
1459*/
1461 Array<OneD, Array<OneD, NekDouble>> &d0factors,
1462 Array<OneD, Array<OneD, NekDouble>> &d1factors,
1463 Array<OneD, Array<OneD, NekDouble>> &d2factors)
1464{
1465 boost::ignore_unused(d2factors); // for 3D shapes
1466 int nquad0 = GetNumPoints(0);
1467 int nquad1 = GetNumPoints(1);
1468
1470 m_metricinfo->GetDerivFactors(GetPointsKeys());
1471
1472 if (d0factors.size() != 3)
1473 {
1474 d0factors = Array<OneD, Array<OneD, NekDouble>>(3);
1475 d1factors = Array<OneD, Array<OneD, NekDouble>>(3);
1476 }
1477
1478 if (d0factors[0].size() != nquad0)
1479 {
1480 d0factors[0] = Array<OneD, NekDouble>(nquad0);
1481 d1factors[0] = Array<OneD, NekDouble>(nquad0);
1482 }
1483
1484 if (d0factors[1].size() != nquad1)
1485 {
1486 d0factors[1] = Array<OneD, NekDouble>(nquad1);
1487 d0factors[2] = Array<OneD, NekDouble>(nquad1);
1488 d1factors[1] = Array<OneD, NekDouble>(nquad1);
1489 d1factors[2] = Array<OneD, NekDouble>(nquad1);
1490 }
1491
1492 // Outwards normals
1494 GetTraceNormal(0);
1496 GetTraceNormal(1);
1498 GetTraceNormal(2);
1499
1500 int ncoords = normal_0.size();
1501
1502 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1503 {
1504
1505 // d xi_2/dx n_x
1506 for (int i = 0; i < nquad0; ++i)
1507 {
1508 d1factors[0][i] = df[1][i] * normal_0[0][i];
1509 }
1510
1511 // d xi_1/dx n_x
1512 for (int i = 0; i < nquad1; ++i)
1513 {
1514 d0factors[1][i] = df[0][(i + 1) * nquad0 - 1] * normal_1[0][i];
1515 d0factors[2][i] = df[0][i * nquad0] * normal_2[0][i];
1516 }
1517
1518 for (int n = 1; n < ncoords; ++n)
1519 {
1520 // d xi_2/dy n_y
1521 // needs checking for 3D coords
1522 for (int i = 0; i < nquad0; ++i)
1523 {
1524 d1factors[0][i] += df[2 * n + 1][i] * normal_0[n][i];
1525 }
1526
1527 // d xi_1/dy n_y
1528 // needs checking for 3D coords
1529 for (int i = 0; i < nquad1; ++i)
1530 {
1531 d0factors[1][i] +=
1532 df[2 * n][(i + 1) * nquad0 - 1] * normal_1[n][i];
1533 d0factors[2][i] += df[2 * n][i * nquad0] * normal_2[n][i];
1534 }
1535 }
1536
1537 // d0 factors
1538 // d xi_1/dx n_x
1539 for (int i = 0; i < nquad0; ++i)
1540 {
1541 d0factors[0][i] = df[0][i] * normal_0[0][i];
1542 }
1543
1544 // d xi_2/dx n_x
1545 for (int i = 0; i < nquad1; ++i)
1546 {
1547 d1factors[1][i] = df[1][(i + 1) * nquad0 - 1] * normal_1[0][i];
1548 d1factors[2][i] = df[1][i * nquad0] * normal_2[0][i];
1549 }
1550
1551 for (int n = 1; n < ncoords; ++n)
1552 {
1553 // d xi_1/dy n_y
1554 // needs checking for 3D coords
1555 for (int i = 0; i < nquad0; ++i)
1556 {
1557 d0factors[0][i] += df[2 * n][i] * normal_0[n][i];
1558 }
1559
1560 // d xi_2/dy n_y
1561 // needs checking for 3D coords
1562 for (int i = 0; i < nquad1; ++i)
1563 {
1564 d1factors[1][i] +=
1565 df[2 * n + 1][(i + 1) * nquad0 - 1] * normal_1[n][i];
1566 d1factors[2][i] += df[2 * n + 1][i * nquad0] * normal_2[n][i];
1567 }
1568 }
1569 }
1570 else
1571 {
1572 // d xi_2/dx n_x
1573 for (int i = 0; i < nquad0; ++i)
1574 {
1575 d1factors[0][i] = df[1][0] * normal_0[0][i];
1576 }
1577
1578 // d xi_1/dx n_x
1579 for (int i = 0; i < nquad1; ++i)
1580 {
1581 d0factors[1][i] = df[0][0] * normal_1[0][i];
1582 d0factors[2][i] = df[0][0] * normal_2[0][i];
1583 }
1584
1585 for (int n = 1; n < ncoords; ++n)
1586 {
1587 // d xi_2/dy n_y
1588 // needs checking for 3D coords
1589 for (int i = 0; i < nquad0; ++i)
1590 {
1591 d1factors[0][i] += df[2 * n + 1][0] * normal_0[n][i];
1592 }
1593
1594 // d xi_1/dy n_y
1595 // needs checking for 3D coords
1596 for (int i = 0; i < nquad1; ++i)
1597 {
1598 d0factors[1][i] += df[2 * n][0] * normal_1[n][i];
1599 d0factors[2][i] += df[2 * n][0] * normal_2[n][i];
1600 }
1601 }
1602
1603 // d1factors
1604 // d xi_1/dx n_x
1605 for (int i = 0; i < nquad0; ++i)
1606 {
1607 d0factors[0][i] = df[0][0] * normal_0[0][i];
1608 }
1609
1610 // d xi_2/dx n_x
1611 for (int i = 0; i < nquad1; ++i)
1612 {
1613 d1factors[1][i] = df[1][0] * normal_1[0][i];
1614 d1factors[2][i] = df[1][0] * normal_2[0][i];
1615 }
1616
1617 for (int n = 1; n < ncoords; ++n)
1618 {
1619 // d xi_1/dy n_y
1620 // needs checking for 3D coords
1621 for (int i = 0; i < nquad0; ++i)
1622 {
1623 d0factors[0][i] += df[2 * n][0] * normal_0[n][i];
1624 }
1625
1626 // d xi_2/dy n_y
1627 // needs checking for 3D coords
1628 for (int i = 0; i < nquad1; ++i)
1629 {
1630 d1factors[1][i] += df[2 * n + 1][0] * normal_1[n][i];
1631 d1factors[2][i] += df[2 * n + 1][0] * normal_2[n][i];
1632 }
1633 }
1634 }
1635}
1636} // namespace LocalRegions
1637} // 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
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:122
PointsKey GetPointsKey() const
Return distribution of points.
Definition: Basis.h:139
Defines a specification for a set of points.
Definition: Points.h:55
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
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
Definition: Expansion.cpp:608
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: Expansion.cpp:535
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_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: TriExp.cpp:436
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: TriExp.h:248
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const override
Definition: TriExp.cpp:645
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: TriExp.cpp:1153
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: TriExp.cpp:1131
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: TriExp.cpp:700
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: TriExp.cpp:402
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
Definition: TriExp.cpp:467
void v_DropLocMatrix(const MatrixKey &mkey) override
Definition: TriExp.cpp:1094
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: TriExp.cpp:409
virtual void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp) override
Definition: TriExp.cpp:1160
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: TriExp.cpp:680
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const MatrixKey &mkey) override
Definition: TriExp.cpp:1099
virtual void v_GetTracePhysMap(const int edge, Array< OneD, int > &outarray) override
Definition: TriExp.cpp:787
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdRegions::StdMatrixKey &mkey) override
Definition: TriExp.cpp:1079
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: TriExp.h:250
virtual void v_IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: TriExp.cpp:532
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: TriExp.cpp:188
virtual StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const override
Definition: TriExp.cpp:652
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: TriExp.cpp:1116
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
Definition: TriExp.cpp:663
virtual DNekScalMatSharedPtr v_GetLocMatrix(const MatrixKey &mkey) override
Definition: TriExp.cpp:1089
virtual StdRegions::Orientation v_GetTraceOrient(int edge) override
Definition: TriExp.cpp:1052
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
Definition: TriExp.cpp:692
virtual void v_ComputeTraceNormal(const int edge) override
Definition: TriExp.cpp:822
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: TriExp.cpp:721
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
Definition: TriExp.cpp:1057
void v_DropLocStaticCondMatrix(const MatrixKey &mkey) override
Definition: TriExp.cpp:1104
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[p]*base1[pq] and put into ou...
Definition: TriExp.cpp:396
virtual void v_GetTraceQFactors(const int edge, Array< OneD, NekDouble > &outarray) override
Definition: TriExp.cpp:780
virtual void v_WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: TriExp.cpp:1139
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: TriExp.cpp:269
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: TriExp.cpp:1372
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: TriExp.cpp:1012
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: TriExp.cpp:1109
virtual void v_NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors) override
: This method gets all of the factors which are required as part of the Gradient Jump Penalty stabili...
Definition: TriExp.cpp:1460
virtual void v_IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Directinoal Derivative in the modal space in the dir direction of varcoeffs.
Definition: TriExp.cpp:544
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey) override
Definition: TriExp.cpp:1428
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: TriExp.cpp:98
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: TriExp.cpp:610
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: TriExp.cpp:253
virtual void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
Definition: TriExp.cpp:1146
TriExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const SpatialDomains::Geometry2DSharedPtr &geom)
Constructor using BasisKey class for quadrature points and order definition.
Definition: TriExp.cpp:49
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrates the specified function over the domain.
Definition: TriExp.cpp:75
virtual void v_ComputeLaplacianMetric() override
Definition: TriExp.cpp:1216
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
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
const LibUtilities::BasisKey GetTraceBasisKey(const int i, int k=-1) const
This function returns the basis key belonging to the i-th trace.
Definition: StdExpansion.h:305
void 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
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:211
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:224
void IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:855
Array< OneD, LibUtilities::BasisSharedPtr > m_base
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 void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:151
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
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ 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< StdExpansion > StdExpansionSharedPtr
std::shared_ptr< StdTriExp > StdTriExpSharedPtr
Definition: StdTriExp.h:232
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