Nektar++
StdExpansion2D.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: StdExpansion2D.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: Daughter of StdExpansion. This class contains routine
32// which are common to 2D expansion. Typically this inolves physiocal
33// space operations.
34//
35///////////////////////////////////////////////////////////////////////////////
36
38
39#ifdef max
40#undef max
41#endif
42
43namespace Nektar::StdRegions
44{
45
47 const LibUtilities::BasisKey &Bb)
48 : StdExpansion(numcoeffs, 2, Ba, Bb)
49{
50}
51
52//----------------------------
53// Differentiation Methods
54//----------------------------
56 const Array<OneD, const NekDouble> &inarray,
57 Array<OneD, NekDouble> &outarray_d0, Array<OneD, NekDouble> &outarray_d1)
58{
59 int nquad0 = m_base[0]->GetNumPoints();
60 int nquad1 = m_base[1]->GetNumPoints();
61
62 if (outarray_d0.size() > 0) // calculate du/dx_0
63 {
64 DNekMatSharedPtr D0 = m_base[0]->GetD();
65 if (inarray.data() == outarray_d0.data())
66 {
67 Array<OneD, NekDouble> wsp(nquad0 * nquad1);
68 Vmath::Vcopy(nquad0 * nquad1, inarray.get(), 1, wsp.get(), 1);
69 Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
70 &(D0->GetPtr())[0], nquad0, &wsp[0], nquad0, 0.0,
71 &outarray_d0[0], nquad0);
72 }
73 else
74 {
75 Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
76 &(D0->GetPtr())[0], nquad0, &inarray[0], nquad0, 0.0,
77 &outarray_d0[0], nquad0);
78 }
79 }
80
81 if (outarray_d1.size() > 0) // calculate du/dx_1
82 {
83 DNekMatSharedPtr D1 = m_base[1]->GetD();
84 if (inarray.data() == outarray_d1.data())
85 {
86 Array<OneD, NekDouble> wsp(nquad0 * nquad1);
87 Vmath::Vcopy(nquad0 * nquad1, inarray.get(), 1, wsp.get(), 1);
88 Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &wsp[0], nquad0,
89 &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0],
90 nquad0);
91 }
92 else
93 {
94 Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &inarray[0],
95 nquad0, &(D1->GetPtr())[0], nquad1, 0.0,
96 &outarray_d1[0], nquad0);
97 }
98 }
99}
100
102 const Array<OneD, const NekDouble> &coords,
103 const Array<OneD, const NekDouble> &physvals)
104{
105 ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
106 ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
107 ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
108 ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
109
111 LocCoordToLocCollapsed(coords, coll);
112
113 const int nq0 = m_base[0]->GetNumPoints();
114 const int nq1 = m_base[1]->GetNumPoints();
115
116 Array<OneD, NekDouble> wsp(nq1);
117 for (int i = 0; i < nq1; ++i)
118 {
119 wsp[i] = StdExpansion::BaryEvaluate<0>(coll[0], &physvals[0] + i * nq0);
120 }
121
122 return StdExpansion::BaryEvaluate<1>(coll[1], &wsp[0]);
123}
124
127 const Array<OneD, const NekDouble> &physvals)
128{
129 NekDouble val;
130 int i;
131 int nq0 = m_base[0]->GetNumPoints();
132 int nq1 = m_base[1]->GetNumPoints();
133 Array<OneD, NekDouble> wsp1(nq1);
134
135 // interpolate first coordinate direction
136 for (i = 0; i < nq1; ++i)
137 {
138 wsp1[i] =
139 Blas::Ddot(nq0, &(I[0]->GetPtr())[0], 1, &physvals[i * nq0], 1);
140 }
141
142 // interpolate in second coordinate direction
143 val = Blas::Ddot(nq1, I[1]->GetPtr(), 1, wsp1, 1);
144
145 return val;
146}
147
149 [[maybe_unused]] const Array<OneD, NekDouble> &coord,
150 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
151 [[maybe_unused]] std::array<NekDouble, 3> &firstOrderDerivs)
152{
153 return 0;
154}
155
156//////////////////////////////
157// Integration Methods
158//////////////////////////////
159
163{
164 int i;
165 NekDouble Int = 0.0;
166 int nquad0 = m_base[0]->GetNumPoints();
167 int nquad1 = m_base[1]->GetNumPoints();
168 Array<OneD, NekDouble> tmp(nquad0 * nquad1);
169
170 // multiply by integration constants
171 for (i = 0; i < nquad1; ++i)
172 {
173 Vmath::Vmul(nquad0, &inarray[0] + i * nquad0, 1, w0.get(), 1,
174 &tmp[0] + i * nquad0, 1);
175 }
176
177 for (i = 0; i < nquad0; ++i)
178 {
179 Vmath::Vmul(nquad1, &tmp[0] + i, nquad0, w1.get(), 1, &tmp[0] + i,
180 nquad0);
181 }
182 Int = Vmath::Vsum(nquad0 * nquad1, tmp, 1);
183
184 return Int;
185}
186
188 const Array<OneD, const NekDouble> &base0,
189 const Array<OneD, const NekDouble> &base1,
190 const Array<OneD, const NekDouble> &inarray,
192 bool doCheckCollDir0, bool doCheckCollDir1)
193{
194 v_BwdTrans_SumFacKernel(base0, base1, inarray, outarray, wsp,
195 doCheckCollDir0, doCheckCollDir1);
196}
197
199 const Array<OneD, const NekDouble> &base0,
200 const Array<OneD, const NekDouble> &base1,
201 const Array<OneD, const NekDouble> &inarray,
203 bool doCheckCollDir0, bool doCheckCollDir1)
204{
205 v_IProductWRTBase_SumFacKernel(base0, base1, inarray, outarray, wsp,
206 doCheckCollDir0, doCheckCollDir1);
207}
208
210{
211 ASSERTL1((dir == 0) || (dir == 1), "Invalid direction.");
212
213 int nquad0 = m_base[0]->GetNumPoints();
214 int nquad1 = m_base[1]->GetNumPoints();
215 int nqtot = nquad0 * nquad1;
216 int nmodes0 = m_base[0]->GetNumModes();
217
218 Array<OneD, NekDouble> tmp1(2 * nqtot + m_ncoeffs + nmodes0 * nquad1, 0.0);
219 Array<OneD, NekDouble> tmp3(tmp1 + 2 * nqtot);
220 Array<OneD, NekDouble> tmp4(tmp1 + 2 * nqtot + m_ncoeffs);
221
222 switch (dir)
223 {
224 case 0:
225 for (int i = 0; i < nqtot; i++)
226 {
227 tmp1[i] = 1.0;
228 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
229 m_base[1]->GetBdata(), tmp1, tmp3,
230 tmp4, false, true);
231 tmp1[i] = 0.0;
232
233 for (int j = 0; j < m_ncoeffs; j++)
234 {
235 (*mat)(j, i) = tmp3[j];
236 }
237 }
238 break;
239 case 1:
240 for (int i = 0; i < nqtot; i++)
241 {
242 tmp1[i] = 1.0;
244 m_base[1]->GetDbdata(), tmp1, tmp3,
245 tmp4, true, false);
246 tmp1[i] = 0.0;
247
248 for (int j = 0; j < m_ncoeffs; j++)
249 {
250 (*mat)(j, i) = tmp3[j];
251 }
252 }
253 break;
254 default:
255 NEKERROR(ErrorUtil::efatal, "Not a 2D expansion.");
256 break;
257 }
258}
259
261 const Array<OneD, const NekDouble> &inarray,
263{
264 if (mkey.GetNVarCoeff() == 0 &&
267 {
268 using std::max;
269
270 // This implementation is only valid when there are no
271 // coefficients associated to the Laplacian operator
272 int nquad0 = m_base[0]->GetNumPoints();
273 int nquad1 = m_base[1]->GetNumPoints();
274 int nqtot = nquad0 * nquad1;
275 int nmodes0 = m_base[0]->GetNumModes();
276 int nmodes1 = m_base[1]->GetNumModes();
277 int wspsize =
278 max(max(max(nqtot, m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
279
280 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
281 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
282
283 // Allocate temporary storage
284 Array<OneD, NekDouble> wsp0(4 * wspsize); // size wspsize
285 Array<OneD, NekDouble> wsp1(wsp0 + wspsize); // size 3*wspsize
286
287 if (!(m_base[0]->Collocation() && m_base[1]->Collocation()))
288 {
289 // LAPLACIAN MATRIX OPERATION
290 // wsp0 = u = B * u_hat
291 // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
292 // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
293 BwdTrans_SumFacKernel(base0, base1, inarray, wsp0, wsp1, true,
294 true);
295 LaplacianMatrixOp_MatFree_Kernel(wsp0, outarray, wsp1);
296 }
297 else
298 {
299 LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp1);
300 }
301 }
302 else
303 {
305 mkey);
306 }
307}
308
310 const Array<OneD, const NekDouble> &inarray,
312{
313 if (mkey.GetNVarCoeff() == 0 &&
316 {
317 using std::max;
318
319 int nquad0 = m_base[0]->GetNumPoints();
320 int nquad1 = m_base[1]->GetNumPoints();
321 int nqtot = nquad0 * nquad1;
322 int nmodes0 = m_base[0]->GetNumModes();
323 int nmodes1 = m_base[1]->GetNumModes();
324 int wspsize =
325 max(max(max(nqtot, m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
327
328 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
329 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
330
331 // Allocate temporary storage
332 Array<OneD, NekDouble> wsp0(5 * wspsize); // size wspsize
333 Array<OneD, NekDouble> wsp1(wsp0 + wspsize); // size wspsize
334 Array<OneD, NekDouble> wsp2(wsp0 + 2 * wspsize); // size 3*wspsize
335
336 if (!(m_base[0]->Collocation() && m_base[1]->Collocation()))
337 {
338 // MASS MATRIX OPERATION
339 // The following is being calculated:
340 // wsp0 = B * u_hat = u
341 // wsp1 = W * wsp0
342 // outarray = B^T * wsp1 = B^T * W * B * u_hat = M * u_hat
343 BwdTrans_SumFacKernel(base0, base1, inarray, wsp0, wsp2, true,
344 true);
345 MultiplyByQuadratureMetric(wsp0, wsp1);
346 IProductWRTBase_SumFacKernel(base0, base1, wsp1, outarray, wsp2,
347 true, true);
348
349 LaplacianMatrixOp_MatFree_Kernel(wsp0, wsp1, wsp2);
350 }
351 else
352 {
353 MultiplyByQuadratureMetric(inarray, outarray);
354 LaplacianMatrixOp_MatFree_Kernel(inarray, wsp1, wsp2);
355 }
356
357 // outarray = lambda * outarray + wsp1
358 // = (lambda * M + L ) * u_hat
359 Vmath::Svtvp(m_ncoeffs, lambda, &outarray[0], 1, &wsp1[0], 1,
360 &outarray[0], 1);
361 }
362 else
363 {
365 mkey);
366 }
367}
368
370 [[maybe_unused]] const unsigned int traceid,
371 [[maybe_unused]] Array<OneD, unsigned int> &maparray)
372{
373 ASSERTL0(false,
374 "This method must be defined at the individual shape level");
375}
376
377/** \brief Determine the mapping to re-orientate the
378 coefficients along the element trace (assumed to align
379 with the standard element) into the orientation of the
380 local trace given by edgeOrient.
381 */
383 const unsigned int eid, Array<OneD, unsigned int> &maparray,
384 Array<OneD, int> &signarray, Orientation edgeOrient, int P,
385 [[maybe_unused]] int Q)
386{
387 unsigned int i;
388
389 int dir;
390 // determine basis direction for edge.
392 {
393 dir = (eid == 0) ? 0 : 1;
394 }
395 else
396 {
397 dir = eid % 2;
398 }
399
400 int numModes = m_base[dir]->GetNumModes();
401
402 // P is the desired length of the map
403 P = (P == -1) ? numModes : P;
404
405 // decalare maparray
406 if (maparray.size() != P)
407 {
408 maparray = Array<OneD, unsigned int>(P);
409 }
410
411 // fill default mapping as increasing index
412 for (i = 0; i < P; ++i)
413 {
414 maparray[i] = i;
415 }
416
417 if (signarray.size() != P)
418 {
419 signarray = Array<OneD, int>(P, 1);
420 }
421 else
422 {
423 std::fill(signarray.get(), signarray.get() + P, 1);
424 }
425
426 // Zero signmap and set maparray to zero if
427 // elemental modes are not as large as trace modes
428 for (i = numModes; i < P; ++i)
429 {
430 signarray[i] = 0.0;
431 maparray[i] = maparray[0];
432 }
433
434 if (edgeOrient == eBackwards)
435 {
436 const LibUtilities::BasisType bType = GetBasisType(dir);
437
438 if ((bType == LibUtilities::eModified_A) ||
439 (bType == LibUtilities::eModified_B))
440 {
441 std::swap(maparray[0], maparray[1]);
442
443 for (i = 3; i < std::min(P, numModes); i += 2)
444 {
445 signarray[i] *= -1;
446 }
447 }
448 else if (bType == LibUtilities::eGLL_Lagrange ||
450 {
451 ASSERTL1(P == numModes, "Different trace space edge dimension "
452 "and element edge dimension not currently "
453 "possible for GLL-Lagrange bases");
454
455 std::reverse(maparray.get(), maparray.get() + P);
456 }
457 else
458 {
459 ASSERTL0(false, "Mapping not defined for this type of basis");
460 }
461 }
462}
463
466 Array<OneD, int> &signarray,
467 Orientation edgeOrient, int P,
468 int Q)
469{
470 Array<OneD, unsigned int> map1, map2;
471 v_GetTraceCoeffMap(eid, map1);
472 v_GetElmtTraceToTraceMap(eid, map2, signarray, edgeOrient, P, Q);
473
474 if (maparray.size() != map2.size())
475 {
476 maparray = Array<OneD, unsigned int>(map2.size());
477 }
478
479 for (int i = 0; i < map2.size(); ++i)
480 {
481 maparray[i] = map1[map2[i]];
482 }
483}
484
485} // namespace Nektar::StdRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
Describes the specification for a Basis.
Definition: Basis.h:45
void v_GenStdMatBwdDeriv(const int dir, DNekMatSharedPtr &mat) override
virtual void v_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, bool doCheckCollDir1)=0
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points.
void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
void v_GetTraceToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient=eForwards, int P=-1, int Q=-1) override
virtual void v_BwdTrans_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, bool doCheckCollDir1)=0
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
void BwdTrans_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)
void v_GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
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)
NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
This function evaluates the expansion at a single (arbitrary) point of the domain.
void v_GetElmtTraceToTraceMap(const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
Determine the mapping to re-orientate the coefficients along the element trace (assumed to align with...
The base class for all shapes.
Definition: StdExpansion.h:65
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
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
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:124
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:133
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:163
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:383
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ P
Monomial polynomials .
Definition: BasisType.h:62
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:57
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
static const NekDouble kNekZeroTol
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
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 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.hpp:396
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.hpp:608
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825