Nektar++
Foundations/Basis.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Basis.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: Basis definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
38
40{
43
44bool operator<(const BasisKey &lhs, const BasisKey &rhs)
45{
46 PointsKey lhsPointsKey = lhs.GetPointsKey();
47 PointsKey rhsPointsKey = rhs.GetPointsKey();
48
49 if (lhsPointsKey < rhsPointsKey)
50 {
51 return true;
52 }
53 if (lhsPointsKey != rhsPointsKey)
54 {
55 return false;
56 }
57
58 if (lhs.m_nummodes < rhs.m_nummodes)
59 {
60 return true;
61 }
62 if (lhs.m_nummodes > rhs.m_nummodes)
63 {
64 return false;
65 }
66
67 return (lhs.m_basistype < rhs.m_basistype);
68}
69
70bool operator>(const BasisKey &lhs, const BasisKey &rhs)
71{
72 return (rhs < lhs);
73}
74
76 const BasisKey &rhs) const
77{
78 return (lhs.m_basistype < rhs.m_basistype);
79}
80
81std::ostream &operator<<(std::ostream &os, const BasisKey &rhs)
82{
83 os << "NumModes: " << rhs.GetNumModes()
84 << " BasisType: " << BasisTypeMap[rhs.GetBasisType()];
85 os << " " << rhs.GetPointsKey() << std::endl;
86
87 return os;
88}
89
91 : m_basisKey(bkey), m_points(PointsManager()[bkey.GetPointsKey()]),
92 m_bdata(bkey.GetTotNumModes() * bkey.GetTotNumPoints()),
93 m_dbdata(bkey.GetTotNumModes() * bkey.GetTotNumPoints())
94{
95 m_InterpManager.RegisterGlobalCreator(
96 std::bind(&Basis::CalculateInterpMatrix, this, std::placeholders::_1));
97}
98
99std::shared_ptr<Basis> Basis::Create(const BasisKey &bkey)
100{
101 std::shared_ptr<Basis> returnval(new Basis(bkey));
102 returnval->Initialize();
103
104 return returnval;
105}
106
108{
109 ASSERTL0(GetNumModes() > 0,
110 "Cannot call Basis initialisation with zero or negative order");
111 ASSERTL0(GetTotNumPoints() > 0, "Cannot call Basis initialisation with "
112 "zero or negative numbers of points");
113
114 GenBasis();
115}
116
117/** \brief Calculate the interpolation Matrix for coefficient from
118 * one base (m_basisKey) to another (tbasis0)
119 */
120std::shared_ptr<NekMatrix<NekDouble>> Basis::CalculateInterpMatrix(
121 const BasisKey &tbasis0)
122{
123 size_t dim = m_basisKey.GetNumModes();
125 BasisKey fbkey(m_basisKey.GetBasisType(), dim, pkey);
126 BasisKey tbkey(tbasis0.GetBasisType(), dim, pkey);
127
128 // "Constructur" of the basis
129 BasisSharedPtr fbasis = BasisManager()[fbkey];
130 BasisSharedPtr tbasis = BasisManager()[tbkey];
131
132 // Get B Matrices
133 Array<OneD, NekDouble> fB_data = fbasis->GetBdata();
134 Array<OneD, NekDouble> tB_data = tbasis->GetBdata();
135
136 // Convert to a NekMatrix
137 NekMatrix<NekDouble> fB(dim, dim, fB_data);
138 NekMatrix<NekDouble> tB(dim, dim, tB_data);
139
140 // Invert the "to" matrix: tu = tB^(-1)*fB fu = ftB fu
141 tB.Invert();
142
143 // Compute transformation matrix
144 Array<OneD, NekDouble> zero1D(dim * dim, 0.0);
145 std::shared_ptr<NekMatrix<NekDouble>> ftB(
146 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(dim, dim,
147 zero1D));
148 (*ftB) = tB * fB;
149
150 return ftB;
151}
152
153// Method used to generate appropriate basis
154/** The following expansions are generated depending on the
155 * enum type defined in \a m_basisKey.m_basistype:
156 *
157 * NOTE: This definition does not follow the order in the
158 * Karniadakis \& Sherwin book since this leads to a more
159 * compact hierarchical pattern for implementation
160 * purposes. The order of these modes dictates the
161 * ordering of the expansion coefficients.
162 *
163 * In the following m_numModes = P
164 *
165 * \a eModified_A:
166 *
167 * m_bdata[i + j*m_numpoints] =
168 * \f$ \phi^a_i(z_j) = \left \{
169 * \begin{array}{ll} \left ( \frac{1-z_j}{2}\right ) & i = 0 \\
170 * \\
171 * \left ( \frac{1+z_j}{2}\right ) & i = 1 \\
172 * \\
173 * \left ( \frac{1-z_j}{2}\right )\left ( \frac{1+z_j}{2}\right )
174 * P^{1,1}_{i-2}(z_j) & 2\leq i < P\\
175 * \end{array} \right . \f$
176 *
177 * \a eModified_B:
178 *
179 * m_bdata[n(i,j) + k*m_numpoints] =
180 * \f$ \phi^b_{ij}(z_k) = \left \{ \begin{array}{lll}
181 * \phi^a_j(z_k) & i = 0, & 0\leq j < P \\
182 * \\
183 * \left ( \frac{1-z_k}{2}\right )^{i} & 1 \leq i < P,& j = 0 \\
184 * \\
185 * \left ( \frac{1-z_k}{2}\right )^{i} \left ( \frac{1+z_k}{2}\right )
186 * P^{2i-1,1}_{j-1}(z_k) & 1 \leq i < P,\ & 1\leq j < P-i\ \\
187 * \end{array} \right . , \f$
188 *
189 * where \f$ n(i,j) \f$ is a consecutive ordering of the
190 * triangular indices \f$ 0 \leq i, i+j < P \f$ where \a j
191 * runs fastest.
192 *
193 *
194 * \a eModified_C:
195 *
196 * m_bdata[n(i,j,k) + l*m_numpoints] =
197 * \f$ \phi^c_{ij,k}(z_l) = \phi^b_{i+j,k}(z_l) =
198 * \left \{ \begin{array}{llll}
199 * \phi^b_{j,k}(z_l) & i = 0, & 0\leq j < P & 0\leq k < P-j\\
200 * \\
201 * \left ( \frac{1-z_l}{2}\right )^{i+j} & 1\leq i < P,\
202 * & 0\leq j < P-i,\ & k = 0 \\
203 * \\
204 * \left ( \frac{1-z_l}{2}\right )^{i+j}
205 * \left ( \frac{1+z_l}{2}\right )
206 * P^{2i+2j-1,1}_{k-1}(z_k) & 1\leq i < P,& 0\leq j < P-i&
207 * 1\leq k < P-i-j \\
208 * \\
209 * \end{array} \right . , \f$
210 *
211 * where \f$ n(i,j,k) \f$ is a consecutive ordering of the
212 * triangular indices \f$ 0 \leq i, i+j, i+j+k < P \f$ where \a k
213 * runs fastest, then \a j and finally \a i.
214 *
215 */
217{
218 size_t i, p, q;
219 NekDouble scal;
220 Array<OneD, NekDouble> modeSharedArray;
221 NekDouble *mode;
224
225 m_points->GetZW(z, w);
226
227 const NekDouble *D = &(m_points->GetD()->GetPtr())[0];
228 size_t numModes = GetNumModes();
229 size_t numPoints = GetNumPoints();
230
231 switch (GetBasisType())
232 {
233
234 /** \brief Orthogonal basis A
235
236 \f$\tilde \psi_p^a (\eta_1) = L_p(\eta_1) = P_p^{0,0}(\eta_1)\f$
237
238 */
239 case eOrtho_A:
240 case eLegendre:
241 {
242 mode = m_bdata.data();
243
244 for (p = 0; p < numModes; ++p, mode += numPoints)
245 {
246 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, p, 0.0,
247 0.0);
248 // normalise
249 scal = sqrt(0.5 * (2.0 * p + 1.0));
250 for (i = 0; i < numPoints; ++i)
251 {
252 mode[i] *= scal;
253 }
254 }
255
256 // Define derivative basis
257 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
258 numPoints, m_bdata.data(), numPoints, 0.0,
259 m_dbdata.data(), numPoints);
260 } // end scope
261 break;
262
263 /** \brief Orthogonal basis B
264
265 \f$\tilde \psi_{pq}^b(\eta_2) = \left ( {1 - \eta_2} \over 2
266 \right)^p
267 P_q^{2p+1,0}(\eta_2)\f$ \\
268
269 */
270
271 // This is tilde psi_pq in Spencer's book, page 105
272 // The 3-dimensional array is laid out in memory such that
273 // 1) Eta_y values are the changing the fastest, then q and p.
274 // 2) q index increases by the stride of numPoints.
275 case eOrtho_B:
276 {
277 mode = m_bdata.data();
278
279 for (size_t p = 0; p < numModes; ++p)
280 {
281 for (size_t q = 0; q < numModes - p; ++q, mode += numPoints)
282 {
283 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, q,
284 2 * p + 1.0, 0.0);
285 for (size_t j = 0; j < numPoints; ++j)
286 {
287 mode[j] *=
288 sqrt(p + q + 1.0) * pow(0.5 * (1.0 - z[j]), p);
289 }
290 }
291 }
292
293 // Define derivative basis
294 Blas::Dgemm('n', 'n', numPoints, numModes * (numModes + 1) / 2,
295 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
296 0.0, m_dbdata.data(), numPoints);
297 }
298 break;
299
300 /** \brief Orthogonal basis C
301
302 \f$\tilde \psi_{pqr}^c = \left ( {1 - \eta_3} \over 2 \right)^{p+q}
303 P_r^{2p+2q+2, 0}(\eta_3)\f$ \ \
304
305 */
306
307 // This is tilde psi_pqr in Spencer's book, page 105
308 // The 4-dimensional array is laid out in memory such that
309 // 1) Eta_z values are the changing the fastest, then r, q, and finally
310 // p. 2) r index increases by the stride of numPoints.
311 case eOrtho_C:
312 {
313 size_t P = numModes - 1, Q = numModes - 1, R = numModes - 1;
314 mode = m_bdata.data();
315
316 for (size_t p = 0; p <= P; ++p)
317 {
318 for (size_t q = 0; q <= Q - p; ++q)
319 {
320 for (size_t r = 0; r <= R - p - q; ++r, mode += numPoints)
321 {
322 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, r,
323 2 * p + 2 * q + 2.0, 0.0);
324 for (size_t k = 0; k < numPoints; ++k)
325 {
326 // Note factor of 0.5 is part of normalisation
327 mode[k] *= pow(0.5 * (1.0 - z[k]), p + q);
328
329 // finish normalisation
330 mode[k] *= sqrt(r + p + q + 1.5);
331 }
332 }
333 }
334 }
335
336 // Define derivative basis
337 Blas::Dgemm('n', 'n', numPoints,
338 numModes * (numModes + 1) * (numModes + 2) / 6,
339 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
340 0.0, m_dbdata.data(), numPoints);
341 }
342 break;
343
344 /** \brief Orthogonal basis C for Pyramid expansion
345 (which is richer than tets)
346
347 \f$\tilde \psi_{pqr}^c = \left ( {1 - \eta_3} \over
348 2\right)^{pq} P_r^{2pq+2, 0}(\eta_3)\f$ \f$ \mbox{where
349 }pq = max(p+q,0) \f$
350
351 This orthogonal expansion has modes that are
352 always in the Cartesian space, however the
353 equivalent ModifiedPyr_C has vertex modes that do
354 not lie in this space. If one chooses \f$pq =
355 max(p+q-1,0)\f$ then the expansion will space the
356 same space as the vertices but the order of the
357 expanion in 'r' is reduced by one.
358
359 1) Eta_z values are the changing the fastest, then
360 r, q, and finally p. 2) r index increases by the
361 stride of numPoints.
362 */
363 case eOrthoPyr_C:
364 {
365 size_t P = numModes - 1, Q = numModes - 1, R = numModes - 1;
366 mode = m_bdata.data();
367
368 for (size_t p = 0; p <= P; ++p)
369 {
370 for (size_t q = 0; q <= Q; ++q)
371 {
372 for (size_t r = 0; r <= R - std::max(p, q);
373 ++r, mode += numPoints)
374 {
375 // this offset allows for orthogonal
376 // expansion to span linear FE space
377 // of modified basis but means that
378 // the cartesian polynomial space
379 // spanned by the expansion is one
380 // order lower.
381 // size_t pq = max(p + q -1,0);
382 size_t pq = std::max(p + q, size_t(0));
383
384 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, r,
385 2 * pq + 2.0, 0.0);
386 for (size_t k = 0; k < numPoints; ++k)
387 {
388 // Note factor of 0.5 is part of normalisation
389 mode[k] *= pow(0.5 * (1.0 - z[k]), pq);
390
391 // finish normalisation
392 mode[k] *= sqrt(r + pq + 1.5);
393 }
394 }
395 }
396 }
397
398 // Define derivative basis
399 Blas::Dgemm('n', 'n', numPoints,
400 numModes * (numModes + 1) * (numModes + 2) / 6,
401 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
402 0.0, m_dbdata.data(), numPoints);
403 }
404 break;
405
406 case eModified_A:
407 {
408 // Note the following packing deviates from the
409 // definition in the Book by Karniadakis in that we
410 // put the vertex degrees of freedom at the lower
411 // index range to follow a more hierarchic structure.
412
413 for (i = 0; i < numPoints; ++i)
414 {
415 m_bdata[i] = 0.5 * (1 - z[i]);
416 m_bdata[numPoints + i] = 0.5 * (1 + z[i]);
417 }
418
419 mode = m_bdata.data() + 2 * numPoints;
420
421 for (p = 2; p < numModes; ++p, mode += numPoints)
422 {
423 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, p - 2, 1.0,
424 1.0);
425
426 for (i = 0; i < numPoints; ++i)
427 {
428 mode[i] *= m_bdata[i] * m_bdata[numPoints + i];
429 }
430 }
431
432 // Define derivative basis
433 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
434 numPoints, m_bdata.data(), numPoints, 0.0,
435 m_dbdata.data(), numPoints);
436 }
437 break;
438
439 case eModified_B:
440 {
441 // Note the following packing deviates from the
442 // definition in the Book by Karniadakis in two
443 // ways. 1) We put the vertex degrees of freedom
444 // at the lower index range to follow a more
445 // hierarchic structure. 2) We do not duplicate
446 // the singular vertex definition so that only a
447 // triangular number (i.e. (modes)*(modes+1)/2) of
448 // modes are required consistent with the
449 // orthogonal basis.
450
451 // In the current structure the q index runs
452 // faster than the p index so that the matrix has
453 // a more compact structure
454
455 const NekDouble *one_m_z_pow, *one_p_z;
456
457 // bdata should be of size order*(order+1)/2*zorder
458
459 // first fow
460 for (i = 0; i < numPoints; ++i)
461 {
462 m_bdata[0 * numPoints + i] = 0.5 * (1 - z[i]);
463 m_bdata[1 * numPoints + i] = 0.5 * (1 + z[i]);
464 }
465
466 mode = m_bdata.data() + 2 * numPoints;
467
468 for (q = 2; q < numModes; ++q, mode += numPoints)
469 {
470 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, q - 2, 1.0,
471 1.0);
472
473 for (i = 0; i < numPoints; ++i)
474 {
475 mode[i] *= m_bdata[i] * m_bdata[numPoints + i];
476 }
477 }
478
479 // second row
480 for (i = 0; i < numPoints; ++i)
481 {
482 mode[i] = 0.5 * (1 - z[i]);
483 }
484
485 mode += numPoints;
486
487 for (q = 2; q < numModes; ++q, mode += numPoints)
488 {
489 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, q - 2, 1.0,
490 1.0);
491
492 for (i = 0; i < numPoints; ++i)
493 {
494 mode[i] *= m_bdata[i] * m_bdata[numPoints + i];
495 }
496 }
497
498 // third and higher rows
499 one_m_z_pow = m_bdata.data();
500 one_p_z = m_bdata.data() + numPoints;
501
502 for (p = 2; p < numModes; ++p)
503 {
504 for (i = 0; i < numPoints; ++i)
505 {
506 mode[i] = m_bdata[i] * one_m_z_pow[i];
507 }
508
509 one_m_z_pow = mode;
510 mode += numPoints;
511
512 for (q = 1; q < numModes - p; ++q, mode += numPoints)
513 {
514 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, q - 1,
515 2 * p - 1, 1.0);
516
517 for (i = 0; i < numPoints; ++i)
518 {
519 mode[i] *= one_m_z_pow[i] * one_p_z[i];
520 }
521 }
522 }
523
524 Blas::Dgemm('n', 'n', numPoints, numModes * (numModes + 1) / 2,
525 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
526 0.0, m_dbdata.data(), numPoints);
527 }
528 break;
529
530 case eModified_C:
531 {
532 // Note the following packing deviates from the
533 // definition in the Book by Karniadakis &
534 // Sherwin in two ways. 1) We put the vertex
535 // degrees of freedom at the lower index range to
536 // follow a more hierarchic structure. 2) We do
537 // not duplicate the singular vertex definition
538 // (or the duplicated face information in the book
539 // ) so that only a tetrahedral number
540 // (i.e. (modes)*(modes+1)*(modes+2)/6) of modes
541 // are required consistent with the orthogonal
542 // basis.
543
544 // In the current structure the r index runs
545 // fastest rollowed by q and than the p index so
546 // that the matrix has a more compact structure
547
548 // Note that eModified_C is a re-organisation/
549 // duplication of eModified_B so will get a
550 // temporary Modified_B expansion and copy the
551 // correct components.
552
553 // Generate Modified_B basis;
556 BasisSharedPtr ModB = BasisManager()[ModBKey];
557
558 Array<OneD, const NekDouble> ModB_data = ModB->GetBdata();
559
560 // Copy Modified_B basis into first
561 // (numModes*(numModes+1)/2)*numPoints entires of
562 // bdata. This fills in the complete (r,p) face.
563
564 // Set up \phi^c_{p,q,r} = \phi^b_{p+q,r}
565
566 size_t N;
567 size_t B_offset = 0;
568 size_t offset = 0;
569 for (p = 0; p < numModes; ++p)
570 {
571 N = numPoints * (numModes - p) * (numModes - p + 1) / 2;
572 Vmath::Vcopy(N, &ModB_data[0] + B_offset, 1,
573 &m_bdata[0] + offset, 1);
574 B_offset += numPoints * (numModes - p);
575 offset += N;
576 }
577
578 // set up derivative of basis.
579 Blas::Dgemm('n', 'n', numPoints,
580 numModes * (numModes + 1) * (numModes + 2) / 6,
581 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
582 0.0, m_dbdata.data(), numPoints);
583 }
584 break;
585
586 case eModifiedPyr_C:
587 {
588 // Note the following packing deviates from the
589 // definition in the Book by Karniadakis &
590 // Sherwin in two ways. 1) We put the vertex
591 // degrees of freedom at the lower index range to
592 // follow a more hierarchic structure. 2) We do
593 // not duplicate the singular vertex definition
594 // so that only a pyramidic number
595 // (i.e. (modes)*(modes+1)*(2*modes+1)/6) of modes
596 // are required consistent with the orthogonal
597 // basis.
598
599 // In the current structure the r index runs
600 // fastest rollowed by q and than the p index so
601 // that the matrix has a more compact structure
602
603 // Generate Modified_B basis;
606 BasisSharedPtr ModB = BasisManager()[ModBKey];
607
608 Array<OneD, const NekDouble> ModB_data = ModB->GetBdata();
609
610 // Copy Modified_B basis into first
611 // (numModes*(numModes+1)/2)*numPoints entires of
612 // bdata.
613
614 size_t N;
615 size_t B_offset = 0;
616 size_t offset = 0;
617
618 // Vertex 0,3,4, edges 3,4,7, face 4
619 N = numPoints * (numModes) * (numModes + 1) / 2;
620 Vmath::Vcopy(N, &ModB_data[0], 1, &m_bdata[0], 1);
621 offset += N;
622
623 B_offset += numPoints * (numModes);
624 // Vertex 1 edges 5
625 N = numPoints * (numModes - 1);
626 Vmath::Vcopy(N, &ModB_data[0] + B_offset, 1, &m_bdata[0] + offset,
627 1);
628 offset += N;
629
630 // Vertex 2 edges 1,6, face 2
631 N = numPoints * (numModes - 1) * (numModes) / 2;
632 Vmath::Vcopy(N, &ModB_data[0] + B_offset, 1, &m_bdata[0] + offset,
633 1);
634 offset += N;
635
636 B_offset += numPoints * (numModes - 1);
637
638 NekDouble *one_m_z_pow, *one_p_z;
639
640 mode = m_bdata.data() + offset;
641
642 for (p = 2; p < numModes; ++p)
643 {
644 // edges 0 2, faces 1 3
645 N = numPoints * (numModes - p);
646 Vmath::Vcopy(N, &ModB_data[0] + B_offset, 1, mode, 1);
647 mode += N;
648 Vmath::Vcopy(N, &ModB_data[0] + B_offset, 1, mode, 1);
649 mode += N;
650 B_offset += N;
651
652 one_p_z = m_bdata.data() + numPoints;
653
654 for (q = 2; q < numModes; ++q)
655 {
656 // face 0
657 for (i = 0; i < numPoints; ++i)
658 {
659 // [(1-z)/2]^{p+q-2} Note in book it
660 // seems to suggest p+q-1 but that
661 // does not seem to give complete
662 // polynomial space for pyramids
663 mode[i] = pow(m_bdata[i], p + q - 2);
664 }
665
666 one_m_z_pow = mode;
667 mode += numPoints;
668
669 // interior
670 for (size_t r = 1; r < numModes - std::max(p, q); ++r)
671 {
672 Polylib::jacobfd(numPoints, z.data(), mode, nullptr,
673 r - 1, 2 * p + 2 * q - 3, 1.0);
674
675 for (i = 0; i < numPoints; ++i)
676 {
677 mode[i] *= one_m_z_pow[i] * one_p_z[i];
678 }
679 mode += numPoints;
680 }
681 }
682 }
683
684 // set up derivative of basis.
685 Blas::Dgemm('n', 'n', numPoints,
686 numModes * (numModes + 1) * (2 * numModes + 1) / 6,
687 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
688 0.0, m_dbdata.data(), numPoints);
689 } // end scope
690 break;
691
692 case eGLL_Lagrange:
693 {
694 mode = m_bdata.data();
695 std::shared_ptr<Points<NekDouble>> m_points =
697 const Array<OneD, const NekDouble> &zp(m_points->GetZ());
698
699 for (p = 0; p < numModes; ++p, mode += numPoints)
700 {
701 for (q = 0; q < numPoints; ++q)
702 {
703 mode[q] =
704 Polylib::hglj(p, z[q], zp.data(), numModes, 0.0, 0.0);
705 }
706 }
707
708 // Define derivative basis
709 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
710 numPoints, m_bdata.data(), numPoints, 0.0,
711 m_dbdata.data(), numPoints);
712
713 } // end scope
714 break;
715
716 case eGauss_Lagrange:
717 {
718 mode = m_bdata.data();
719 std::shared_ptr<Points<NekDouble>> m_points =
721 const Array<OneD, const NekDouble> &zp(m_points->GetZ());
722
723 for (p = 0; p < numModes; ++p, mode += numPoints)
724 {
725 for (q = 0; q < numPoints; ++q)
726 {
727 mode[q] =
728 Polylib::hgj(p, z[q], zp.data(), numModes, 0.0, 0.0);
729 }
730 }
731
732 // Define derivative basis
733 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
734 numPoints, m_bdata.data(), numPoints, 0.0,
735 m_dbdata.data(), numPoints);
736
737 } // end scope
738 break;
739
740 case eFourier:
741 {
742 ASSERTL0(numModes % 2 == 0,
743 "Fourier modes should be a factor of 2");
744
745 for (i = 0; i < numPoints; ++i)
746 {
747 m_bdata[i] = 1.0;
748 m_bdata[numPoints + i] = 0.0;
749
750 m_dbdata[i] = m_dbdata[numPoints + i] = 0.0;
751 }
752
753 for (p = 1; p < numModes / 2; ++p)
754 {
755 for (i = 0; i < numPoints; ++i)
756 {
757 m_bdata[2 * p * numPoints + i] = cos(p * M_PI * (z[i] + 1));
758 m_bdata[(2 * p + 1) * numPoints + i] =
759 -sin(p * M_PI * (z[i] + 1));
760
761 m_dbdata[2 * p * numPoints + i] =
762 -p * M_PI * sin(p * M_PI * (z[i] + 1));
763 m_dbdata[(2 * p + 1) * numPoints + i] =
764 -p * M_PI * cos(p * M_PI * (z[i] + 1));
765 }
766 }
767 } // end scope
768 break;
769
770 // Fourier Single Mode (1st mode)
772 {
773 for (i = 0; i < numPoints; ++i)
774 {
775 m_bdata[i] = cos(M_PI * (z[i] + 1));
776 m_bdata[numPoints + i] = -sin(M_PI * (z[i] + 1));
777
778 m_dbdata[i] = -M_PI * sin(M_PI * (z[i] + 1));
779 m_dbdata[numPoints + i] = -M_PI * cos(M_PI * (z[i] + 1));
780 }
781
782 for (p = 1; p < numModes / 2; ++p)
783 {
784 for (i = 0; i < numPoints; ++i)
785 {
786 m_bdata[2 * p * numPoints + i] = 0.;
787 m_bdata[(2 * p + 1) * numPoints + i] = 0.;
788
789 m_dbdata[2 * p * numPoints + i] = 0.;
790 m_dbdata[(2 * p + 1) * numPoints + i] = 0.;
791 }
792 }
793 } // end scope
794 break;
795
796 // Fourier Real Half Mode
798 {
799 m_bdata[0] = cos(M_PI * z[0]);
800 m_dbdata[0] = -M_PI * sin(M_PI * z[0]);
801 } // end scope
802 break;
803
804 // Fourier Imaginary Half Mode
806 {
807 m_bdata[0] = -sin(M_PI * z[0]);
808 m_dbdata[0] = -M_PI * cos(M_PI * z[0]);
809 } // end scope
810 break;
811
812 case eChebyshev:
813 {
814 mode = m_bdata.data();
815
816 for (p = 0, scal = 1; p < numModes; ++p, mode += numPoints)
817 {
818 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, p, -0.5,
819 -0.5);
820
821 for (i = 0; i < numPoints; ++i)
822 {
823 mode[i] *= scal;
824 }
825
826 scal *= 4 * (p + 1) * (p + 1) / (2 * p + 2) / (2 * p + 1);
827 }
828
829 // Define derivative basis
830 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
831 numPoints, m_bdata.data(), numPoints, 0.0,
832 m_dbdata.data(), numPoints);
833 } // end scope
834 break;
835
836 case eMonomial:
837 {
838 mode = m_bdata.data();
839
840 for (size_t p = 0; p < numModes; ++p, mode += numPoints)
841 {
842 for (size_t i = 0; i < numPoints; ++i)
843 {
844 mode[i] = pow(z[i], p);
845 }
846 }
847
848 // Define derivative basis
849 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
850 numPoints, m_bdata.data(), numPoints, 0.0,
851 m_dbdata.data(), numPoints);
852 } // end scope
853 break;
854
855 default:
856 NEKERROR(ErrorUtil::efatal, "Basis Type not known or "
857 "not implemented at this time.");
858 }
859}
860
861/** \brief Determine if polynomial basis can be exactly integrated
862 * with itself
863 */
865{
866 bool returnval = false;
867
868 switch (GetPointsType())
869 {
874 case eGaussKronrodLegendre:
877 returnval = (GetNumPoints() >= GetNumModes());
878 break;
879
880 default:
881 break;
882 }
883
884 return returnval;
885}
886
887/** \brief Determine if basis has collocation property,
888 * i.e. GLL_Lagrange with Lobatto integration of appropriate order,
889 * Gauss_Lagrange with Gauss integration of appropriate order.
890 */
892{
893 return ((m_basistype == eGLL_Lagrange &&
895 GetNumModes() == GetNumPoints()) ||
898 GetNumModes() == GetNumPoints()));
899}
900
901// BasisKey compared to BasisKey
902bool operator==(const BasisKey &x, const BasisKey &y)
903{
904 return (x.GetPointsKey() == y.GetPointsKey() &&
905 x.m_basistype == y.m_basistype &&
906 x.GetNumModes() == y.GetNumModes());
907}
908
909// BasisKey* compared to BasisKey
910bool operator==(const BasisKey *x, const BasisKey &y)
911{
912 return (*x == y);
913}
914
915// \brief BasisKey compared to BasisKey*
916bool operator==(const BasisKey &x, const BasisKey *y)
917{
918 return (x == *y);
919}
920
921// \brief BasisKey compared to BasisKey
922bool operator!=(const BasisKey &x, const BasisKey &y)
923{
924 return (!(x == y));
925}
926
927// BasisKey* compared to BasisKey
928bool operator!=(const BasisKey *x, const BasisKey &y)
929{
930 return (!(*x == y));
931}
932
933// BasisKey compared to BasisKey*
934bool operator!=(const BasisKey &x, const BasisKey *y)
935{
936 return (!(x == *y));
937}
938
939} // namespace Nektar::LibUtilities
#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
int GetNumModes() const
Return order of basis from the basis specification.
Definition: Basis.h:207
Array< OneD, NekDouble > m_dbdata
Derivative Basis definition.
Definition: Basis.h:324
Basis()=delete
Private default constructor.
static std::shared_ptr< Basis > Create(const BasisKey &bkey)
Returns a new instance of a Basis with given BasisKey.
PointsSharedPtr m_points
Set of points.
Definition: Basis.h:322
int GetTotNumPoints() const
Return total number of points from the basis specification.
Definition: Basis.h:225
BasisKey m_basisKey
Basis specification.
Definition: Basis.h:321
BasisType GetBasisType() const
Return the type of expansion basis.
Definition: Basis.h:231
int GetNumPoints() const
Return the number of points from the basis specification.
Definition: Basis.h:219
std::shared_ptr< NekMatrix< NekDouble > > CalculateInterpMatrix(const BasisKey &tbasis0)
Calculate the interpolation Matrix for coefficient from one base (m_basisKey) to another (tbasis0)
NekManager< BasisKey, NekMatrix< NekDouble >, BasisKey::opLess > m_InterpManager
Definition: Basis.h:326
void GenBasis()
Generate appropriate basis and their derivatives.
Array< OneD, NekDouble > m_bdata
Basis definition.
Definition: Basis.h:323
static bool initBasisManager
Definition: Basis.h:329
Describes the specification for a Basis.
Definition: Basis.h:45
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:120
BasisType GetBasisType() const
Return type of expansion basis.
Definition: Basis.h:131
PointsKey GetPointsKey() const
Return distribution of points.
Definition: Basis.h:137
size_t m_nummodes
Expansion order.
Definition: Basis.h:185
int GetTotNumPoints() const
Definition: Basis.h:125
int GetTotNumModes() const
Definition: Basis.h:80
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:74
PointsType GetPointsType() const
Return type of quadrature.
Definition: Basis.h:143
BasisType m_basistype
Expansion type.
Definition: Basis.h:186
bool Collocation() const
Determine if basis has collocation properties.
bool ExactIprodInt() const
Determine if basis has exact integration for inner product.
bool RegisterGlobalCreator(const CreateFuncType &createFunc)
Register the Global Create Function. The return value is just to facilitate calling statically.
Definition: NekManager.hpp:177
Defines a specification for a set of points.
Definition: Points.h:50
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
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
BasisManagerT & BasisManager(void)
bool operator==(const BasisKey &x, const BasisKey &y)
std::shared_ptr< Basis > BasisSharedPtr
const char *const BasisTypeMap[]
Definition: Foundations.hpp:44
PointsManagerT & PointsManager(void)
bool operator<(const BasisKey &lhs, const BasisKey &rhs)
bool operator>(const BasisKey &lhs, const BasisKey &rhs)
bool operator!=(const BasisKey &x, const BasisKey &y)
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:47
@ eGaussRadauKronrodMLegendre
1D Gauss-Radau-Kronrod-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:67
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eGaussLobattoKronrodLegendre
1D Lobatto Kronrod quadrature points
Definition: PointsType.h:72
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:46
@ eGaussRadauPLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition: PointsType.h:49
@ 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
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:50
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
Definition: BasisType.h:64
@ eChebyshev
Chebyshev Polynomials.
Definition: BasisType.h:61
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eModifiedPyr_C
Principle Modified Functions.
Definition: BasisType.h:53
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:68
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:66
@ eOrthoPyr_C
Principle Orthogonal Functions .
Definition: BasisType.h:51
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::vector< double > w(NPUPPER)
std::vector< double > z(NPUPPER)
std::vector< double > q(NPUPPER *NPUPPER)
double NekDouble
double hgj(const int i, const double z, const double *zgj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Jacobi points zgj at the ar...
Definition: Polylib.cpp:914
double hglj(const int i, const double z, const double *zglj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Lobatto-Jacobi points zgrj ...
Definition: Polylib.cpp:1025
void jacobfd(const int np, const double *z, double *poly_in, double *polyd, const int n, const double alpha, const double beta)
Routine to calculate Jacobi polynomials, , and their first derivative, .
Definition: Polylib.cpp:1248
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294
bool operator()(const BasisKey &lhs, const BasisKey &rhs) const