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 [[maybe_unused]] int numcoeffs,
48 [[maybe_unused]] const LibUtilities::BasisKey &Ba,
49 [[maybe_unused]] const LibUtilities::BasisKey &Bb)
50{
51}
52
53//----------------------------
54// Differentiation Methods
55//----------------------------
57 const Array<OneD, const NekDouble> &inarray,
58 Array<OneD, NekDouble> &outarray_d0, Array<OneD, NekDouble> &outarray_d1)
59{
60 int nquad0 = m_base[0]->GetNumPoints();
61 int nquad1 = m_base[1]->GetNumPoints();
62
63 if (outarray_d0.size() > 0) // calculate du/dx_0
64 {
65 DNekMatSharedPtr D0 = m_base[0]->GetD();
66 if (inarray.data() == outarray_d0.data())
67 {
68 Array<OneD, NekDouble> wsp(nquad0 * nquad1);
69 Vmath::Vcopy(nquad0 * nquad1, inarray.get(), 1, wsp.get(), 1);
70 Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
71 &(D0->GetPtr())[0], nquad0, &wsp[0], nquad0, 0.0,
72 &outarray_d0[0], nquad0);
73 }
74 else
75 {
76 Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
77 &(D0->GetPtr())[0], nquad0, &inarray[0], nquad0, 0.0,
78 &outarray_d0[0], nquad0);
79 }
80 }
81
82 if (outarray_d1.size() > 0) // calculate du/dx_1
83 {
84 DNekMatSharedPtr D1 = m_base[1]->GetD();
85 if (inarray.data() == outarray_d1.data())
86 {
87 Array<OneD, NekDouble> wsp(nquad0 * nquad1);
88 Vmath::Vcopy(nquad0 * nquad1, inarray.get(), 1, wsp.get(), 1);
89 Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &wsp[0], nquad0,
90 &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0],
91 nquad0);
92 }
93 else
94 {
95 Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &inarray[0],
96 nquad0, &(D1->GetPtr())[0], nquad1, 0.0,
97 &outarray_d1[0], nquad0);
98 }
99 }
100}
101
103 const Array<OneD, const NekDouble> &coords,
104 const Array<OneD, const NekDouble> &physvals)
105{
106 ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
107 ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
108 ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
109 ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
110
112 LocCoordToLocCollapsed(coords, coll);
113
114 const int nq0 = m_base[0]->GetNumPoints();
115 const int nq1 = m_base[1]->GetNumPoints();
116
117 Array<OneD, NekDouble> wsp(nq1);
118 for (int i = 0; i < nq1; ++i)
119 {
120 wsp[i] = StdExpansion::BaryEvaluate<0>(coll[0], &physvals[0] + i * nq0);
121 }
122
123 return StdExpansion::BaryEvaluate<1>(coll[1], &wsp[0]);
124}
125
128 const Array<OneD, const NekDouble> &physvals)
129{
130 NekDouble val;
131 int i;
132 int nq0 = m_base[0]->GetNumPoints();
133 int nq1 = m_base[1]->GetNumPoints();
134 Array<OneD, NekDouble> wsp1(nq1);
135
136 // interpolate first coordinate direction
137 for (i = 0; i < nq1; ++i)
138 {
139 wsp1[i] =
140 Blas::Ddot(nq0, &(I[0]->GetPtr())[0], 1, &physvals[i * nq0], 1);
141 }
142
143 // interpolate in second coordinate direction
144 val = Blas::Ddot(nq1, I[1]->GetPtr(), 1, wsp1, 1);
145
146 return val;
147}
148
150 [[maybe_unused]] const Array<OneD, NekDouble> &coord,
151 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
152 [[maybe_unused]] std::array<NekDouble, 3> &firstOrderDerivs)
153{
154 return 0;
155}
156
157//////////////////////////////
158// Integration Methods
159//////////////////////////////
160
164{
165 int i;
166 NekDouble Int = 0.0;
167 int nquad0 = m_base[0]->GetNumPoints();
168 int nquad1 = m_base[1]->GetNumPoints();
169 Array<OneD, NekDouble> tmp(nquad0 * nquad1);
170
171 // multiply by integration constants
172 for (i = 0; i < nquad1; ++i)
173 {
174 Vmath::Vmul(nquad0, &inarray[0] + i * nquad0, 1, w0.get(), 1,
175 &tmp[0] + i * nquad0, 1);
176 }
177
178 for (i = 0; i < nquad0; ++i)
179 {
180 Vmath::Vmul(nquad1, &tmp[0] + i, nquad0, w1.get(), 1, &tmp[0] + i,
181 nquad0);
182 }
183 Int = Vmath::Vsum(nquad0 * nquad1, tmp, 1);
184
185 return Int;
186}
187
189 const Array<OneD, const NekDouble> &base0,
190 const Array<OneD, const NekDouble> &base1,
191 const Array<OneD, const NekDouble> &inarray,
193 bool doCheckCollDir0, bool doCheckCollDir1)
194{
195 v_BwdTrans_SumFacKernel(base0, base1, inarray, outarray, wsp,
196 doCheckCollDir0, doCheckCollDir1);
197}
198
200 const Array<OneD, const NekDouble> &base0,
201 const Array<OneD, const NekDouble> &base1,
202 const Array<OneD, const NekDouble> &inarray,
204 bool doCheckCollDir0, bool doCheckCollDir1)
205{
206 v_IProductWRTBase_SumFacKernel(base0, base1, inarray, outarray, wsp,
207 doCheckCollDir0, doCheckCollDir1);
208}
209
211{
212 ASSERTL1((dir == 0) || (dir == 1), "Invalid direction.");
213
214 int nquad0 = m_base[0]->GetNumPoints();
215 int nquad1 = m_base[1]->GetNumPoints();
216 int nqtot = nquad0 * nquad1;
217 int nmodes0 = m_base[0]->GetNumModes();
218
219 Array<OneD, NekDouble> tmp1(2 * nqtot + m_ncoeffs + nmodes0 * nquad1, 0.0);
220 Array<OneD, NekDouble> tmp3(tmp1 + 2 * nqtot);
221 Array<OneD, NekDouble> tmp4(tmp1 + 2 * nqtot + m_ncoeffs);
222
223 switch (dir)
224 {
225 case 0:
226 for (int i = 0; i < nqtot; i++)
227 {
228 tmp1[i] = 1.0;
229 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
230 m_base[1]->GetBdata(), tmp1, tmp3,
231 tmp4, false, true);
232 tmp1[i] = 0.0;
233
234 for (int j = 0; j < m_ncoeffs; j++)
235 {
236 (*mat)(j, i) = tmp3[j];
237 }
238 }
239 break;
240 case 1:
241 for (int i = 0; i < nqtot; i++)
242 {
243 tmp1[i] = 1.0;
245 m_base[1]->GetDbdata(), tmp1, tmp3,
246 tmp4, true, false);
247 tmp1[i] = 0.0;
248
249 for (int j = 0; j < m_ncoeffs; j++)
250 {
251 (*mat)(j, i) = tmp3[j];
252 }
253 }
254 break;
255 default:
256 NEKERROR(ErrorUtil::efatal, "Not a 2D expansion.");
257 break;
258 }
259}
260
262 const Array<OneD, const NekDouble> &inarray,
264{
265 if (mkey.GetNVarCoeff() == 0 &&
268 {
269 using std::max;
270
271 // This implementation is only valid when there are no
272 // coefficients associated to the Laplacian operator
273 int nquad0 = m_base[0]->GetNumPoints();
274 int nquad1 = m_base[1]->GetNumPoints();
275 int nqtot = nquad0 * nquad1;
276 int nmodes0 = m_base[0]->GetNumModes();
277 int nmodes1 = m_base[1]->GetNumModes();
278 int wspsize =
279 max(max(max(nqtot, m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
280
281 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
282 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
283
284 // Allocate temporary storage
285 Array<OneD, NekDouble> wsp0(4 * wspsize); // size wspsize
286 Array<OneD, NekDouble> wsp1(wsp0 + wspsize); // size 3*wspsize
287
288 if (!(m_base[0]->Collocation() && m_base[1]->Collocation()))
289 {
290 // LAPLACIAN MATRIX OPERATION
291 // wsp0 = u = B * u_hat
292 // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
293 // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
294 BwdTrans_SumFacKernel(base0, base1, inarray, wsp0, wsp1, true,
295 true);
296 LaplacianMatrixOp_MatFree_Kernel(wsp0, outarray, wsp1);
297 }
298 else
299 {
300 LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp1);
301 }
302 }
303 else
304 {
306 mkey);
307 }
308}
309
311 const Array<OneD, const NekDouble> &inarray,
313{
314 if (mkey.GetNVarCoeff() == 0 &&
317 {
318 using std::max;
319
320 int nquad0 = m_base[0]->GetNumPoints();
321 int nquad1 = m_base[1]->GetNumPoints();
322 int nqtot = nquad0 * nquad1;
323 int nmodes0 = m_base[0]->GetNumModes();
324 int nmodes1 = m_base[1]->GetNumModes();
325 int wspsize =
326 max(max(max(nqtot, m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
328
329 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
330 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
331
332 // Allocate temporary storage
333 Array<OneD, NekDouble> wsp0(5 * wspsize); // size wspsize
334 Array<OneD, NekDouble> wsp1(wsp0 + wspsize); // size wspsize
335 Array<OneD, NekDouble> wsp2(wsp0 + 2 * wspsize); // size 3*wspsize
336
337 if (!(m_base[0]->Collocation() && m_base[1]->Collocation()))
338 {
339 // MASS MATRIX OPERATION
340 // The following is being calculated:
341 // wsp0 = B * u_hat = u
342 // wsp1 = W * wsp0
343 // outarray = B^T * wsp1 = B^T * W * B * u_hat = M * u_hat
344 BwdTrans_SumFacKernel(base0, base1, inarray, wsp0, wsp2, true,
345 true);
346 MultiplyByQuadratureMetric(wsp0, wsp1);
347 IProductWRTBase_SumFacKernel(base0, base1, wsp1, outarray, wsp2,
348 true, true);
349
350 LaplacianMatrixOp_MatFree_Kernel(wsp0, wsp1, wsp2);
351 }
352 else
353 {
354 MultiplyByQuadratureMetric(inarray, outarray);
355 LaplacianMatrixOp_MatFree_Kernel(inarray, wsp1, wsp2);
356 }
357
358 // outarray = lambda * outarray + wsp1
359 // = (lambda * M + L ) * u_hat
360 Vmath::Svtvp(m_ncoeffs, lambda, &outarray[0], 1, &wsp1[0], 1,
361 &outarray[0], 1);
362 }
363 else
364 {
366 mkey);
367 }
368}
369
371 [[maybe_unused]] const unsigned int traceid,
372 [[maybe_unused]] Array<OneD, unsigned int> &maparray)
373{
374 ASSERTL0(false,
375 "This method must be defined at the individual shape level");
376}
377
378/** \brief Determine the mapping to re-orientate the
379 coefficients along the element trace (assumed to align
380 with the standard element) into the orientation of the
381 local trace given by edgeOrient.
382 */
384 const unsigned int eid, Array<OneD, unsigned int> &maparray,
385 Array<OneD, int> &signarray, Orientation edgeOrient, int P,
386 [[maybe_unused]] int Q)
387{
388 unsigned int i;
389
390 int dir;
391 // determine basis direction for edge.
393 {
394 dir = (eid == 0) ? 0 : 1;
395 }
396 else
397 {
398 dir = eid % 2;
399 }
400
401 int numModes = m_base[dir]->GetNumModes();
402
403 // P is the desired length of the map
404 P = (P == -1) ? numModes : P;
405
406 // decalare maparray
407 if (maparray.size() != P)
408 {
409 maparray = Array<OneD, unsigned int>(P);
410 }
411
412 // fill default mapping as increasing index
413 for (i = 0; i < P; ++i)
414 {
415 maparray[i] = i;
416 }
417
418 if (signarray.size() != P)
419 {
420 signarray = Array<OneD, int>(P, 1);
421 }
422 else
423 {
424 std::fill(signarray.get(), signarray.get() + P, 1);
425 }
426
427 // Zero signmap and set maparray to zero if
428 // elemental modes are not as large as trace modes
429 for (i = numModes; i < P; ++i)
430 {
431 signarray[i] = 0.0;
432 maparray[i] = maparray[0];
433 }
434
435 if (edgeOrient == eBackwards)
436 {
437 const LibUtilities::BasisType bType = GetBasisType(dir);
438
439 if ((bType == LibUtilities::eModified_A) ||
440 (bType == LibUtilities::eModified_B))
441 {
442 std::swap(maparray[0], maparray[1]);
443
444 for (i = 3; i < std::min(P, numModes); i += 2)
445 {
446 signarray[i] *= -1;
447 }
448 }
449 else if (bType == LibUtilities::eGLL_Lagrange ||
451 {
452 ASSERTL1(P == numModes, "Different trace space edge dimension "
453 "and element edge dimension not currently "
454 "possible for GLL-Lagrange bases");
455
456 std::reverse(maparray.get(), maparray.get() + P);
457 }
458 else
459 {
460 ASSERTL0(false, "Mapping not defined for this type of basis");
461 }
462 }
463}
464
467 Array<OneD, int> &signarray,
468 Orientation edgeOrient, int P,
469 int Q)
470{
471 Array<OneD, unsigned int> map1, map2;
472 v_GetTraceCoeffMap(eid, map1);
473 v_GetElmtTraceToTraceMap(eid, map2, signarray, edgeOrient, P, Q);
474
475 if (maparray.size() != map2.size())
476 {
477 maparray = Array<OneD, unsigned int>(map2.size());
478 }
479
480 for (int i = 0; i < map2.size(); ++i)
481 {
482 maparray[i] = map1[map2[i]];
483 }
484}
485
486} // 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...
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