HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
Loading...
Searching...
No Matches
basis.hpp
Go to the documentation of this file.
1// Core data structures and methods required to implement polytopal methods in 3D
2//
3// Provides:
4// - Full and partial polynomial spaces on the element, faces, and edges
5// - Generic routines to create quadrature nodes over cells and faces of the mesh
6//
7// Authors: Daniele Di Pietro (daniele.di-pietro@umontpellier.fr)
8// Jerome Droniou (jerome.droniou@monash.edu)
9
10/*
11 *
12 * This library was developed around HHO methods, although some parts of it have a more
13 * general purpose. If you use this code or part of it in a scientific publication,
14 * please mention the following book as a reference for the underlying principles
15 * of HHO schemes:
16 *
17 * The Hybrid High-Order Method for Polytopal Meshes: Design, Analysis, and Applications.
18 * D. A. Di Pietro and J. Droniou. Modeling, Simulation and Applications, vol. 19.
19 * Springer International Publishing, 2020, xxxi + 525p. doi: 10.1007/978-3-030-37203-3.
20 * url: https://hal.archives-ouvertes.fr/hal-02151813.
21 *
22 */
23
24#ifndef BASIS_HPP
25#define BASIS_HPP
26
27#include <boost/multi_array.hpp>
28
29#include <mesh.hpp>
30#include <iostream>
31
33
34#include <quadraturerule.hpp>
35
36namespace HArDCore3D
37{
38
50 constexpr int dimspace = 3;
51 typedef Eigen::Matrix3d MatrixRd;
52 static const std::vector<VectorRd> basisRd = { VectorRd(1., 0., 0.), VectorRd(0., 1., 0.), VectorRd(0., 0., 1.) };
53
54 template <typename T>
55 using FType = std::function<T(const VectorRd&)>;
56
57 template <typename T>
58 using CellFType = std::function<T(const VectorRd&, const Cell *)>;
59
60 template <typename T>
61 using BasisQuad = boost::multi_array<T, 2>;
62
64 {
65 Scalar = 0,
66 Vector = 1,
67 Matrix = 2
68 };
69
70 //-----------------------------------------------------------------------------
71 // Powers of monomial basis
72 //-----------------------------------------------------------------------------
73
75 template<typename GeometricSupport>
77 {
78 // Only specializations are relevant
79 };
80
81 template<>
82 struct MonomialPowers<Cell>
83 {
84 // Powers of homogeneous polynomials of degree L
85 static std::vector<VectorZd> homogeneous(const size_t L)
86 {
87 std::vector<VectorZd> powers;
89 for (size_t i = 0; i <= L; i++)
90 {
91 for (size_t j = 0; i + j <= L; j++)
92 {
93 powers.push_back(VectorZd(i, j, L - i - j));
94 } // for j
95 } // for i
96 return powers;
97 }
98
99 // Powers for degrees from 0 to degree
100 static std::vector<VectorZd> complete(const size_t degree)
101 {
102 std::vector<VectorZd> powers;
103 powers.reserve( PolynomialSpaceDimension<Cell>::Poly(degree) );
104 for (size_t l = 0; l <= degree; l++)
105 {
106 std::vector<VectorZd> pow_hom = homogeneous(l);
107 for (VectorZd p : pow_hom){
108 powers.push_back(p);
109 }
110 } // for l
111 return powers;
112 }
113 };
114
115 template<>
116 struct MonomialPowers<Face>
117 {
118 // Powers of homogeneous polynomials of degree L
119 static std::vector<Eigen::Vector2i> homogeneous(const size_t L)
120 {
121 std::vector<Eigen::Vector2i> powers;
123 for (size_t i = 0; i <= L; i++)
124 {
125 powers.push_back(Eigen::Vector2i(i, L - i));
126 } // for i
127 return powers;
128 }
129
130 // Powers for degrees from 0 to degree
131 static std::vector<Eigen::Vector2i> complete(size_t degree)
132 {
133 std::vector<Eigen::Vector2i> powers;
134 powers.reserve( PolynomialSpaceDimension<Face>::Poly(degree) );
135 for (size_t l = 0; l <= degree; l++)
136 {
137 std::vector<Eigen::Vector2i> pow_hom = homogeneous(l);
138 for (Eigen::Vector2i p : pow_hom){
139 powers.push_back(p);
140 }
141 } // for l
142
143 return powers;
144 }
145 };
146
147 //----------------------------------------------------------------------
148 //----------------------------------------------------------------------
149 // SCALAR MONOMIAL BASIS
150 //----------------------------------------------------------------------
151 //----------------------------------------------------------------------
152
155 {
156 public:
157 typedef double FunctionValue;
158 typedef VectorRd GradientValue;
159 typedef VectorRd CurlValue;
160 typedef double DivergenceValue;
162
163 typedef Cell GeometricSupport;
164
165 constexpr static const TensorRankE tensorRank = Scalar;
166 constexpr static const bool hasAncestor = false;
167 static const bool hasFunction = true;
168 static const bool hasGradient = true;
169 static const bool hasCurl = false;
170 static const bool hasDivergence = false;
171 static const bool hasHessian = true;
172 static const bool hasCurlCurl = false;
173
176 const Cell &T,
177 size_t degree
178 );
179
181 inline size_t dimension() const
182 {
183 return (m_degree >= 0 ? (m_degree * (m_degree + 1) * (2 * m_degree + 1) + 9 * m_degree * (m_degree + 1) + 12 * (m_degree + 1)) / 12 : 0);
184 }
185
187 FunctionValue function(size_t i, const VectorRd &x) const;
188
190 GradientValue gradient(size_t i, const VectorRd &x) const;
191
193 HessianValue hessian(size_t i, const VectorRd &x) const;
194
196 inline size_t max_degree() const
197 {
198 return m_degree;
199 }
200
202 inline VectorZd powers(size_t i) const
203 {
204 return m_powers[i];
205 }
206
207 private:
209 inline VectorRd _coordinate_transform(const VectorRd &x) const
210 {
211 return (x - m_xT) / m_hT;
212 }
213
214 Cell m_T;
215 size_t m_degree;
216 VectorRd m_xT;
217 double m_hT;
218 std::vector<VectorZd> m_powers;
219 };
220
221
222 //------------------------------------------------------------------------------
223
226 {
227 public:
228 typedef double FunctionValue;
229 typedef VectorRd GradientValue;
230 typedef VectorRd CurlValue;
231 typedef double DivergenceValue;
233
234 typedef Face GeometricSupport;
235
236 typedef Eigen::Matrix<double, 2, dimspace> JacobianType;
237
238 constexpr static const TensorRankE tensorRank = Scalar;
239 constexpr static const bool hasAncestor = false;
240 static const bool hasFunction = true;
241 static const bool hasGradient = true;
242 static const bool hasCurl = true;
243 static const bool hasDivergence = false;
244 static const bool hasHessian = true;
245 static const bool hasCurlCurl = false;
246
249 const Face &F,
250 size_t degree
251 );
252
254 inline size_t dimension() const
255 {
256 return (m_degree + 1) * (m_degree + 2) / 2;
257 }
258
260 FunctionValue function(size_t i, const VectorRd &x) const;
261
263 GradientValue gradient(size_t i, const VectorRd &x) const;
264
266 CurlValue curl(size_t i, const VectorRd &x) const;
267
269 HessianValue hessian(size_t i, const VectorRd &x) const;
270
272 inline size_t max_degree() const
273 {
274 return m_degree;
275 }
276
278 inline const VectorRd &normal() const
279 {
280 return m_nF;
281 }
282
284 inline const JacobianType &jacobian() const
285 {
286 return m_jacobian;
287 }
288
290 inline const JacobianType coordinates_system() const
291 {
292 return m_jacobian * m_hF;
293 }
294
296 inline Eigen::Vector2i powers(size_t i) const
297 {
298 return m_powers[i];
299 }
300
301 private:
303 inline Eigen::Vector2d _coordinate_transform(const VectorRd &x) const
304 {
305 return m_jacobian * (x - m_xF);
306 }
307
308 size_t m_degree;
309 VectorRd m_xF;
310 double m_hF;
311 VectorRd m_nF;
312 JacobianType m_jacobian;
313 std::vector<Eigen::Vector2i> m_powers;
314 };
315
316 //------------------------------------------------------------------------------
317
320 {
321 public:
322 typedef double FunctionValue;
323 typedef VectorRd GradientValue;
324 typedef VectorRd CurlValue;
325 typedef double DivergenceValue;
326 typedef double HessianValue;
327
328 typedef Edge GeometricSupport;
329
330 constexpr static const TensorRankE tensorRank = Scalar;
331 constexpr static const bool hasAncestor = false;
332 static const bool hasFunction = true;
333 static const bool hasGradient = true;
334 static const bool hasCurl = false;
335 static const bool hasDivergence = false;
336 static const bool hasHessian = false;
337 static const bool hasCurlCurl = false;
338
341 const Edge &E,
342 size_t degree
343 );
344
346 inline size_t dimension() const
347 {
348 return m_degree + 1;
349 }
350
352 FunctionValue function(size_t i, const VectorRd &x) const;
353
355 GradientValue gradient(size_t i, const VectorRd &x) const;
356
358 inline size_t max_degree() const
359 {
360 return m_degree;
361 }
362
363 private:
364 inline double _coordinate_transform(const VectorRd &x) const
365 {
366 return (x - m_xE).dot(m_tE) / m_hE;
367 }
368
369 size_t m_degree;
370 VectorRd m_xE;
371 double m_hE;
372 VectorRd m_tE;
373 };
374
375 //------------------------------------------------------------------------------
376 //------------------------------------------------------------------------------
377 // DERIVED BASIS
378 //------------------------------------------------------------------------------
379 //------------------------------------------------------------------------------
380
381 //----------------------------FAMILY------------------------------------------
382 //
384
387 template <typename BasisType>
388 class Family
389 {
390 public:
391 typedef typename BasisType::FunctionValue FunctionValue;
392 typedef typename BasisType::GradientValue GradientValue;
393 typedef typename BasisType::CurlValue CurlValue;
394 typedef typename BasisType::DivergenceValue DivergenceValue;
395 typedef typename BasisType::HessianValue HessianValue;
396
397 typedef typename BasisType::GeometricSupport GeometricSupport;
398
399 constexpr static const TensorRankE tensorRank = BasisType::tensorRank;
400 constexpr static const bool hasAncestor = true;
401 static const bool hasFunction = BasisType::hasFunction;
402 static const bool hasGradient = BasisType::hasGradient;
403 static const bool hasCurl = BasisType::hasCurl;
404 static const bool hasDivergence = BasisType::hasDivergence;
405 static const bool hasHessian = BasisType::hasHessian;
406 static const bool hasCurlCurl = BasisType::hasCurlCurl;
407
409
412 const BasisType &basis,
413 const Eigen::MatrixXd &matrix
414 )
415 : m_basis(basis),
416 m_matrix(matrix)
417 {
418 assert((size_t)matrix.cols() == basis.dimension() || "Inconsistent family initialization");
419 }
420
422 inline size_t dimension() const
423 {
424 return m_matrix.rows();
425 }
426
428 FunctionValue function(size_t i, const VectorRd &x) const
429 {
430 static_assert(hasFunction, "Call to function() not available");
431
432 FunctionValue f = m_matrix(i, 0) * m_basis.function(0, x);
433 for (auto j = 1; j < m_matrix.cols(); j++)
434 {
435 f += m_matrix(i, j) * m_basis.function(j, x);
436 } // for j
437 return f;
438 }
439
441 FunctionValue function(size_t i, size_t iqn, const boost::multi_array<FunctionValue, 2> &ancestor_value_quad) const
442 {
443 static_assert(hasFunction, "Call to function() not available");
444
445 FunctionValue f = m_matrix(i, 0) * ancestor_value_quad[0][iqn];
446 for (auto j = 1; j < m_matrix.cols(); j++)
447 {
448 f += m_matrix(i, j) * ancestor_value_quad[j][iqn];
449 } // for j
450 return f;
451 }
452
453
455 GradientValue gradient(size_t i, const VectorRd &x) const
456 {
457 static_assert(hasGradient, "Call to gradient() not available");
458
459 GradientValue G = m_matrix(i, 0) * m_basis.gradient(0, x);
460 for (auto j = 1; j < m_matrix.cols(); j++)
461 {
462 G += m_matrix(i, j) * m_basis.gradient(j, x);
463 } // for j
464 return G;
465 }
466
468 GradientValue gradient(size_t i, size_t iqn, const boost::multi_array<GradientValue, 2> &ancestor_gradient_quad) const
469 {
470 static_assert(hasGradient, "Call to gradient() not available");
471
472 GradientValue G = m_matrix(i, 0) * ancestor_gradient_quad[0][iqn];
473 for (auto j = 1; j < m_matrix.cols(); j++)
474 {
475 G += m_matrix(i, j) * ancestor_gradient_quad[j][iqn];
476 } // for j
477 return G;
478 }
479
481 CurlValue curl(size_t i, const VectorRd &x) const
482 {
483 static_assert(hasCurl, "Call to curl() not available");
484
485 CurlValue C = m_matrix(i, 0) * m_basis.curl(0, x);
486 for (auto j = 1; j < m_matrix.cols(); j++)
487 {
488 C += m_matrix(i, j) * m_basis.curl(j, x);
489 } // for j
490 return C;
491 }
492
494 CurlValue curl(size_t i, size_t iqn, const boost::multi_array<CurlValue, 2> &ancestor_curl_quad) const
495 {
496 static_assert(hasCurl, "Call to curl() not available");
497
498 CurlValue C = m_matrix(i, 0) * ancestor_curl_quad[0][iqn];
499 for (auto j = 1; j < m_matrix.cols(); j++)
500 {
501 C += m_matrix(i, j) * ancestor_curl_quad[j][iqn];
502 } // for j
503 return C;
504 }
505
507 DivergenceValue divergence(size_t i, const VectorRd &x) const
508 {
509 static_assert(hasDivergence, "Call to divergence() not available");
510
511 DivergenceValue D = m_matrix(i, 0) * m_basis.divergence(0, x);
512 for (auto j = 1; j < m_matrix.cols(); j++)
513 {
514 D += m_matrix(i, j) * m_basis.divergence(j, x);
515 } // for j
516 return D;
517 }
518
520 DivergenceValue divergence(size_t i, size_t iqn, const boost::multi_array<DivergenceValue, 2> &ancestor_divergence_quad) const
521 {
522 static_assert(hasDivergence, "Call to divergence() not available");
523
524 DivergenceValue D = m_matrix(i, 0) * ancestor_divergence_quad[0][iqn];
525 for (auto j = 1; j < m_matrix.cols(); j++)
526 {
527 D += m_matrix(i, j) * ancestor_divergence_quad[j][iqn];
528 } // for j
529 return D;
530 }
531
533 HessianValue hessian(size_t i, const VectorRd &x) const
534 {
535 static_assert(hasHessian, "Call to hessian() not available");
536
537 HessianValue H = m_matrix(i, 0) * m_basis.hessian(0, x);
538 for (auto j = 1; j < m_matrix.cols(); j++)
539 {
540 H += m_matrix(i, j) * m_basis.hessian(j, x);
541 } // for j
542 return H;
543 }
544
546 HessianValue hessian(size_t i, size_t iqn, const boost::multi_array<HessianValue, 2> &ancestor_hessian_quad) const
547 {
548 static_assert(hasHessian, "Call to hessian() not available");
549
550 HessianValue H = m_matrix(i, 0) * ancestor_hessian_quad[0][iqn];
551 for (auto j = 1; j < m_matrix.cols(); j++)
552 {
553 H += m_matrix(i, j) * ancestor_hessian_quad[j][iqn];
554 } // for j
555 return H;
556 }
557
559 inline const Eigen::MatrixXd & matrix() const
560 {
561 return m_matrix;
562 }
563
565 constexpr inline const BasisType &ancestor() const
566 {
567 return m_basis;
568 }
569
571 inline size_t max_degree() const
572 {
573 return m_basis.max_degree();
574 }
575
576 private:
577 BasisType m_basis;
578 Eigen::MatrixXd m_matrix;
579 };
580
581 //----------------------TENSORIZED--------------------------------------------------------
582
584
608 template <typename ScalarFamilyType, size_t N>
610 {
611 public:
612 typedef typename Eigen::Matrix<double, N, 1> FunctionValue;
613 typedef typename Eigen::Matrix<double, N, dimspace> GradientValue;
614 typedef VectorRd CurlValue;
615 typedef double DivergenceValue;
616 typedef void HessianValue;
617
618 typedef typename ScalarFamilyType::GeometricSupport GeometricSupport;
619
620 constexpr static const TensorRankE tensorRank = Vector;
621 constexpr static const bool hasAncestor = true;
622 static const bool hasFunction = ScalarFamilyType::hasFunction;
623 static const bool hasGradient = ScalarFamilyType::hasGradient;
624 // We know how to compute the curl and divergence if gradient is available, and
625 // if we tensorize at the dimension of the space
626 static const bool hasDivergence = (ScalarFamilyType::hasGradient && N == dimspace);
627 static const bool hasCurl = (ScalarFamilyType::hasGradient && N == dimspace);
628 static const bool hasHessian = false;
629 static const bool hasCurlCurl = (ScalarFamilyType::hasHessian && N == dimspace);
630
632
634 : m_scalar_family(scalar_family)
635 {
636 static_assert(ScalarFamilyType::tensorRank == Scalar,
637 "Vector family can only be constructed from scalar families");
638 }
639
641 inline size_t dimension() const
642 {
643 return m_scalar_family.dimension() * N;
644 }
645
647 FunctionValue function(size_t i, const VectorRd &x) const
648 {
649 static_assert(hasFunction, "Call to function() not available");
650
651 FunctionValue ek = Eigen::Matrix<double, N, 1>::Zero();
652 ek(i / m_scalar_family.dimension()) = 1.;
653 return ek * m_scalar_family.function(i % m_scalar_family.dimension(), x);
654 }
655
657 FunctionValue function(size_t i, size_t iqn, const boost::multi_array<double, 2> &ancestor_value_quad) const
658 {
659 static_assert(hasFunction, "Call to function() not available");
660
661 FunctionValue ek = Eigen::Matrix<double, N, 1>::Zero();
662 ek(i / m_scalar_family.dimension()) = 1.;
663 return ek * ancestor_value_quad[i % m_scalar_family.dimension()][iqn];
664 }
665
667 GradientValue gradient(size_t i, const VectorRd &x) const
668 {
669 static_assert(hasGradient, "Call to gradient() not available");
670
671 GradientValue G = Eigen::Matrix<double, N, dimspace>::Zero();
672 G.row(i / m_scalar_family.dimension()) = m_scalar_family.gradient(i % m_scalar_family.dimension(), x);
673 return G;
674 }
675
677 GradientValue gradient(size_t i, size_t iqn, const boost::multi_array<VectorRd, 2> &ancestor_gradient_quad) const
678 {
679 static_assert(hasGradient, "Call to gradient() not available");
680
681 GradientValue G = Eigen::Matrix<double, N, dimspace>::Zero();
682 G.row(i / m_scalar_family.dimension()) = ancestor_gradient_quad[i % m_scalar_family.dimension()][iqn];
683 return G;
684 }
685
687 CurlValue curl(size_t i, const VectorRd &x) const
688 {
689 static_assert(hasCurl, "Call to curl() not available");
690
691 VectorRd ek = VectorRd::Zero();
692 ek(i / m_scalar_family.dimension()) = 1.;
693 return m_scalar_family.gradient(i % m_scalar_family.dimension(), x).cross(ek);
694 }
695
697 CurlValue curl(size_t i, size_t iqn, const boost::multi_array<VectorRd, 2> &ancestor_gradient_quad) const
698 {
699 static_assert(hasCurl, "Call to curl() not available");
700
701 VectorRd ek = VectorRd::Zero();
702 ek(i / m_scalar_family.dimension()) = 1.;
703 return ancestor_gradient_quad[i % m_scalar_family.dimension()][iqn].cross(ek);
704 }
705
707 DivergenceValue divergence(size_t i, const VectorRd &x) const
708 {
709 static_assert(hasDivergence, "Call to divergence() not available");
710
711 return m_scalar_family.gradient(i % m_scalar_family.dimension(), x)(i / m_scalar_family.dimension());
712 }
713
715 DivergenceValue divergence(size_t i, size_t iqn, const boost::multi_array<VectorRd, 2> &ancestor_gradient_quad) const
716 {
717 static_assert(hasDivergence, "Call to divergence() not available");
718
719 return ancestor_gradient_quad[i % m_scalar_family.dimension()][iqn](i / m_scalar_family.dimension());
720 }
721
723 CurlValue curlcurl(size_t i, const VectorRd &x) const
724 {
725 static_assert(hasCurlCurl, "Call to curlcurl() not available");
726
727 FunctionValue ek = Eigen::Matrix<double, N, 1>::Zero();
728 int p = i / m_scalar_family.dimension();
729 ek(p) = 1.;
730 typename ScalarFamilyType::HessianValue hess = m_scalar_family.hessian(i % m_scalar_family.dimension(), x);
731 return - hess.trace() * ek + hess.col(p);
732 }
733
735 CurlValue curlcurl(size_t i, size_t iqn, const boost::multi_array<typename ScalarFamilyType::HessianValue, 2> &ancestor_hessian_quad) const
736 {
737 static_assert(hasCurlCurl, "Call to curlcurl() not available");
738
739 FunctionValue ek = Eigen::Matrix<double, N, 1>::Zero();
740 int p = i / m_scalar_family.dimension();
741 ek(p) = 1.;
742 typename ScalarFamilyType::HessianValue hess = ancestor_hessian_quad[i % m_scalar_family.dimension()][iqn];
743 return - hess.trace() * ek + hess.col(p);
744 }
745
746
748 constexpr inline const ScalarFamilyType &ancestor() const
749 {
750 return m_scalar_family;
751 }
752
754 inline size_t max_degree() const
755 {
756 return m_scalar_family.max_degree();
757 }
758
759 private:
760 ScalarFamilyType m_scalar_family;
761 };
762
763
764 //----------------------MATRIX FAMILY--------------------------------------------------------
765
767
772 // Useful formulas to manipulate this basis (all indices starting from 0):
773 // the i-th element in the basis is f_{i%r} E_{i/r}
774 // the matrix E_m has 1 at position (m%N, m/N)
775 // the element f_aE_b is the i-th in the family with i = b*r + a
776 template<typename ScalarFamilyType, size_t N>
778 {
779 public:
780 typedef typename Eigen::Matrix<double, N, N> FunctionValue;
781 typedef void GradientValue;
782 typedef void CurlValue;
783 typedef typename Eigen::Matrix<double, N, 1> DivergenceValue;
784 typedef void HessianValue;
785
786 typedef typename ScalarFamilyType::GeometricSupport GeometricSupport;
787
788 constexpr static const TensorRankE tensorRank = Matrix;
789 constexpr static const bool hasAncestor = true;
790 static const bool hasFunction = ScalarFamilyType::hasFunction;
791 static const bool hasGradient = false;
792 static const bool hasDivergence = ( ScalarFamilyType::hasGradient && N==dimspace );
793 static const bool hasCurl = false;
794 static const bool hasHessian = false;
795 static const bool hasCurlCurl = false;
796
798
800 : m_scalar_family(scalar_family),
801 m_E(N*N),
802 m_transposeOperator(Eigen::MatrixXd::Zero(scalar_family.dimension()*N*N, scalar_family.dimension()*N*N))
803 {
804 static_assert(ScalarFamilyType::tensorRank == Scalar,
805 "Vector family can only be constructed from scalar families");
806
807 // Construct the basis for NxN matrices
808 for (size_t j = 0; j < N; j++){
809 for (size_t i = 0; i < N; i++){
810 m_E[j*N + i] = Eigen::Matrix<double, N, N>::Zero();
811 m_E[j*N + i](i,j) = 1.;
812 }
813 }
814
815 // Construct the transpose operator
816 for (size_t i = 0; i < dimension(); i++){
817 // The i-th basis element is sent to the j-th basis element s.t. j%r=i%r, (j/r)%N=(i/r)/N, (j/r)/N=(i/r)%N
818 size_t r = m_scalar_family.dimension();
819 size_t j = ( ( (i/r)%N )*N + (i/r)/N ) * r + (i%r);
820 m_transposeOperator(i,j) = 1.;
821 }
822 }
823
825 inline size_t dimension() const
826 {
827 return m_scalar_family.dimension() * N * N;
828 }
829
831 FunctionValue function(size_t i, const VectorRd & x) const
832 {
833 static_assert(hasFunction, "Call to function() not available");
834
835 return m_scalar_family.function(i % m_scalar_family.dimension(), x) * m_E[i / m_scalar_family.dimension()];
836 }
837
839 FunctionValue function(size_t i, size_t iqn, const boost::multi_array<double, 2> &ancestor_value_quad) const
840 {
841 static_assert(hasFunction, "Call to function() not available");
842
843 return ancestor_value_quad[i % m_scalar_family.dimension()][iqn] * m_E[i / m_scalar_family.dimension()];
844 }
845
847 DivergenceValue divergence(size_t i, const VectorRd & x) const
848 {
849 static_assert(hasDivergence, "Call to divergence() not available");
850
851 Eigen::Matrix<double, N, 1> V = m_E[i/m_scalar_family.dimension()].col( (i/m_scalar_family.dimension()) /N);
852 return m_scalar_family.gradient(i % m_scalar_family.dimension(), x)( (i / m_scalar_family.dimension()) / N) * V;
853 }
854
856 DivergenceValue divergence(size_t i, size_t iqn, const boost::multi_array<VectorRd, 2> &ancestor_gradient_quad) const
857 {
858 static_assert(hasDivergence, "Call to divergence() not available");
859
860 Eigen::Matrix<double, N, 1> V = m_E[i/m_scalar_family.dimension()].col( (i/m_scalar_family.dimension()) /N);
861 return ancestor_gradient_quad[i % m_scalar_family.dimension()][iqn]( (i / m_scalar_family.dimension()) / N) * V;
862 }
863
865 inline const ScalarFamilyType &ancestor() const
866 {
867 return m_scalar_family;
868 }
869
871 inline const size_t matrixSize() const
872 {
873 return N;
874 }
875
877 inline const Eigen::MatrixXd transposeOperator() const
878 {
879 return m_transposeOperator;
880 }
881
883 inline const Eigen::MatrixXd symmetriseOperator() const
884 {
885 return (Eigen::MatrixXd::Identity(dimension(), dimension()) + m_transposeOperator)/2;
886 }
887
889 inline const Eigen::MatrixXd symmetricBasis() const
890 {
891 size_t dim_scalar = m_scalar_family.dimension();
892 Eigen::MatrixXd symOpe = symmetriseOperator();
893
894 Eigen::MatrixXd SB = Eigen::MatrixXd::Zero(dim_scalar * N*(N+1)/2, dimension());
895
896 // The matrix consists in selecting certain rows of symmetriseOperator, corresponding to the lower trianglular part
897 // To select these rows, we loop over the columns of the underlying matrices, and get the rows in symmetriseOperator()
898 // corresponding to the lower triangular part of each column
899 size_t position = 0;
900 for (size_t i = 0; i<N; i++){
901 // Skip the first i*dim_scalar elements in the column (upper triangle) in symOpe, take the remaining (N-i)*dim_scalar
902 SB.middleRows(position, (N-i)*dim_scalar) = symOpe.middleRows(i*(N+1)*dim_scalar, (N-i)*dim_scalar);
903 // update current position
904 position += (N-i)*dim_scalar;
905 }
906
907 return SB;
908 }
909
911
914 inline const Eigen::MatrixXd traceOperator() const
915 {
916 size_t r = m_scalar_family.dimension();
917
918 Eigen::MatrixXd Tr = Eigen::MatrixXd::Zero(r, dimension());
919
920 /* According to the formulas described at the start of this class, the l-th element of the
921 matrix family is f_{l%r} E_{l/r}.
922 Given the construction of E_m=E_{l/r}, the trace is either 0, or 1 if m%N=m/N, that is,
923 (l/r)%N = (l/r)/N.
924 So Tr_{kl}=1 only if k=l%r and (l/r)%N = (l/r)/N
925 */
926 for (size_t l=0; l < dimension(); l++){
927 if ( size_t(l/r)%N == size_t(size_t(l/r)/N) ){
928 size_t k = l%r;
929 Tr(k,l) = 1.;
930 }
931 }
932
933 return Tr;
934 }
935
936 private:
937 ScalarFamilyType m_scalar_family;
938 std::vector<Eigen::Matrix<double, N, N>> m_E;
939 Eigen::MatrixXd m_transposeOperator;
940 };
941
942
943 //----------------------TANGENT FAMILY--------------------------------------------------------
944
946 template <typename ScalarFamilyType>
948 {
949 public:
950 typedef VectorRd FunctionValue;
951 typedef VectorRd GradientValue;
952 typedef VectorRd CurlValue;
953 typedef double DivergenceValue;
954 typedef Eigen::Matrix2d HessianValue;
955
956 typedef Face GeometricSupport;
957
958 constexpr static const TensorRankE tensorRank = Vector;
959 constexpr static const bool hasAncestor = true;
960 static const bool hasFunction = true;
961 static const bool hasGradient = false;
962 static const bool hasCurl = false;
963 static const bool hasDivergence = ScalarFamilyType::hasGradient;
964 static const bool hasHessian = false;
965 static const bool hasCurlCurl = ScalarFamilyType::hasHessian;
966
968
972 const Eigen::Matrix<double, 2, dimspace> &generators
973 )
974 : m_scalar_family(scalar_family),
975 m_generators(generators),
976 m_normal( ( generators.row(0).cross(generators.row(1)) ).normalized() )
977 {
978 static_assert(ScalarFamilyType::hasFunction, "Call to function() not available");
979 static_assert(std::is_same<typename ScalarFamilyType::GeometricSupport, Face>::value,
980 "Tangent families can only be defined on faces");
981 }
982
984 inline size_t dimension() const
985 {
986 return m_scalar_family.dimension() * 2;
987 }
988
990 inline FunctionValue function(size_t i, const VectorRd &x) const
991 {
992 return m_generators.row(i / m_scalar_family.dimension()) * m_scalar_family.function(i % m_scalar_family.dimension(), x);
993 }
994
996 inline DivergenceValue divergence(size_t i, const VectorRd &x) const
997 {
998 return m_generators.row(i / m_scalar_family.dimension()).dot(m_scalar_family.gradient(i % m_scalar_family.dimension(), x));
999 }
1000
1002 inline CurlValue curlcurl(size_t i, const VectorRd &x) const
1003 {
1004 return m_normal.cross(
1005 m_scalar_family.hessian(i % m_scalar_family.dimension(), x) *
1006 (m_normal.cross(m_generators.row(i / m_scalar_family.dimension())) ) );
1007 }
1008
1010 constexpr inline const ScalarFamilyType &ancestor() const
1011 {
1012 return m_scalar_family;
1013 }
1014
1016 inline size_t max_degree() const
1017 {
1018 return m_scalar_family.max_degree();
1019 }
1020
1022 inline Eigen::Matrix<double, 2, dimspace> generators() const
1023 {
1024 return m_generators;
1025 }
1026
1027
1028 private:
1029 ScalarFamilyType m_scalar_family;
1030 Eigen::Matrix<double, 2, dimspace> m_generators;
1031 VectorRd m_normal; // unit normal
1032 };
1033
1034 //--------------------SHIFTED BASIS----------------------------------------------------------
1035
1037
1038 template <typename BasisType>
1040 {
1041 public:
1042 typedef typename BasisType::FunctionValue FunctionValue;
1043 typedef typename BasisType::GradientValue GradientValue;
1044 typedef typename BasisType::CurlValue CurlValue;
1045 typedef typename BasisType::DivergenceValue DivergenceValue;
1046 typedef typename BasisType::HessianValue HessianValue;
1047
1048 typedef typename BasisType::GeometricSupport GeometricSupport;
1049
1050 constexpr static const TensorRankE tensorRank = BasisType::tensorRank;
1051 constexpr static const bool hasAncestor = true;
1052 static const bool hasFunction = BasisType::hasFunction;
1053 static const bool hasGradient = BasisType::hasGradient;
1054 static const bool hasCurl = BasisType::hasCurl;
1055 static const bool hasDivergence = BasisType::hasDivergence;
1056 static const bool hasHessian = BasisType::hasHessian;
1057 static const bool hasCurlCurl = BasisType::hasCurlCurl;
1058
1060
1063 const BasisType &basis,
1064 const int shift
1065 )
1066 : m_basis(basis),
1067 m_shift(shift)
1068 {
1069 // Do nothing
1070 }
1071
1073 inline size_t dimension() const
1074 {
1075 return m_basis.dimension() - m_shift;
1076 }
1077
1079 inline size_t shift() const
1080 {
1081 return m_shift;
1082 }
1083
1085 constexpr inline const BasisType &ancestor() const
1086 {
1087 return m_basis;
1088 }
1089
1091 inline size_t max_degree() const
1092 {
1093 return m_basis.max_degree();
1094 }
1095
1097 inline FunctionValue function(size_t i, const VectorRd &x) const
1098 {
1099 static_assert(hasFunction, "Call to function() not available");
1100
1101 return m_basis.function(i + m_shift, x);
1102 }
1103
1105 inline GradientValue gradient(size_t i, const VectorRd &x) const
1106 {
1107 static_assert(hasGradient, "Call to gradient() not available");
1108
1109 return m_basis.gradient(i + m_shift, x);
1110 }
1111
1113 inline CurlValue curl(size_t i, const VectorRd &x) const
1114 {
1115 static_assert(hasCurl, "Call to curl() not available");
1116
1117 return m_basis.curl(i + m_shift, x);
1118 }
1119
1121 inline DivergenceValue divergence(size_t i, const VectorRd &x) const
1122 {
1123 static_assert(hasDivergence, "Call to divergence() not available");
1124
1125 return m_basis.divergence(i + m_shift, x);
1126 }
1127
1129 inline HessianValue hessian(size_t i, const VectorRd &x) const
1130 {
1131 static_assert(hasHessian, "Call to hessian() not available");
1132
1133 return m_basis.hessian(i + m_shift, x);
1134 }
1135
1136 private:
1137 BasisType m_basis;
1138 int m_shift;
1139 };
1140
1141 //-----------------------RESTRICTED BASIS-------------------------------------------------------
1142
1144
1145 template <typename BasisType>
1147 {
1148 public:
1149 typedef typename BasisType::FunctionValue FunctionValue;
1150 typedef typename BasisType::GradientValue GradientValue;
1151 typedef typename BasisType::CurlValue CurlValue;
1152 typedef typename BasisType::DivergenceValue DivergenceValue;
1153 typedef typename BasisType::HessianValue HessianValue;
1154
1155 typedef typename BasisType::GeometricSupport GeometricSupport;
1156
1157 constexpr static const TensorRankE tensorRank = BasisType::tensorRank;
1158 constexpr static const bool hasAncestor = true;
1159 static const bool hasFunction = BasisType::hasFunction;
1160 static const bool hasGradient = BasisType::hasGradient;
1161 static const bool hasCurl = BasisType::hasCurl;
1162 static const bool hasDivergence = BasisType::hasDivergence;
1163 static const bool hasHessian = BasisType::hasHessian;
1164 static const bool hasCurlCurl = BasisType::hasCurlCurl;
1165
1167
1170 const BasisType &basis,
1171 const size_t &dimension
1172 )
1173 : m_basis(basis),
1174 m_dimension(dimension)
1175 {
1176 // Make sure that the provided dimension is smaller than the one
1177 // of the provided basis
1178 assert(dimension <= basis.dimension());
1179 }
1180
1182 inline size_t dimension() const
1183 {
1184 return m_dimension;
1185 }
1186
1188 constexpr inline const BasisType &ancestor() const
1189 {
1190 return m_basis;
1191 }
1192
1194 inline size_t max_degree() const
1195 {
1196 // We need to find the degree based on the dimension, assumes the basis is hierarchical
1197 size_t deg = 0;
1198 while (PolynomialSpaceDimension<GeometricSupport>::Poly(deg) < m_dimension){
1199 deg++;
1200 }
1201 return deg;
1202 }
1203
1205 inline FunctionValue function(size_t i, const VectorRd &x) const
1206 {
1207 static_assert(hasFunction, "Call to function() not available");
1208
1209 return m_basis.function(i, x);
1210 }
1211
1213 GradientValue gradient(size_t i, const VectorRd &x) const
1214 {
1215 static_assert(hasGradient, "Call to gradient() not available");
1216
1217 return m_basis.gradient(i, x);
1218 }
1219
1221 CurlValue curl(size_t i, const VectorRd &x) const
1222 {
1223 static_assert(hasCurl, "Call to curl() not available");
1224
1225 return m_basis.curl(i, x);
1226 }
1227
1229 DivergenceValue divergence(size_t i, const VectorRd &x) const
1230 {
1231 static_assert(hasDivergence, "Call to divergence() not available");
1232
1233 return m_basis.divergence(i, x);
1234 }
1235
1237 HessianValue hessian(size_t i, const VectorRd &x) const
1238 {
1239 static_assert(hasHessian, "Call to hessian() not available");
1240
1241 return m_basis.hessian(i, x);
1242 }
1243
1244 private:
1245 BasisType m_basis;
1246 size_t m_dimension;
1247 };
1248
1249 //--------------------GRADIENT BASIS----------------------------------------------------------
1250
1252
1254 template <typename BasisType>
1256 {
1257 public:
1258 typedef typename BasisType::GradientValue FunctionValue;
1259 typedef void GradientValue;
1260 typedef void CurlValue;
1261 typedef void DivergenceValue;
1262 typedef void HessianValue;
1263
1264 typedef typename BasisType::GeometricSupport GeometricSupport;
1265
1266 constexpr static const TensorRankE tensorRank = (BasisType::tensorRank==Scalar ? Vector : Matrix);
1267 constexpr static const bool hasAncestor = true;
1268 static const bool hasFunction = true;
1269 static const bool hasGradient = false;
1270 static const bool hasCurl = false;
1271 static const bool hasDivergence = false;
1272 static const bool hasHessian = false;
1273 static const bool hasCurlCurl = false;
1274
1276
1279 : m_scalar_basis(basis)
1280 {
1281 static_assert(BasisType::tensorRank == Scalar || BasisType::tensorRank == Vector,
1282 "Gradient basis can only be constructed starting from scalar or vector bases");
1283 static_assert(BasisType::hasGradient,
1284 "Gradient basis requires gradient() for the original basis to be available");
1285 // Do nothing
1286 }
1287
1289 inline size_t dimension() const
1290 {
1291 return m_scalar_basis.dimension();
1292 }
1293
1295 inline FunctionValue function(size_t i, const VectorRd &x) const
1296 {
1297 return m_scalar_basis.gradient(i, x);
1298 }
1299
1301 constexpr inline const BasisType &ancestor() const
1302 {
1303 return m_scalar_basis;
1304 }
1305
1306
1307 private:
1308 BasisType m_scalar_basis;
1309 };
1310
1311 //-------------------CURL BASIS-----------------------------------------------------------
1312
1314
1316 template <typename BasisType>
1318 {
1319 public:
1320 typedef VectorRd FunctionValue;
1321 typedef Eigen::Matrix<double, dimspace, dimspace> GradientValue;
1322 typedef Eigen::Matrix<double, dimspace, dimspace> CurlValue;
1323 typedef double DivergenceValue;
1324 typedef Eigen::Matrix<double, dimspace, dimspace*dimspace> HessianValue;
1325
1326 typedef typename BasisType::GeometricSupport GeometricSupport;
1327
1328 constexpr static const TensorRankE tensorRank = Vector;
1329 constexpr static const bool hasAncestor = true;
1330 static const bool hasFunction = true;
1331 static const bool hasGradient = false;
1332 static const bool hasCurl = false;
1333 static const bool hasDivergence = false;
1334 static const bool hasHessian = false;
1335 static const bool hasCurlCurl = false;
1336
1338
1340 CurlBasis(const BasisType &basis)
1341 : m_basis(basis)
1342 {
1343 static_assert((BasisType::tensorRank == Vector && std::is_same<typename BasisType::GeometricSupport, Cell>::value) ||
1344 (BasisType::tensorRank == Scalar && std::is_same<typename BasisType::GeometricSupport, Face>::value),
1345 "Curl basis can only be constructed starting from vector bases on elements or scalar bases on faces");
1346 static_assert(BasisType::hasCurl,
1347 "Curl basis requires curl() for the original basis to be available");
1348 }
1349
1351 inline size_t dimension() const
1352 {
1353 return m_basis.dimension();
1354 }
1355
1357 inline FunctionValue function(size_t i, const VectorRd &x) const
1358 {
1359 return m_basis.curl(i, x);
1360 }
1361
1363 constexpr inline const BasisType &ancestor() const
1364 {
1365 return m_basis;
1366 }
1367
1368
1369 private:
1370 size_t m_degree;
1371 BasisType m_basis;
1372 };
1373
1374
1375 //--------------------DIVERGENCE BASIS----------------------------------------------------------
1376
1378
1379 template <typename BasisType>
1381 {
1382 public:
1383 typedef typename BasisType::DivergenceValue FunctionValue;
1384 typedef void GradientValue;
1385 typedef void CurlValue;
1386 typedef void DivergenceValue;
1387 typedef void HessianValue;
1388
1389 typedef typename BasisType::GeometricSupport GeometricSupport;
1390
1391 constexpr static const TensorRankE tensorRank = Scalar;
1392 constexpr static const bool hasAncestor = true;
1393 static const bool hasFunction = true;
1394 static const bool hasGradient = false;
1395 static const bool hasCurl = false;
1396 static const bool hasDivergence = false;
1397 static const bool hasHessian = false;
1398 static const bool hasCurlCurl = false;
1399
1401
1404 : m_vector_basis(basis)
1405 {
1406 static_assert(BasisType::tensorRank == Vector || BasisType::tensorRank == Matrix,
1407 "Divergence basis can only be constructed starting from vector bases");
1408 static_assert(BasisType::hasDivergence,
1409 "Divergence basis requires divergence() for the original basis to be available");
1410 // Do nothing
1411 }
1412
1414 inline size_t dimension() const
1415 {
1416 return m_vector_basis.dimension();
1417 }
1418
1420 inline FunctionValue function(size_t i, const VectorRd &x) const
1421 {
1422 return m_vector_basis.divergence(i, x);
1423 }
1424
1426 constexpr inline const BasisType &ancestor() const
1427 {
1428 return m_vector_basis;
1429 }
1430
1431
1432 private:
1433 BasisType m_vector_basis;
1434 };
1435
1436
1438 /* It is created as the DivergenceBasis of a TangentFamily obtained rotating the generators by -pi/2 in the oriented face */
1439 template <typename ScalarFamilyType>
1441 {
1442 // Take the generator and rotates them
1443 Eigen::MatrixXd gen = tangent_family.generators();
1444 VectorRd nF = F.normal();
1445 Eigen::MatrixXd rotated_gen = Eigen::MatrixXd::Zero(2, 3);
1446 rotated_gen.row(0) = VectorRd(gen.row(0)).cross(nF);
1447 rotated_gen.row(1) = VectorRd(gen.row(1)).cross(nF);
1448
1450
1452 }
1453
1454 //---------------------------------------------------------------------
1455 // BASES FOR THE KOSZUL COMPLEMENTS OF G^k, R^k
1456 //---------------------------------------------------------------------
1457
1460 {
1461 public:
1462 typedef VectorRd FunctionValue;
1463 typedef Eigen::Matrix<double, dimspace, dimspace> GradientValue;
1464 typedef Eigen::Matrix<double, dimspace, dimspace> CurlValue;
1465 typedef double DivergenceValue;
1466 typedef Eigen::Matrix<double, dimspace, dimspace*dimspace> HessianValue;
1467
1468 typedef Cell GeometricSupport;
1469
1470 constexpr static const TensorRankE tensorRank = Vector;
1471 constexpr static const bool hasAncestor = false;
1472 static const bool hasFunction = true;
1473 static const bool hasGradient = false;
1474 static const bool hasCurl = false;
1475 static const bool hasDivergence = true;
1476 static const bool hasHessian = false;
1477 static const bool hasCurlCurl = false;
1478
1481 const Cell &T,
1482 size_t degree
1483 );
1484
1486 inline size_t dimension() const
1487 {
1489 }
1490
1492 FunctionValue function(size_t i, const VectorRd &x) const;
1493
1495 DivergenceValue divergence(size_t i, const VectorRd &x) const;
1496
1498 inline size_t max_degree() const
1499 {
1500 return m_degree;
1501 }
1502
1504 inline VectorZd powers(size_t i) const
1505 {
1506 return m_powers[i];
1507 }
1508
1509 private:
1511 inline VectorRd _coordinate_transform(const VectorRd &x) const
1512 {
1513 return (x - m_xT) / m_hT;
1514 }
1515
1516 size_t m_degree;
1517 VectorRd m_xT;
1518 double m_hT;
1519 std::vector<VectorZd> m_powers;
1520 };
1521
1523 // This basis is obtained by translation and scaling of the following basis of x \times (P^{k-1}(R^3))^3
1524 // (suggested by Daniel Mathews (daniel.mathews@monash.edu)):
1525 // {scalar monomials in (x1,x2,x3) of degree <= k-1} * {(0, x3, -x2), (-x3, 0, x1)}
1526 // and
1527 // {scalar monomials in (x1,x2) of degree <=k-1} * (x2, -x1, 0)
1528 // The vectors (0,x3,-x2), (-x3,0,x1) and (x2,-x1,0) are the "directions", the scalar factors in P^{k-1}
1529 // are given by m_powers
1531 {
1532 public:
1533 typedef VectorRd FunctionValue;
1534 typedef Eigen::Matrix<double, dimspace, dimspace> GradientValue;
1535 typedef VectorRd CurlValue;
1536 typedef double DivergenceValue;
1537 typedef Eigen::Matrix<double, dimspace, dimspace*dimspace> HessianValue;
1538
1539 typedef Cell GeometricSupport;
1540
1541 constexpr static const TensorRankE tensorRank = Vector;
1542 constexpr static const bool hasAncestor = false;
1543 static const bool hasFunction = true;
1544 static const bool hasGradient = false;
1545 static const bool hasCurl = true;
1546 static const bool hasDivergence = false;
1547 static const bool hasHessian = false;
1548 static const bool hasCurlCurl = false;
1549
1552 const Cell &T,
1553 size_t degree
1554 );
1555
1557 inline size_t dimension() const
1558 {
1560 }
1561
1563 FunctionValue function(size_t i, const VectorRd &x) const;
1564
1566 CurlValue curl(size_t i, const VectorRd &x) const;
1567
1569 inline size_t max_degree() const
1570 {
1571 return m_degree;
1572 }
1573
1575 inline VectorZd powers(size_t i) const
1576 {
1577 return m_powers[i];
1578 }
1579
1581 inline size_t dimPkmo() const
1582 {
1583 return m_dimPkmo3D;
1584 }
1585
1586 private:
1588 inline VectorRd _coordinate_transform(const VectorRd &x) const
1589 {
1590 return (x - m_xT) / m_hT;
1591 }
1592
1593 // To evaluate the value and curl of the "direction"
1594 VectorRd direction_value(size_t i, const VectorRd &x) const;
1595 VectorRd direction_curl(size_t i, const VectorRd &x) const;
1596
1597 size_t m_degree; // Degree k of the basis
1598 VectorRd m_xT; // center of cell
1599 double m_hT; // diameter of cell
1600 size_t m_dimPkmo3D; // shorthand for dimension of P^{k-1}(R^3)
1601 size_t m_dimPkmo2D; // shorthand for dimension of P^{k-1}(R^2)
1602 std::vector<VectorZd> m_powers; // vector listing the powers of the scalar factors in the basis
1603
1604 };
1605
1608 {
1609 public:
1610 typedef VectorRd FunctionValue;
1611 typedef Eigen::Matrix<double, dimspace, dimspace> GradientValue;
1612 typedef Eigen::Matrix<double, dimspace, dimspace> CurlValue;
1613 typedef double DivergenceValue;
1614 typedef Eigen::Matrix<double, dimspace, dimspace*dimspace> HessianValue;
1615
1616 typedef Face GeometricSupport;
1617
1618 typedef Eigen::Matrix<double, 2, dimspace> JacobianType;
1619
1620 constexpr static const TensorRankE tensorRank = Vector;
1621 constexpr static const bool hasAncestor = false;
1622 static const bool hasFunction = true;
1623 static const bool hasGradient = false;
1624 static const bool hasCurl = false;
1625 static const bool hasDivergence = true;
1626 static const bool hasHessian = false;
1627 static const bool hasCurlCurl = false;
1628
1631 const Face &F,
1632 size_t degree
1633 );
1634
1636 inline size_t dimension() const
1637 {
1639 }
1640
1642 FunctionValue function(size_t i, const VectorRd &x) const;
1643
1645 DivergenceValue divergence(size_t i, const VectorRd &x) const;
1646
1648 inline const JacobianType &jacobian() const
1649 {
1650 return m_jacobian;
1651 }
1652
1654 inline size_t max_degree() const
1655 {
1656 return m_degree;
1657 }
1658
1660 inline Eigen::Vector2i powers(size_t i) const
1661 {
1662 return m_powers[i];
1663 }
1664
1665 private:
1667 inline Eigen::Vector2d _coordinate_transform(const VectorRd &x) const
1668 {
1669 return m_jacobian * (x - m_xF);
1670 }
1671
1672 size_t m_degree;
1673 VectorRd m_xF;
1674 double m_hF;
1675 JacobianType m_jacobian;
1676 std::vector<Eigen::Vector2i> m_powers;
1677 };
1678
1679
1682 {
1683 public:
1684 typedef VectorRd FunctionValue;
1685 typedef Eigen::Matrix<double, dimspace, dimspace> GradientValue;
1686 typedef Eigen::Matrix<double, dimspace, dimspace> CurlValue;
1687 typedef double DivergenceValue;
1688 typedef Eigen::Matrix<double, dimspace, dimspace*dimspace> HessianValue;
1689
1690 typedef Face GeometricSupport;
1691
1692 typedef Eigen::Matrix<double, 2, dimspace> JacobianType;
1693
1694 constexpr static const TensorRankE tensorRank = Vector;
1695 constexpr static const bool hasAncestor = false;
1696 static const bool hasFunction = true;
1697 static const bool hasGradient = false;
1698 static const bool hasCurl = false;
1699 static const bool hasDivergence = false;
1700 static const bool hasHessian = false;
1701 static const bool hasCurlCurl = false;
1702
1705 const Face &F,
1706 size_t degree
1707 );
1708
1710 inline size_t dimension() const
1711 {
1713 }
1714
1716 FunctionValue function(size_t i, const VectorRd &x) const;
1717
1719 inline const VectorRd &normal() const
1720 {
1721 return m_nF;
1722 }
1723
1725 inline const std::shared_ptr<RolyComplBasisFace> &rck() const
1726 {
1727 return m_Rck_basis;
1728 }
1729
1730 private:
1731 size_t m_degree;
1732 VectorRd m_nF;
1733 std::shared_ptr<RolyComplBasisFace> m_Rck_basis;
1734 };
1735
1736 //------------------------------------------------------------------------------
1737 // Free functions
1738 //------------------------------------------------------------------------------
1739
1749
1751 template <typename outValue, typename inValue, typename FunctionType>
1752 boost::multi_array<outValue, 2> transform_values_quad(
1753 const boost::multi_array<inValue, 2> &B_quad,
1754 const FunctionType &F
1755 )
1756 {
1757 boost::multi_array<outValue, 2> transformed_B_quad( boost::extents[B_quad.shape()[0]][B_quad.shape()[1]] );
1758
1759 std::transform( B_quad.origin(), B_quad.origin() + B_quad.num_elements(), transformed_B_quad.origin(), F);
1760
1761 return transformed_B_quad;
1762 }
1763
1765
1766 template <typename ScalarBasisType, size_t N>
1768 const ScalarBasisType & B,
1769 const std::vector<Eigen::VectorXd> & v
1770 )
1771 {
1772 size_t r = B.dimension();
1773 size_t k = v.size();
1774
1775 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(r*k, r*N);
1776
1777 for (size_t i=0; i < k; i++){
1778 // Check the dimension of vector v[i]
1779 assert(v[i].size() == N);
1780
1781 // Fill in M
1782 for (size_t j=0; j < N; j++){
1783 M.block(i*r, j*r, r, r) = v[i](j) * Eigen::MatrixXd::Identity(r,r);
1784 }
1785 }
1786
1789 }
1790
1792
1795 Eigen::MatrixXd PermuteTensorization( const size_t a,
1796 const size_t b
1797 );
1798
1800 inline static std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)>
1801 symmetrise_matrix = [](const Eigen::MatrixXd & x)->Eigen::MatrixXd { return 0.5*(x+x.transpose());};
1802
1804 inline static std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)>
1805 skew_symmetrise_matrix = [](const Eigen::MatrixXd & x)->Eigen::MatrixXd { return 0.5*(x-x.transpose());};
1806
1808
1809 template <typename ScalarBasisType, size_t N>
1811 const ScalarBasisType & B
1812 )
1813 {
1814 size_t r = B.dimension();
1815 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(r, r*N*N);
1816
1817 for (size_t k=0; k < N; k++){
1818 M.block(0, k*r*(N+1), r, r) = Eigen::MatrixXd::Identity(r, r);
1819 }
1820
1823 }
1824
1825 //------------------------------------------------------------------------------
1826 // BASIS EVALUATIONS
1827 //------------------------------------------------------------------------------
1828
1829 //-----------------------DETAILS FOR EVALUATIONS----------------------
1830
1831 namespace detail
1832 {
1834
1835 template <typename BasisType, BasisFunctionE BasisFunction>
1837 {
1838 };
1839
1840 // Evaluate the function value at x
1841 template <typename BasisType>
1843 {
1844 static_assert(BasisType::hasFunction, "Call to function not available");
1845 typedef typename BasisType::FunctionValue ReturnValue;
1846 static inline ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
1847 {
1848 return basis.function(i, x);
1849 }
1850
1851 // Computes function value at quadrature node iqn, knowing values of ancestor basis at quadrature nodes
1852 static inline ReturnValue evaluate(
1853 const BasisType &basis,
1854 size_t i,
1855 size_t iqn,
1856 const boost::multi_array<ReturnValue, 2> &ancestor_value_quad
1857 )
1858 {
1859 return basis.function(i, iqn, ancestor_value_quad);
1860 }
1861
1862 };
1863
1864 // Evaluate the gradient value at x
1865 template <typename BasisType>
1867 {
1868 static_assert(BasisType::hasGradient, "Call to gradient not available");
1869 typedef typename BasisType::GradientValue ReturnValue;
1870 static inline ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
1871 {
1872 return basis.gradient(i, x);
1873 }
1874
1875 // Computes gradient value at quadrature node iqn, knowing gradients of ancestor basis at quadrature nodes
1876 static inline ReturnValue evaluate(
1877 const BasisType &basis,
1878 size_t i,
1879 size_t iqn,
1880 const boost::multi_array<ReturnValue, 2> &ancestor_gradient_quad
1881 )
1882 {
1883 return basis.gradient(i, iqn, ancestor_gradient_quad);
1884 }
1885 };
1886
1887 // Evaluate the curl value at x
1888 template <typename BasisType>
1890 {
1891 static_assert(BasisType::hasCurl, "Call to curl not available");
1892 typedef typename BasisType::CurlValue ReturnValue;
1893 static inline ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
1894 {
1895 return basis.curl(i, x);
1896 }
1897
1898 // Computes curl value at quadrature node iqn, knowing curls of ancestor basis at quadrature nodes
1899 static inline ReturnValue evaluate(
1900 const BasisType &basis,
1901 size_t i,
1902 size_t iqn,
1903 const boost::multi_array<ReturnValue, 2> &ancestor_curl_quad
1904 )
1905 {
1906 return basis.curl(i, iqn, ancestor_curl_quad);
1907 }
1908 };
1909
1910 // Evaluate the divergence value at x
1911 template <typename BasisType>
1913 {
1914 static_assert(BasisType::hasDivergence, "Call to divergence not available");
1915 typedef typename BasisType::DivergenceValue ReturnValue;
1916 static inline ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
1917 {
1918 return basis.divergence(i, x);
1919 }
1920
1921 // Computes divergence value at quadrature node iqn, knowing divergences of ancestor basis at quadrature nodes
1922 static inline ReturnValue evaluate(
1923 const BasisType &basis,
1924 size_t i,
1925 size_t iqn,
1926 const boost::multi_array<ReturnValue, 2> &ancestor_divergence_quad
1927 )
1928 {
1929 return basis.divergence(i, iqn, ancestor_divergence_quad);
1930 }
1931
1932 };
1933
1934 // Evaluate the hessian value at x
1935 template <typename BasisType>
1937 {
1938 static_assert(BasisType::hasHessian, "Call to hessian not available");
1939 typedef typename BasisType::HessianValue ReturnValue;
1940 static inline ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
1941 {
1942 return basis.hessian(i, x);
1943 }
1944
1945 // Computes hessian value at quadrature node iqn, knowing hessian of ancestor basis at quadrature nodes
1946 static inline ReturnValue evaluate(
1947 const BasisType &basis,
1948 size_t i,
1949 size_t iqn,
1950 const boost::multi_array<ReturnValue, 2> &ancestor_hessian_quad
1951 )
1952 {
1953 return basis.hessian(i, iqn, ancestor_hessian_quad);
1954 }
1955
1956 };
1957
1958 // Evaluate the curl curl value at x
1959 template <typename BasisType>
1961 {
1962 static_assert(BasisType::hasCurlCurl, "Call to curl curl not available");
1963 typedef typename BasisType::CurlValue ReturnValue;
1964 static inline ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
1965 {
1966 return basis.curlcurl(i, x);
1967 }
1968
1969 // Computes hessian value at quadrature node iqn, knowing hessian of ancestor basis at quadrature nodes
1970 static inline ReturnValue evaluate(
1971 const BasisType &basis,
1972 size_t i,
1973 size_t iqn,
1974 const boost::multi_array<ReturnValue, 2> &ancestor_curlcurl_quad
1975 )
1976 {
1977 return basis.curlcurl(i, iqn, ancestor_curlcurl_quad);
1978 }
1979
1980 };
1981
1982 //---------- TENSORIZED FAMILY EVALUATION TRAITS -----------------------------------------
1983 // Evaluate the value at x of a TensorizedVectorFamily (includes information on the ancestor basis, to optimise eval_quad for tensorized bases)
1984 template <typename ScalarBasisType, size_t N>
1986 {
1987 static_assert(TensorizedVectorFamily<ScalarBasisType, N>::hasFunction, "Call to function not available");
1988
1990 static const BasisFunctionE AncestorBasisFunction = Function; // Type of values needed from ancestor basis
1992
1993 // Computes function values at x
1994 static inline ReturnValue evaluate(
1996 size_t i,
1997 const VectorRd &x
1998 )
1999 {
2000 return basis.function(i, x);
2001 }
2002
2003 // Computes function value at quadrature node iqn, knowing ancestor basis at quadrature nodes
2004 static inline ReturnValue evaluate(
2006 size_t i,
2007 size_t iqn,
2008 const boost::multi_array<double, 2> &ancestor_basis_quad
2009 )
2010 {
2011 return basis.function(i, iqn, ancestor_basis_quad);
2012 }
2013 };
2014
2015 // Evaluate the gradient at x of a TensorizedVectorFamily (includes information on the ancestor basis, to optimise eval_quad for tensorized bases)
2016 template <typename ScalarBasisType, size_t N>
2018 {
2019 static_assert(TensorizedVectorFamily<ScalarBasisType, N>::hasGradient, "Call to gradient not available");
2020
2022 static const BasisFunctionE AncestorBasisFunction = Gradient; // Type of values needed from ancestor basis
2024
2025 // Computes gradient value at x
2026 static inline ReturnValue evaluate(
2028 size_t i,
2029 const VectorRd &x
2030 )
2031 {
2032 return basis.gradient(i, x);
2033 }
2034
2035 // Computes gradient value at quadrature node iqn, knowing ancestor basis at quadrature nodes
2036 static inline ReturnValue evaluate(
2038 size_t i,
2039 size_t iqn,
2040 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2041 )
2042 {
2043 return basis.gradient(i, iqn, ancestor_basis_quad);
2044 }
2045 };
2046
2047 // Evaluate the curl at x of a TensorizedVectorFamily (includes information on the ancestor basis, to optimise eval_quad for tensorized bases)
2048 template <typename ScalarBasisType, size_t N>
2050 {
2051 static_assert(TensorizedVectorFamily<ScalarBasisType, N>::hasCurl, "Call to curl not available");
2052
2054 static const BasisFunctionE AncestorBasisFunction = Gradient; // Type of values needed from ancestor basis
2056
2057 // Computes curl values at x
2058 static inline ReturnValue evaluate(
2060 size_t i,
2061 const VectorRd &x
2062 )
2063 {
2064 return basis.curl(i, x);
2065 }
2066
2067 // Computes curl values at quadrature node iqn, knowing ancestor basis at quadrature nodes
2068 static inline ReturnValue evaluate(
2070 size_t i,
2071 size_t iqn,
2072 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2073 )
2074 {
2075 return basis.curl(i, iqn, ancestor_basis_quad);
2076 }
2077
2078 };
2079
2080 // Evaluate the divergence at x of a TensorizedVectorFamily (includes information on the ancestor basis, to optimise eval_quad for tensorized bases)
2081 template <typename ScalarBasisType, size_t N>
2083 {
2084 static_assert(TensorizedVectorFamily<ScalarBasisType, N>::hasDivergence, "Call to divergence not available");
2085
2087 static const BasisFunctionE AncestorBasisFunction = Gradient; // Type of values needed from ancestor basis
2089
2090 // Computes divergence value at x
2091 static inline ReturnValue evaluate(
2093 size_t i,
2094 const VectorRd &x
2095 )
2096 {
2097 return basis.divergence(i, x);
2098 }
2099
2100 // Computes divergence value at quadrature node iqn, knowing ancestor basis at quadrature nodes
2101 static inline ReturnValue evaluate(
2103 size_t i,
2104 size_t iqn,
2105 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2106 )
2107 {
2108 return basis.divergence(i, iqn, ancestor_basis_quad);
2109 }
2110
2111 };
2112
2113 // Evaluate the curl curl at x of a TensorizedVectorFamily (includes information on the ancestor basis, to optimise eval_quad for tensorized bases)
2114 template <typename ScalarBasisType, size_t N>
2116 {
2117 static_assert(TensorizedVectorFamily<ScalarBasisType, N>::hasCurlCurl, "Call to curl curl not available");
2118
2120 static const BasisFunctionE AncestorBasisFunction = Hessian; // Type of values needed from ancestor basis
2122
2123 // Computes curl curl value at x
2124 static inline ReturnValue evaluate(
2126 size_t i,
2127 const VectorRd &x
2128 )
2129 {
2130 return basis.curlcurl(i, x);
2131 }
2132
2133 // Computes curl curl value at quadrature node iqn, knowing ancestor basis at quadrature nodes
2134 static inline ReturnValue evaluate(
2136 size_t i,
2137 size_t iqn,
2138 const boost::multi_array<MatrixRd, 2> &ancestor_basis_quad
2139 )
2140 {
2141 return basis.curlcurl(i, iqn, ancestor_basis_quad);
2142 }
2143
2144 };
2145
2146 //---------- MATRIX FAMILY EVALUATION TRAITS -----------------------------------------
2147 // Evaluate the value at x of a MatrixFamily (includes information on the ancestor basis, to optimise eval_quad for matrix families)
2148 template <typename ScalarBasisType, size_t N>
2150 {
2151 static_assert(MatrixFamily<ScalarBasisType, N>::hasFunction, "Call to function not available");
2152
2154 static const BasisFunctionE AncestorBasisFunction = Function; // Type of values needed from ancestor basis
2156
2157 // Computes function values at x
2158 static inline ReturnValue evaluate(
2160 size_t i,
2161 const VectorRd &x
2162 )
2163 {
2164 return basis.function(i, x);
2165 }
2166
2167 // Computes function value at quadrature node iqn, knowing ancestor basis at quadrature nodes
2168 static inline ReturnValue evaluate(
2170 size_t i,
2171 size_t iqn,
2172 const boost::multi_array<double, 2> &ancestor_basis_quad
2173 )
2174 {
2175 return basis.function(i, iqn, ancestor_basis_quad);
2176 }
2177 };
2178
2179 // Evaluate the divergence at x of a MatrixFamily (includes information on the ancestor basis, to optimise eval_quad for matrix families)
2180 template <typename ScalarBasisType, size_t N>
2182 {
2183 static_assert(MatrixFamily<ScalarBasisType, N>::hasDivergence, "Call to divergence not available");
2184
2186 static const BasisFunctionE AncestorBasisFunction = Gradient; // Type of values needed from ancestor basis
2188
2189 // Computes divergence value at x
2190 static inline ReturnValue evaluate(
2192 size_t i,
2193 const VectorRd &x
2194 )
2195 {
2196 return basis.divergence(i, x);
2197 }
2198
2199 // Computes divergence value at quadrature node iqn, knowing ancestor basis at quadrature nodes
2200 static inline ReturnValue evaluate(
2202 size_t i,
2203 size_t iqn,
2204 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2205 )
2206 {
2207 return basis.divergence(i, iqn, ancestor_basis_quad);
2208 }
2209
2210 };
2211
2212 } // end of namespace detail
2213
2214
2215 //-----------------------------------EVALUATE_QUAD--------------------------------------
2216
2218 template <BasisFunctionE BasisFunction>
2220 {
2222 template <typename BasisType>
2223 static boost::multi_array<typename detail::basis_evaluation_traits<BasisType, BasisFunction>::ReturnValue, 2>
2225 const BasisType &basis,
2226 const QuadratureRule &quad
2227 )
2228 {
2230
2231 boost::multi_array<typename traits::ReturnValue, 2>
2232 basis_quad(boost::extents[basis.dimension()][quad.size()]);
2233
2234 for (size_t i = 0; i < basis.dimension(); i++)
2235 {
2236 for (size_t iqn = 0; iqn < quad.size(); iqn++)
2237 {
2238 basis_quad[i][iqn] = traits::evaluate(basis, i, quad[iqn].vector());
2239 } // for iqn
2240 } // for i
2241
2242 return basis_quad;
2243 }
2244
2246 template <typename BasisType>
2247 static boost::multi_array<typename detail::basis_evaluation_traits<Family<BasisType>, BasisFunction>::ReturnValue, 2>
2249 const Family<BasisType> &basis,
2250 const QuadratureRule &quad
2251 )
2252 {
2254
2255 boost::multi_array<typename traits::ReturnValue, 2>
2256 basis_quad(boost::extents[basis.dimension()][quad.size()]);
2257
2258 // Compute values at quadrature note on ancestor basis, and then do a transformation
2259 const auto &ancestor_basis = basis.ancestor();
2260 boost::multi_array<typename traits::ReturnValue, 2>
2262
2263 for (size_t iqn = 0; iqn < quad.size(); iqn++)
2264 {
2265 for (size_t i = 0; i < size_t(basis.matrix().rows()); i++){
2266 basis_quad[i][iqn] = traits::evaluate(basis, i, iqn, ancestor_basis_quad);
2267 }
2268 }
2269 return basis_quad;
2270 }
2271
2273 template <typename BasisType, size_t N>
2274 static boost::multi_array<typename detail::basis_evaluation_traits<TensorizedVectorFamily<BasisType, N>, BasisFunction>::ReturnValue, 2>
2277 const QuadratureRule &quad
2278 )
2279 {
2281
2282 boost::multi_array<typename traits::ReturnValue, 2>
2283 basis_quad(boost::extents[basis.dimension()][quad.size()]);
2284
2285 const auto &ancestor_basis = basis.ancestor();
2286 const boost::multi_array<typename traits::AncestorBasisFunctionValue, 2> ancestor_basis_quad
2288
2289 for (size_t i = 0; i < basis.dimension(); i++){
2290 for (size_t iqn = 0; iqn < quad.size(); iqn++){
2291 basis_quad[i][iqn] = traits::evaluate(basis, i, iqn, ancestor_basis_quad);
2292 }
2293 }
2294 return basis_quad;
2295 }
2296
2298 template <typename BasisType, size_t N>
2299 static boost::multi_array<typename detail::basis_evaluation_traits<MatrixFamily<BasisType, N>, BasisFunction>::ReturnValue, 2>
2301 const MatrixFamily<BasisType, N> &basis,
2302 const QuadratureRule &quad
2303 )
2304 {
2306
2307 boost::multi_array<typename traits::ReturnValue, 2>
2308 basis_quad(boost::extents[basis.dimension()][quad.size()]);
2309
2310 const auto &ancestor_basis = basis.ancestor();
2311 const boost::multi_array<typename traits::AncestorBasisFunctionValue, 2> ancestor_basis_quad
2313
2314 for (size_t i = 0; i < basis.dimension(); i++){
2315 for (size_t iqn = 0; iqn < quad.size(); iqn++){
2316 basis_quad[i][iqn] = traits::evaluate(basis, i, iqn, ancestor_basis_quad);
2317 }
2318 }
2319 return basis_quad;
2320 }
2321
2322 };
2323
2324
2325 //------------------------------------------------------------------------------
2326 // ORTHONORMALISATION
2327 //------------------------------------------------------------------------------
2328
2337 template <typename T>
2338 Eigen::MatrixXd gram_schmidt(
2339 boost::multi_array<T, 2> &basis_eval,
2340 const std::function<double(size_t, size_t)> &inner_product
2341 )
2342 {
2343 auto induced_norm = [inner_product](size_t i) -> double {
2344 return std::sqrt(inner_product(i, i));
2345 };
2346
2347 // Number of basis functions
2348 size_t Nb = basis_eval.shape()[0];
2349 // Number of quadrature nodes
2350 size_t Nqn = basis_eval.shape()[1];
2351
2352 Eigen::MatrixXd B = Eigen::MatrixXd::Zero(Nb, Nb);
2353
2354 // Normalise the first element
2355 double norm = induced_norm(0);
2356 for (size_t iqn = 0; iqn < Nqn; iqn++)
2357 {
2358 basis_eval[0][iqn] /= norm;
2359 } // for iqn
2360 B(0, 0) = 1. / norm;
2361
2362 for (size_t ib = 1; ib < Nb; ib++)
2363 {
2364 // 'coeffs' represents the coefficients of the ib-th orthogonal function on the ON basis functions 0 to ib-1
2365 Eigen::RowVectorXd coeffs = Eigen::RowVectorXd::Zero(ib);
2366 for (size_t pb = 0; pb < ib; pb++)
2367 {
2368 coeffs(pb) = -inner_product(ib, pb);
2369 } // for pb
2370
2371 // store the values of the orthogonalised ib-th basis function
2372 for (size_t pb = 0; pb < ib; pb++)
2373 {
2374 for (size_t iqn = 0; iqn < Nqn; iqn++)
2375 {
2377 } // for iqn
2378 } // for pb
2379
2380 // normalise ib-th basis function
2381 double norm = induced_norm(ib);
2382 for (size_t iqn = 0; iqn < Nqn; iqn++)
2383 {
2384 basis_eval[ib][iqn] /= norm;
2385 } // for iqn
2386 coeffs /= norm;
2387 // Compute ib-th row of B.
2388 // B.topLeftCorner(ib, ib) contains the rows that represent the ON basis functions 0 to ib-1 on
2389 // the original basis. Multiplying on the left by coeffs gives the coefficients of the ON basis
2390 // functions on the original basis functions from 0 to ib-1.
2391 B.block(ib, 0, 1, ib) = coeffs * B.topLeftCorner(ib, ib);
2392 B(ib, ib) = 1. / norm;
2393 }
2394 return B;
2395 }
2396
2397 //------------------------------------------------------------------------------
2398
2400 double scalar_product(const double &x, const double &y);
2401
2403 double scalar_product(const double &x, const Eigen::Matrix<double, 1, 1> &y);
2404
2406 double scalar_product(const VectorRd &x, const VectorRd &y);
2407
2409 template<int N>
2410 double scalar_product(const Eigen::Matrix<double, N, N> & x, const Eigen::Matrix<double, N, N> & y)
2411 {
2412 return (x.transpose() * y).trace();
2413 }
2414
2416 template <typename Value>
2417 boost::multi_array<double, 2> scalar_product(
2418 const boost::multi_array<Value, 2> &basis_quad,
2419 const Value &v
2420 )
2421 {
2422 boost::multi_array<double, 2> basis_dot_v_quad(boost::extents[basis_quad.shape()[0]][basis_quad.shape()[1]]);
2423 std::transform(basis_quad.origin(), basis_quad.origin() + basis_quad.num_elements(),
2424 basis_dot_v_quad.origin(), [&v](const Value &x) -> double { return scalar_product(x,v); });
2425 return basis_dot_v_quad;
2426 }
2427
2428
2430 boost::multi_array<VectorRd, 2>
2432 const boost::multi_array<VectorRd, 2> &basis_quad,
2433 const VectorRd &v
2434 );
2435
2437 template <typename BasisType>
2439 const BasisType &basis,
2440 const QuadratureRule &qr,
2441 boost::multi_array<typename BasisType::FunctionValue, 2> &basis_quad
2442 )
2443 {
2444 // Check that the basis evaluation and quadrature rule are coherent
2445 assert(basis.dimension() == basis_quad.shape()[0] && qr.size() == basis_quad.shape()[1]);
2446
2447 // The inner product between the i-th and j-th basis vectors
2448 std::function<double(size_t, size_t)> inner_product = [&basis_quad, &qr](size_t i, size_t j) -> double {
2449 double r = 0.;
2450 for (size_t iqn = 0; iqn < qr.size(); iqn++)
2451 {
2453 } // for iqn
2454 return r;
2455 };
2456
2457 Eigen::MatrixXd B = gram_schmidt(basis_quad, inner_product);
2458
2459 return Family<BasisType>(basis, B);
2460 }
2461
2463 /* This method is more robust that the previous one based on values at quadrature nodes */
2464 template <typename BasisType>
2466 const BasisType &basis,
2467 const Eigen::MatrixXd & GM
2468 )
2469 {
2470 // Check that the basis and Gram matrix are coherent
2471 assert(basis.dimension() == size_t(GM.rows()) && GM.rows() == GM.cols());
2472
2473 Eigen::MatrixXd L = GM.llt().matrixL();
2474
2475 return Family<BasisType>(basis, L.inverse());
2476 }
2477
2478 //------------------------------------------------------------------------------
2479 // Gram matrices
2480 //------------------------------------------------------------------------------
2481
2485 template <typename FunctionValue>
2486 Eigen::MatrixXd compute_gram_matrix(const boost::multi_array<FunctionValue, 2> &B1,
2487 const boost::multi_array<FunctionValue, 2> &B2,
2488 const QuadratureRule &qr,
2489 const size_t nrows,
2490 const size_t ncols,
2491 const std::string sym = "nonsym"
2492 )
2493 {
2494 // Check that the basis evaluation and quadrature rule are coherent
2495 assert(qr.size() == B1.shape()[1] && qr.size() == B2.shape()[1]);
2496 // Check that we don't ask for more members of family than available
2497 assert(nrows <= B1.shape()[0] && ncols <= B2.shape()[0]);
2498
2499 // Re-cast quadrature weights into ArrayXd to make computations faster
2500 Eigen::ArrayXd qr_weights = Eigen::ArrayXd::Zero(qr.size());
2501 for (size_t iqn = 0; iqn < qr.size(); iqn++)
2502 {
2503 qr_weights(iqn) = qr[iqn].w;
2504 }
2505
2506 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(nrows, ncols);
2507 for (size_t i = 0; i < nrows; i++)
2508 {
2509 size_t jcut = 0;
2510 if (sym == "sym")
2511 jcut = i;
2512 for (size_t j = 0; j < jcut; j++)
2513 {
2514 M(i, j) = M(j, i);
2515 }
2516 for (size_t j = jcut; j < ncols; j++)
2517 {
2518 std::vector<double> tmp(B1.shape()[1]);
2519 // Extract values at quadrature nodes for elements i of B1 and j of B2
2520 auto B1i = B1[boost::indices[i][boost::multi_array_types::index_range(0, B1.shape()[1])]];
2521 auto B2j = B2[boost::indices[j][boost::multi_array_types::index_range(0, B1.shape()[1])]];
2522 // Compute scalar product of B1i and B2j, and recast it as ArrayXd
2523 std::transform(B1i.begin(), B1i.end(), B2j.begin(), tmp.begin(), [](FunctionValue a, FunctionValue b) -> double { return scalar_product(a, b); });
2524 Eigen::ArrayXd tmp_array = Eigen::Map<Eigen::ArrayXd, Eigen::Unaligned>(tmp.data(), tmp.size());
2525 // Multiply by quadrature weights and sum (using .sum() of ArrayXd makes this step faster than a loop)
2526 M(i, j) = (qr_weights * tmp_array).sum();
2527
2528 // Simple version with loop (replaces everything above after the for on j)
2529 /*
2530 for (size_t iqn = 0; iqn < qr.size(); iqn++) {
2531 M(i,j) += qr[iqn].w * scalar_product(B1[i][iqn], B2[j][iqn]);
2532 } // for iqn
2533 */
2534
2535 } // for j
2536 } // for i
2537 return M;
2538 }
2539
2543 template <typename FunctionValue>
2544 Eigen::MatrixXd compute_gram_matrix(const boost::multi_array<FunctionValue, 2> &B1,
2545 const boost::multi_array<FunctionValue, 2> &B2,
2546 const QuadratureRule &qr,
2547 const std::string sym = "nonsym"
2548 )
2549 {
2550 return compute_gram_matrix<FunctionValue>(B1, B2, qr, B1.shape()[0], B2.shape()[0], sym);
2551 }
2552
2554 template <typename FunctionValue>
2555 inline Eigen::MatrixXd compute_gram_matrix(const boost::multi_array<FunctionValue, 2> &B,
2556 const QuadratureRule &qr
2557 )
2558 {
2559 return compute_gram_matrix<FunctionValue>(B, B, qr, "sym");
2560 }
2561
2563 Eigen::MatrixXd compute_gram_matrix(const boost::multi_array<VectorRd, 2> &B1,
2564 const boost::multi_array<double, 2> &B2,
2565 const QuadratureRule &qr
2566 );
2567
2571 Eigen::MatrixXd compute_gram_matrix(const boost::multi_array<double, 2> &B1,
2572 const boost::multi_array<double, 2> &B2,
2573 const QuadratureRule &qr,
2574 const size_t nrows,
2575 const size_t ncols,
2576 const std::string sym = "nonsym"
2577 );
2578
2582 Eigen::MatrixXd compute_gram_matrix(const boost::multi_array<double, 2> &B1,
2583 const boost::multi_array<double, 2> &B2,
2584 const QuadratureRule &qr,
2585 const std::string sym = "nonsym"
2586 );
2587
2591 Eigen::MatrixXd compute_gram_matrix(const boost::multi_array<VectorRd, 2> &B1,
2592 const boost::multi_array<VectorRd, 2> &B2,
2593 const QuadratureRule &qr,
2594 const size_t nrows,
2595 const size_t ncols,
2596 const std::string sym = "nonsym"
2597 );
2598
2600 Eigen::MatrixXd compute_gram_matrix(const boost::multi_array<VectorRd, 2> &B1,
2601 const boost::multi_array<VectorRd, 2> &B2,
2602 const QuadratureRule &qr,
2603 const std::string sym = "nonsym"
2604 );
2605
2607 template<typename ScalarFamilyType, size_t N>
2609 const boost::multi_array<double, 2> & scalar_family_quad,
2610 const QuadratureRule & qr
2611 )
2612 {
2613 Eigen::MatrixXd Gram = Eigen::MatrixXd::Zero(MatFam.dimension(), MatFam.dimension());
2614
2615 // Gram matrix of the scalar family
2617 for (size_t i = 0; i < MatFam.dimension(); i = i+scalarGram.rows()){
2618 Gram.block(i, i, scalarGram.rows(), scalarGram.cols()) = scalarGram;
2619 }
2620
2621 return Gram;
2622 }
2623
2625 template <typename T>
2626 Eigen::VectorXd integrate(
2627 const FType<T> &f,
2628 const BasisQuad<T> &B,
2629 const QuadratureRule &qr,
2630 size_t n_rows = 0
2631 )
2632 {
2633 // If default, set n_rows to size of family
2634 if (n_rows == 0)
2635 {
2636 n_rows = B.shape()[0];
2637 }
2638
2639 // Number of quadrature nodes
2640 const size_t num_quads = qr.size();
2641
2642 // Check number of quadrature nodes is compatible with B
2643 assert(num_quads == B.shape()[1]);
2644
2645 // Check that we don't ask for more members of family than available
2646 assert(n_rows <= B.shape()[0]);
2647
2648 Eigen::VectorXd V = Eigen::VectorXd::Zero(n_rows);
2649
2650 for (size_t iqn = 0; iqn < num_quads; iqn++)
2651 {
2652 double qr_weight = qr[iqn].w;
2653 T f_on_qr = f(qr[iqn].vector());
2654 for (size_t i = 0; i < n_rows; i++)
2655 {
2656 V(i) += qr_weight * scalar_product(B[i][iqn], f_on_qr);
2657 }
2658 }
2659
2660 return V;
2661 }
2662
2664 template <typename T, typename U>
2666 const FType<U> &f,
2667 const BasisQuad<T> &B1,
2668 const BasisQuad<T> &B2,
2669 const QuadratureRule &qr,
2670 size_t n_rows = 0,
2671 size_t n_cols = 0,
2672 const std::string sym = "nonsym"
2673 )
2674 {
2675 // If default, set n_rows and n_cols to size of families
2676 if (n_rows == 0 && n_cols == 0)
2677 {
2678 n_rows = B1.shape()[0];
2679 n_cols = B2.shape()[0];
2680 }
2681
2682 // Number of quadrature nodes
2683 const size_t num_quads = qr.size();
2684 // Check number of quadrature nodes is compatible with B1 and B2
2685 assert(num_quads == B1.shape()[1] && num_quads == B2.shape()[1]);
2686 // Check that we don't ask for more members of family than available
2687 assert(n_rows <= B1.shape()[0] && n_cols <= B2.shape()[0]);
2688
2689 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(n_rows, n_cols);
2690 for (size_t iqn = 0; iqn < num_quads; iqn++)
2691 {
2692 double qr_weight = qr[iqn].w;
2693 U f_on_qr = f(qr[iqn].vector());
2694 for (size_t i = 0; i < n_rows; i++)
2695 {
2696 T f_B1 = f_on_qr * B1[i][iqn];
2697 size_t jcut = 0;
2698 if (sym == "sym")
2699 jcut = i;
2700 for (size_t j = 0; j < jcut; j++)
2701 {
2702 M(i, j) = M(j, i);
2703 }
2704 for (size_t j = jcut; j < n_cols; j++)
2705 {
2706 M(i, j) += qr_weight * scalar_product(f_B1, B2[j][iqn]);
2707 }
2708 }
2709 }
2710 return M;
2711 }
2712
2714 template <typename T, typename U>
2716 const FType<U> &f,
2717 const BasisQuad<T> &B1,
2718 const BasisQuad<T> &B2,
2719 const QuadratureRule &qr,
2720 const std::string sym
2721 )
2722 {
2723 return compute_weighted_gram_matrix(f, B1, B2, qr, B1.shape()[0], B2.shape()[0], sym);
2724 }
2725
2727 Eigen::MatrixXd compute_weighted_gram_matrix(
2728 const FType<VectorRd> &f,
2729 const BasisQuad<VectorRd> &B1,
2730 const BasisQuad<double> &B2,
2731 const QuadratureRule &qr,
2732 size_t n_rows = 0,
2733 size_t n_cols = 0
2734 );
2735
2737 Eigen::MatrixXd compute_weighted_gram_matrix(
2738 const FType<VectorRd> &f,
2739 const BasisQuad<double> &B1,
2740 const BasisQuad<VectorRd> &B2,
2741 const QuadratureRule &qr,
2742 size_t n_rows = 0,
2743 size_t n_cols = 0
2744 );
2745
2746 //------------------------------------------------------------------------------
2747 // L2 projection of a function
2748 //------------------------------------------------------------------------------
2749
2751 template <typename BasisType>
2752 Eigen::VectorXd l2_projection(
2753 const std::function<typename BasisType::FunctionValue(const VectorRd &)> &f,
2754 const BasisType &basis,
2756 const boost::multi_array<typename BasisType::FunctionValue, 2> &basis_quad,
2757 const Eigen::MatrixXd &mass_basis = Eigen::MatrixXd::Zero(1,1)
2758 )
2759 {
2760 Eigen::MatrixXd Mass = mass_basis;
2761 if (Mass.norm() < 1e-13){
2763 }
2764 Eigen::LDLT<Eigen::MatrixXd> cholesky_mass(Mass);
2765 Eigen::VectorXd b = Eigen::VectorXd::Zero(basis.dimension());
2766 for (size_t i = 0; i < basis.dimension(); i++)
2767 {
2768 for (size_t iqn = 0; iqn < quad.size(); iqn++)
2769 {
2770 b(i) += quad[iqn].w * scalar_product(f(quad[iqn].vector()), basis_quad[i][iqn]);
2771 } // for iqn
2772 } // for i
2773 return cholesky_mass.solve(b);
2774 }
2775
2776 //------------------------------------------------------------------------------
2777 // Decomposition on a polynomial basis
2778 //------------------------------------------------------------------------------
2779
2781
2786 template <typename BasisType>
2788 {
2790
2796 const Face &F,
2797 const BasisType &basis
2798 ):
2799 m_dim(basis.dimension()),
2800 m_basis(basis),
2801 m_on_basis(basis, Eigen::MatrixXd::Identity(m_dim,m_dim)) // Need to initialise before we compute the real ON basis
2802 {
2803 /* The decomposition will be performed using an orthonormalised version of m_basis, for stability.
2804 The scalar product for which we orthonormalise is a minimal one (much less expensive than integrating
2805 over the face): \sum_i phi(x_i) psi(x_i) where (x_i)_i are nodes that are sufficient to determine
2806 entirely polynomials up to the considered degree.
2807 The nodes we build here correspond to P^k nodes on a triangle on a face (the triangle being centered
2808 at the face center) */
2809
2810 // Create nodes
2811 std::vector<Eigen::Vector2i> indices = MonomialPowers<Face>::complete(basis.max_degree());
2812 VectorRd tF = F.edge(0)->tangent();
2813 VectorRd nF = F.edge_normal(0);
2814 VectorRd x0 = F.center_mass() - F.diam()*(tF+nF)/3.0;
2815
2816 // Nodes
2817 m_nb_nodes = indices.size();
2818 m_nodes.reserve(m_nb_nodes);
2819 for (size_t i=0; i<m_nb_nodes; i++){
2820 VectorRd xi = x0 + F.diam()*(double(indices[i](0)) * tF + double(indices[i](1)) * nF)/std::max(1.,double(basis.max_degree()));
2821 m_nodes.emplace_back(xi.x(), xi.y(), xi.z(), 1.);
2822 }
2823
2824 // We then orthonormalise Orthonormalise basis for simple scalar product
2825 m_on_basis_nodes.resize(boost::extents[m_dim][m_nb_nodes]);
2828 };
2829
2832 const Cell &T,
2833 const BasisType &basis
2834 ):
2835 m_dim(basis.dimension()),
2836 m_basis(basis),
2837 m_on_basis(basis, Eigen::MatrixXd::Identity(m_dim,m_dim)) // Need to initialise before we compute the real ON basis
2838 {
2839 /* The decomposition will be performed using an orthonormalised version of m_basis, for stability.
2840 The scalar product for which we orthonormalise is a minimal one (much less expensive than integrating
2841 over the cell): \sum_i phi(x_i) psi(x_i) where (x_i)_i are nodes that are sufficient to determine
2842 entirely polynomials up to the considered degree.
2843 The nodes we build here correspond to P^k nodes on a tetrahedron */
2844
2845 // Create nodes
2846 std::vector<Eigen::Vector3i> indices = MonomialPowers<Cell>::complete(basis.max_degree());
2847 VectorRd e0(1., 0., 0.);
2848 VectorRd e1(0., 1., 0.);
2849 VectorRd e2(0., 0., 1.);
2850 VectorRd xT = T.center_mass() - T.diam()*(e0+e1+e2)/3.0;
2851
2852 // Nodes
2853 m_nb_nodes = indices.size();
2854 m_nodes.reserve(m_nb_nodes);
2855 for (size_t i=0; i<m_nb_nodes; i++){
2856 VectorRd xi = xT + T.diam()*(double(indices[i](0)) * e0 + double(indices[i](1)) * e1 + double(indices[i](2)) * e2)/std::max(1.,double(basis.max_degree()));
2857 m_nodes.emplace_back(xi.x(), xi.y(), xi.z(), 1.);
2858 }
2859
2860 // We then orthonormalise Orthonormalise basis for simple scalar product
2861 m_on_basis_nodes.resize(boost::extents[m_dim][m_nb_nodes]);
2864 };
2865
2868 {
2869 return m_nodes;
2870 };
2871
2874 boost::multi_array<typename BasisType::FunctionValue, 2> &values
2875 )
2876 {
2877 // Gram matrix of the polynomials and the ON basis
2879 // L=matrix of P as a family of the ON basis. In theory, Gram_onbasis is identity, but due to rounding errors it
2880 // can be a bit different. Computing L as if it's not the case improves stability
2882 Eigen::MatrixXd L = Gram_onbasis.ldlt().solve(Gram_P_onbasis.transpose());
2883 return Family<BasisType>(m_basis, L.transpose() * m_on_basis.matrix());
2884 };
2885
2886
2887 // Member variables
2888 size_t m_dim; // dimension of basis
2889 BasisType m_basis; // Basis on which we decompose
2890 Family<BasisType> m_on_basis; // Orthonormalised basis (for simple scalar product)
2891 boost::multi_array<typename BasisType::FunctionValue, 2> m_on_basis_nodes; // Values of ON basis at the nodes
2892 size_t m_nb_nodes; // nb of nodes
2894
2895 };
2896
2897
2898
2900
2901} // end of namespace HArDCore3D
2902
2903#endif
Basis for the space of curls of polynomials.
Definition basis.hpp:1318
Basis (or rather family) of divergence of an existing basis.
Definition basis.hpp:1381
Family of functions expressed as linear combination of the functions of a given basis.
Definition basis.hpp:389
Basis for the complement G^{c,k}(T) in P^k(T)^3 of the range of grad.
Definition basis.hpp:1531
Basis for the complement G^{c,k}(F) in P^k(F)^2 of the range of the gradient on a face.
Definition basis.hpp:1682
Basis for the space of gradients of polynomials.
Definition basis.hpp:1256
Matrix family obtained from a scalar family.
Definition basis.hpp:778
Scalar monomial basis on a cell.
Definition basis.hpp:155
Scalar monomial basis on an edge.
Definition basis.hpp:320
Scalar monomial basis on a face.
Definition basis.hpp:226
Generate a basis restricted to the first "dimension" functions.
Definition basis.hpp:1147
Basis for the complement R^{c,k}(T) in P^k(T)^3 of the range of curl.
Definition basis.hpp:1460
Basis for the complement R^{c,k}(F) in P^k(F)^2 of the range of the vectorial rotational on a face.
Definition basis.hpp:1608
Generate a basis where the function indices are shifted.
Definition basis.hpp:1040
Vector family for polynomial functions that are tangent to a certain place (determined by the generat...
Definition basis.hpp:948
Vector family obtained by tensorization of a scalar family.
Definition basis.hpp:610
static constexpr const bool hasAncestor
Definition basis.hpp:400
Eigen::Matrix< double, dimspace, dimspace *dimspace > HessianValue
Definition basis.hpp:1466
static constexpr const bool hasAncestor
Definition basis.hpp:239
void GradientValue
Definition basis.hpp:1259
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:1229
double DivergenceValue
Definition basis.hpp:1687
double DivergenceValue
Definition basis.hpp:1536
VectorRd FunctionValue
Definition basis.hpp:1462
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:358
static const bool hasCurlCurl
Definition basis.hpp:1335
Family< BasisType > l2_orthonormalize(const BasisType &basis, const QuadratureRule &qr, boost::multi_array< typename BasisType::FunctionValue, 2 > &basis_quad)
-orthonormalization: simply consists in using gram_schmidt() with the specific l2 inner product
Definition basis.hpp:2438
VectorRd CurlValue
Definition basis.hpp:1535
const Eigen::MatrixXd symmetriseOperator() const
Return the symmetrisation operator, the rN^2 square matrix that to a given vector of coefficients on ...
Definition basis.hpp:883
Eigen::Matrix< double, dimspace, dimspace *dimspace > HessianValue
Definition basis.hpp:1614
DivergenceValue divergence(size_t i, size_t iqn, const boost::multi_array< DivergenceValue, 2 > &ancestor_divergence_quad) const
Evaluate the divergence of the i-th function at a quadrature point iqn, knowing all the divergences o...
Definition basis.hpp:520
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_hessian_quad)
Definition basis.hpp:1946
static const bool hasGradient
Definition basis.hpp:1269
static const bool hasHessian
Definition basis.hpp:1547
size_t dimension() const
Return the dimension of the family.
Definition basis.hpp:984
Eigen::Matrix< double, 2, dimspace > JacobianType
Definition basis.hpp:1692
size_t dimension() const
Return the dimension of the basis.
Definition basis.hpp:1073
static const bool hasCurl
Definition basis.hpp:1474
static constexpr const bool hasAncestor
Definition basis.hpp:1051
double DivergenceValue
Definition basis.hpp:160
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2058
VectorRd FunctionValue
Definition basis.hpp:1610
constexpr const BasisType & ancestor() const
Return the ancestor.
Definition basis.hpp:565
MatrixRd HessianValue
Definition basis.hpp:161
VectorRd FunctionValue
Definition basis.hpp:950
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< double, 2 > &ancestor_basis_quad)
Definition basis.hpp:2168
static const bool hasHessian
Definition basis.hpp:1056
static const bool hasDivergence
Definition basis.hpp:1699
static const bool hasHessian
Definition basis.hpp:1626
void CurlValue
Definition basis.hpp:782
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1301
VectorRd CurlValue
Definition basis.hpp:230
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1157
Eigen::Matrix< double, dimspace, dimspace *dimspace > HessianValue
Definition basis.hpp:1688
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2200
BasisType::FunctionValue FunctionValue
Definition basis.hpp:1149
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the two-dimensional curl of the i-th basis function at point x.
Definition basis.cpp:91
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1870
Face GeometricSupport
Definition basis.hpp:956
static const bool hasHessian
Definition basis.hpp:1476
static const bool hasGradient
Definition basis.hpp:333
static const bool hasCurlCurl
Definition basis.hpp:1164
size_t dimPkmo() const
Returns the dimension of P^{k-1}(R^3)
Definition basis.hpp:1581
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< MatrixRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2134
BasisType::DivergenceValue ReturnValue
Definition basis.hpp:1915
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:1105
static boost::multi_array< typename detail::basis_evaluation_traits< MatrixFamily< BasisType, N >, BasisFunction >::ReturnValue, 2 > compute(const MatrixFamily< BasisType, N > &basis, const QuadratureRule &quad)
Evaluate a Matrix family at quadrature nodes (optimised compared the generic basis evaluation,...
Definition basis.hpp:2300
size_t dimension() const
Return the dimension of the family.
Definition basis.hpp:641
const std::shared_ptr< RolyComplBasisFace > & rck() const
Return the Rck basis.
Definition basis.hpp:1725
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1295
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1050
boost::multi_array< T, 2 > BasisQuad
type for bases evaluated on quadrature nodes
Definition basis.hpp:61
static const bool hasDivergence
Definition basis.hpp:335
constexpr int dimspace
Dimension, and generic types for vector in correct dimension (makes it easier to translate a code bet...
Definition basis.hpp:50
TensorizedVectorFamily< ScalarBasisType, N >::GradientValue ReturnValue
Definition basis.hpp:2021
static const bool hasGradient
Definition basis.hpp:961
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.cpp:81
constexpr const BasisType & ancestor() const
Return the underlying complete basis.
Definition basis.hpp:1085
BasisType AncestorType
Definition basis.hpp:1166
Family< TensorizedVectorFamily< ScalarBasisType, N > > GenericTensorization(const ScalarBasisType &B, const std::vector< Eigen::VectorXd > &v)
From a scalar family B=(B_1..B_r) and vectors (v_1..v_k) in R^N, constructs a "Family" of "Tensorized...
Definition basis.hpp:1767
static std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> symmetrise_matrix
Function to symmetrise a matrix (useful together with transform_values_quad)
Definition basis.hpp:1801
static const bool hasHessian
Definition basis.hpp:1334
DivergenceValue divergence(size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_gradient_quad) const
Evaluate the divergence of the i-th basis function at a quadrature point iqn, knowing all the gradien...
Definition basis.hpp:856
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1498
FunctionValue function(size_t i, size_t iqn, const boost::multi_array< double, 2 > &ancestor_value_quad) const
Evaluate the i-th basis function at a quadrature point iqn, knowing all the values of ancestor basis ...
Definition basis.hpp:657
VectorRd CurlValue
Definition basis.hpp:614
void GradientValue
Definition basis.hpp:1384
Eigen::Matrix< double, dimspace, dimspace *dimspace > HessianValue
Definition basis.hpp:1537
FunctionValue function(size_t i, size_t iqn, const boost::multi_array< FunctionValue, 2 > &ancestor_value_quad) const
Evaluate the i-th function at a quadrature point iqn, knowing all the values of ancestor basis functi...
Definition basis.hpp:441
BasisType::FunctionValue FunctionValue
Definition basis.hpp:391
CurlValue curlcurl(size_t i, size_t iqn, const boost::multi_array< typename ScalarFamilyType::HessianValue, 2 > &ancestor_hessian_quad) const
Evaluate the curl curl of the i-th basis function at a quadrature point iqn, knowing all the values o...
Definition basis.hpp:735
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th basis function at point x.
Definition basis.hpp:687
static const bool hasDivergence
Definition basis.hpp:1396
DivergenceBasis< TangentFamily< ScalarFamilyType > > ScalarRotFamily(const TangentFamily< ScalarFamilyType > &tangent_family, const Face &F)
The following function creates the "scalar rot" basis of a TangentFamily on a face.
Definition basis.hpp:1440
TensorizedVectorFamily(const ScalarFamilyType &scalar_family)
Definition basis.hpp:633
static const bool hasFunction
Definition basis.hpp:332
static std::vector< VectorZd > homogeneous(const size_t L)
Definition basis.hpp:85
static boost::multi_array< typename detail::basis_evaluation_traits< TensorizedVectorFamily< BasisType, N >, BasisFunction >::ReturnValue, 2 > compute(const TensorizedVectorFamily< BasisType, N > &basis, const QuadratureRule &quad)
Evaluate a tensorized family at quadrature nodes (optimised compared the generic basis evaluation,...
Definition basis.hpp:2275
const Eigen::MatrixXd transposeOperator() const
Return the transpose operator, the rN^2 square matrix that to a given vector of coefficients on the M...
Definition basis.hpp:877
static const bool hasDivergence
Definition basis.hpp:170
double DivergenceValue
Definition basis.hpp:231
static constexpr const TensorRankE tensorRank
Definition basis.hpp:330
Eigen::Matrix< double, dimspace, dimspace > CurlValue
Definition basis.hpp:1464
ScalarFamilyType AncestorType
Definition basis.hpp:967
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1266
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1326
static const bool hasFunction
Definition basis.hpp:1622
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th function at point x.
Definition basis.hpp:428
double DivergenceValue
Definition basis.hpp:615
static const bool hasCurlCurl
Definition basis.hpp:172
static boost::multi_array< typename detail::basis_evaluation_traits< BasisType, BasisFunction >::ReturnValue, 2 > compute(const BasisType &basis, const QuadratureRule &quad)
Generic basis evaluation.
Definition basis.hpp:2224
VectorZd powers(size_t i) const
Returns the powers of the i-th basis function (its degree can be found using powers(i)....
Definition basis.hpp:202
static std::vector< Eigen::Vector2i > homogeneous(const size_t L)
Definition basis.hpp:119
static const bool hasCurl
Definition basis.hpp:334
static constexpr const bool hasAncestor
Definition basis.hpp:789
static const bool hasDivergence
Definition basis.hpp:1271
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1557
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.cpp:271
BasisType::DivergenceValue FunctionValue
Definition basis.hpp:1383
static const bool hasGradient
Definition basis.hpp:402
size_t dimension() const
Dimension of the basis.
Definition basis.hpp:254
static const bool hasGradient
Definition basis.hpp:1331
static constexpr const TensorRankE tensorRank
Definition basis.hpp:958
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:667
static const bool hasHessian
Definition basis.hpp:1163
static const bool hasCurlCurl
Definition basis.hpp:629
VectorRd GradientValue
Definition basis.hpp:229
static const bool hasCurl
Definition basis.hpp:1624
const JacobianType & jacobian() const
Return the Jacobian of the coordinate system transformation.
Definition basis.hpp:1648
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1091
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:831
DecomposePoly(const Face &F, const BasisType &basis)
Constructor for face.
Definition basis.hpp:2795
double scalar_product(const double &x, const double &y)
Scalar product between two reals.
Definition basis.cpp:308
VectorRd GradientValue
Definition basis.hpp:951
static constexpr const bool hasAncestor
Definition basis.hpp:1471
static const bool hasFunction
Definition basis.hpp:1472
static const bool hasCurl
Definition basis.hpp:169
CurlValue curlcurl(size_t i, const VectorRd &x) const
Evaluate the curl curl of the i-th basis function at point x.
Definition basis.hpp:1002
Face GeometricSupport
Definition basis.hpp:234
static constexpr const TensorRankE tensorRank
Definition basis.hpp:165
Eigen::Matrix< double, N, 1 > DivergenceValue
Definition basis.hpp:783
VectorRd CurlValue
Definition basis.hpp:159
std::function< T(const VectorRd &, const Cell *)> CellFType
type for function of point. T is the return type of the function
Definition basis.hpp:58
static const bool hasHessian
Definition basis.hpp:628
Eigen::MatrixXd gram_schmidt(boost::multi_array< T, 2 > &basis_eval, const std::function< double(size_t, size_t)> &inner_product)
Definition basis.hpp:2338
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1357
Eigen::Matrix< double, 2, dimspace > generators() const
Returns the generators of the basis functions.
Definition basis.hpp:1022
BasisType::HessianValue ReturnValue
Definition basis.hpp:1939
static const bool hasFunction
Definition basis.hpp:1543
Eigen::Matrix< double, dimspace, dimspace > CurlValue
Definition basis.hpp:1322
VectorRd CurlValue
Definition basis.hpp:324
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1321
BasisType::CurlValue CurlValue
Definition basis.hpp:1151
Eigen::Matrix< double, N, 1 > FunctionValue
Definition basis.hpp:612
static const bool hasDivergence
Definition basis.hpp:1333
std::function< T(const VectorRd &)> FType
type for function of point. T is the return type of the function
Definition basis.hpp:55
BasisType::HessianValue HessianValue
Definition basis.hpp:1046
double DivergenceValue
Definition basis.hpp:1323
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.cpp:74
const Eigen::MatrixXd & matrix() const
Return the coefficient matrix.
Definition basis.hpp:559
BasisType::HessianValue HessianValue
Definition basis.hpp:395
CurlValue curl(size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_gradient_quad) const
Evaluate the curl of the i-th basis function at a quadrature point iqn, knowing all the gradients of ...
Definition basis.hpp:697
VectorZd powers(size_t i) const
Returns the powers of the i-th basis function (not including the vector part)
Definition basis.hpp:1575
const JacobianType coordinates_system() const
Return the system of coordinates (basis in rows) on the face.
Definition basis.hpp:290
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2124
TensorizedVectorFamily< ScalarBasisType, N >::CurlValue ReturnValue
Definition basis.hpp:2119
void DivergenceValue
Definition basis.hpp:1386
const ScalarFamilyType & ancestor() const
Return the ancestor (family that has been tensorized)
Definition basis.hpp:865
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:754
Eigen::VectorXd l2_projection(const std::function< typename BasisType::FunctionValue(const VectorRd &)> &f, const BasisType &basis, QuadratureRule &quad, const boost::multi_array< typename BasisType::FunctionValue, 2 > &basis_quad, const Eigen::MatrixXd &mass_basis=Eigen::MatrixXd::Zero(1, 1))
Compute the L2-projection of a function.
Definition basis.hpp:2752
static constexpr const bool hasAncestor
Definition basis.hpp:1158
MatrixFamily(const ScalarFamilyType &scalar_family)
Definition basis.hpp:799
static const bool hasGradient
Definition basis.hpp:1160
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:1121
double DivergenceValue
Definition basis.hpp:1613
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1097
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1048
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.cpp:221
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th basis function at point x.
Definition basis.hpp:1221
Eigen::VectorXd integrate(const FType< T > &f, const BasisQuad< T > &B, const QuadratureRule &qr, size_t n_rows=0)
Computes the vector of integrals (f, phi_i)
Definition basis.hpp:2626
static const bool hasCurl
Definition basis.hpp:403
static const bool hasGradient
Definition basis.hpp:1394
size_t dimension() const
Dimension of the basis.
Definition basis.hpp:1710
boost::multi_array< typename BasisType::FunctionValue, 2 > m_on_basis_nodes
Definition basis.hpp:2891
double DivergenceValue
Definition basis.hpp:1465
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1940
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2158
VectorRd GradientValue
Definition basis.hpp:323
static const bool hasGradient
Definition basis.hpp:1473
VectorZd powers(size_t i) const
Returns the powers of the i-th basis function (not including the vector part)
Definition basis.hpp:1504
static const bool hasCurlCurl
Definition basis.hpp:406
MatrixFamily< ScalarBasisType, N >::DivergenceValue ReturnValue
Definition basis.hpp:2185
HessianValue hessian(size_t i, const VectorRd &x) const
Evaluate the Hessian of the i-th basis function at point x.
Definition basis.cpp:96
double FunctionValue
Definition basis.hpp:228
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1389
static const bool hasCurlCurl
Definition basis.hpp:965
static constexpr const bool hasAncestor
Definition basis.hpp:1542
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1611
ScalarFamilyType AncestorType
Definition basis.hpp:631
static constexpr const bool hasAncestor
Definition basis.hpp:1329
Eigen::Matrix< double, dimspace, dimspace > CurlValue
Definition basis.hpp:1612
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:990
GradientValue gradient(size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_gradient_quad) const
Evaluate the gradient of the i-th basis function at a quadrature point iqn, knowing all the gradients...
Definition basis.hpp:677
BasisType::HessianValue HessianValue
Definition basis.hpp:1153
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1426
static const bool hasHessian
Definition basis.hpp:1700
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:996
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1470
static const bool hasFunction
Definition basis.hpp:790
BasisFunctionE
Definition basis.hpp:1741
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:707
static const bool hasFunction
Definition basis.hpp:1393
BasisType AncestorType
Definition basis.hpp:1337
static const bool hasGradient
Definition basis.hpp:1053
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2190
MatrixFamily< ScalarBasisType, N >::FunctionValue ReturnValue
Definition basis.hpp:2153
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2026
VectorRd FunctionValue
Definition basis.hpp:1320
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th function at point x.
Definition basis.hpp:455
double FunctionValue
Definition basis.hpp:157
BasisType::GradientValue GradientValue
Definition basis.hpp:1150
VectorRd FunctionValue
Definition basis.hpp:1533
static const bool hasFunction
Definition basis.hpp:240
TensorizedVectorFamily< ScalarBasisType, N >::DivergenceValue ReturnValue
Definition basis.hpp:2086
BasisType m_basis
Definition basis.hpp:2889
boost::multi_array< VectorRd, 2 > vector_product(const boost::multi_array< VectorRd, 2 > &basis_quad, const VectorRd &v)
Compute the vector (cross) product between the evaluation of a basis and a constant vector.
Definition basis.cpp:324
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1893
static const bool hasDivergence
Definition basis.hpp:1475
static const bool hasCurlCurl
Definition basis.hpp:1273
DivergenceBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1403
TensorRankE
Definition basis.hpp:64
static std::vector< Eigen::Vector2i > complete(size_t degree)
Definition basis.hpp:131
static const bool hasCurlCurl
Definition basis.hpp:1627
VectorRd CurlValue
Definition basis.hpp:952
static const bool hasFunction
Definition basis.hpp:1696
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:847
static const bool hasCurl
Definition basis.hpp:1545
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th function at point x.
Definition basis.hpp:507
Eigen::Matrix< double, N, dimspace > GradientValue
Definition basis.hpp:613
CurlValue curlcurl(size_t i, const VectorRd &x) const
Evaluate the curl curl of the i-th basis function at point x. We use the formula curl curl = - Laplac...
Definition basis.hpp:723
static constexpr const bool hasAncestor
Definition basis.hpp:1695
void CurlValue
Definition basis.hpp:1385
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2036
CurlValue curl(size_t i, size_t iqn, const boost::multi_array< CurlValue, 2 > &ancestor_curl_quad) const
Evaluate the curl of the i-th function at a quadrature point iqn, knowing all the curls of ancestor b...
Definition basis.hpp:494
double DivergenceValue
Definition basis.hpp:325
static constexpr const bool hasAncestor
Definition basis.hpp:1267
static constexpr const bool hasAncestor
Definition basis.hpp:621
BasisType::DivergenceValue DivergenceValue
Definition basis.hpp:394
static constexpr const TensorRankE tensorRank
Definition basis.hpp:788
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1964
Eigen::Matrix2d HessianValue
Definition basis.hpp:954
BasisType::CurlValue ReturnValue
Definition basis.hpp:1892
Family< BasisType > m_on_basis
Definition basis.hpp:2890
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1391
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1486
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1155
Eigen::MatrixXd PermuteTensorization(const size_t a, const size_t b)
Returns the matrix giving the permutation of the tensorization of a family of size a with a family of...
Definition basis.cpp:338
QuadratureRule get_nodes() const
Return the set of nodes (useful to compute value of polynomial to decompose via evaluate_quad)
Definition basis.hpp:2867
static constexpr const bool hasAncestor
Definition basis.hpp:166
static const bool hasCurlCurl
Definition basis.hpp:1477
static const bool hasCurl
Definition basis.hpp:242
BasisType::CurlValue ReturnValue
Definition basis.hpp:1963
void CurlValue
Definition basis.hpp:1260
static const bool hasFunction
Definition basis.hpp:622
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1541
Eigen::Matrix< double, 2, dimspace > JacobianType
Definition basis.hpp:1618
Eigen::Vector2i powers(size_t i) const
Returns the powers of the i-th basis function.
Definition basis.hpp:1660
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:397
size_t dimension() const
Dimension of the basis.
Definition basis.hpp:346
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th basis function at point x.
Definition basis.cpp:229
double FunctionValue
Definition basis.hpp:322
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_value_quad)
Definition basis.hpp:1852
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_divergence_quad)
Definition basis.hpp:1922
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.cpp:151
Family(const BasisType &basis, const Eigen::MatrixXd &matrix)
Constructor.
Definition basis.hpp:411
static constexpr const TensorRankE tensorRank
Definition basis.hpp:238
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.cpp:278
boost::multi_array< outValue, 2 > transform_values_quad(const boost::multi_array< inValue, 2 > &B_quad, const FunctionType &F)
Takes an array B_quad of values at quadrature nodes and applies the function F to all of them....
Definition basis.hpp:1752
static const bool hasGradient
Definition basis.hpp:1623
void DivergenceValue
Definition basis.hpp:1261
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< double, 2 > &ancestor_basis_quad)
Definition basis.hpp:2004
const JacobianType & jacobian() const
Return the Jacobian of the coordinate system transformation.
Definition basis.hpp:284
static const bool hasFunction
Definition basis.hpp:1330
BasisType::GradientValue GradientValue
Definition basis.hpp:392
void HessianValue
Definition basis.hpp:1387
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2091
size_t dimension() const
Dimension of the family. This is actually the number of functions in the family, not necessarily line...
Definition basis.hpp:422
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:196
void HessianValue
Definition basis.hpp:616
TensorizedVectorFamily< ScalarBasisType, N >::CurlValue ReturnValue
Definition basis.hpp:2053
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:181
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_gradient_quad)
Definition basis.hpp:1876
const size_t matrixSize() const
Return the dimension of the matrices in the family.
Definition basis.hpp:871
static const bool hasGradient
Definition basis.hpp:168
HessianValue hessian(size_t i, const VectorRd &x) const
Evaluate the hessian of the i-th basis function at point x.
Definition basis.cpp:40
static const bool hasCurl
Definition basis.hpp:793
Family< MatrixFamily< ScalarBasisType, N > > IsotropicMatrixFamily(const ScalarBasisType &B)
From a scalar family B, constructs a "Family" of "MatrixFamily" (built on B, of size NxN) that repres...
Definition basis.hpp:1810
static const bool hasDivergence
Definition basis.hpp:1546
BasisType AncestorType
Definition basis.hpp:1275
ScalarFamilyType::GeometricSupport GeometricSupport
Definition basis.hpp:786
DivergenceValue divergence(size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_gradient_quad) const
Evaluate the divergence of the i-th basis function at a quadrature point iqn, knowing all the gradien...
Definition basis.hpp:715
static const bool hasFunction
Definition basis.hpp:401
Family< BasisType > family(boost::multi_array< typename BasisType::FunctionValue, 2 > &values)
Returns the decomposed polynomials as a Family of the provided basis.
Definition basis.hpp:2873
size_t m_dim
Definition basis.hpp:2888
static const bool hasCurl
Definition basis.hpp:1270
BasisType::DivergenceValue DivergenceValue
Definition basis.hpp:1152
static constexpr const TensorRankE tensorRank
Definition basis.hpp:399
static const bool hasGradient
Definition basis.hpp:791
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1363
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1463
const VectorRd & normal() const
Return the normal to the face used in the computation of the curl.
Definition basis.hpp:1719
const Eigen::MatrixXd traceOperator() const
Returns the matrix corresponding to the trace, expressed as a linear operator between that MatrixFami...
Definition basis.hpp:914
static const bool hasCurl
Definition basis.hpp:1332
Eigen::Matrix3d MatrixRd
Definition basis.hpp:51
BasisType::GradientValue FunctionValue
Definition basis.hpp:1258
static const bool hasHessian
Definition basis.hpp:171
static const bool hasDivergence
Definition basis.hpp:1055
static const bool hasDivergence
Definition basis.hpp:243
static const bool hasGradient
Definition basis.hpp:1697
Eigen::Vector2i powers(size_t i) const
Returns the powers of the i-th basis function (its degree can be found using powers(i)....
Definition basis.hpp:296
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1289
static const bool hasHessian
Definition basis.hpp:964
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:272
static const bool hasCurl
Definition basis.hpp:1161
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:647
static const bool hasFunction
Definition basis.hpp:1052
size_t m_nb_nodes
Definition basis.hpp:2892
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1328
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1569
void GradientValue
Definition basis.hpp:781
static const bool hasCurl
Definition basis.hpp:627
Edge GeometricSupport
Definition basis.hpp:328
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_curlcurl_quad)
Definition basis.hpp:1970
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1534
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:571
VectorRd GradientValue
Definition basis.hpp:158
static const bool hasFunction
Definition basis.hpp:1159
static const bool hasFunction
Definition basis.hpp:1268
ScalarFamilyType::GeometricSupport GeometricSupport
Definition basis.hpp:618
size_t dimension() const
Dimension of the basis.
Definition basis.hpp:1636
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1620
static const bool hasDivergence
Definition basis.hpp:404
RestrictedBasis(const BasisType &basis, const size_t &dimension)
Constructor.
Definition basis.hpp:1169
HessianValue hessian(size_t i, const VectorRd &x) const
Evaluate the hessian of the i-th function at point x.
Definition basis.hpp:533
static const bool hasFunction
Definition basis.hpp:167
size_t shift() const
Return the shift.
Definition basis.hpp:1079
static std::vector< VectorZd > complete(const size_t degree)
Definition basis.hpp:100
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1205
Eigen::Matrix< double, dimspace, dimspace *dimspace > HessianValue
Definition basis.hpp:1324
BasisType::DivergenceValue DivergenceValue
Definition basis.hpp:1045
static const bool hasDivergence
Definition basis.hpp:963
static const bool hasCurl
Definition basis.hpp:1698
Eigen::MatrixXd compute_weighted_gram_matrix(const FType< VectorRd > &f, const BasisQuad< VectorRd > &B1, const BasisQuad< double > &B2, const QuadratureRule &qr, size_t n_rows, size_t n_cols)
Computes the Gram-like matrix of integrals (f dot phi_i, phi_j)
Definition basis.cpp:491
HessianValue hessian(size_t i, const VectorRd &x) const
Evaluate the hessian of the i-th basis function at point x.
Definition basis.hpp:1129
BasisType::GradientValue GradientValue
Definition basis.hpp:1043
static const bool hasGradient
Definition basis.hpp:241
TensorizedVectorFamily< ScalarBasisType, N >::FunctionValue ReturnValue
Definition basis.hpp:1989
static const bool hasCurl
Definition basis.hpp:962
DecomposePoly(const Cell &T, const BasisType &basis)
Constructor for cell.
Definition basis.hpp:2831
static const bool hasDivergence
Definition basis.hpp:1625
static constexpr const bool hasAncestor
Definition basis.hpp:959
static const bool hasFunction
Definition basis.hpp:960
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_curl_quad)
Definition basis.hpp:1899
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1685
static const bool hasHessian
Definition basis.hpp:405
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1016
FunctionValue function(size_t i, size_t iqn, const boost::multi_array< double, 2 > &ancestor_value_quad) const
Evaluate the i-th basis function at a quadrature point iqn, knowing all the values of ancestor basis ...
Definition basis.hpp:839
static const bool hasHessian
Definition basis.hpp:336
BasisType AncestorType
Definition basis.hpp:408
static const bool hasHessian
Definition basis.hpp:794
constexpr const BasisType & ancestor() const
Return the underlying complete basis.
Definition basis.hpp:1188
static const bool hasGradient
Definition basis.hpp:623
static const bool hasCurlCurl
Definition basis.hpp:337
Face GeometricSupport
Definition basis.hpp:1616
static const bool hasCurlCurl
Definition basis.hpp:245
BasisType::FunctionValue ReturnValue
Definition basis.hpp:1845
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2068
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1264
static constexpr const bool hasAncestor
Definition basis.hpp:331
static const bool hasGradient
Definition basis.hpp:1544
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th basis function at point x.
Definition basis.hpp:1113
static constexpr const bool hasAncestor
Definition basis.hpp:1392
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1994
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1414
TangentFamily(const ScalarFamilyType &scalar_family, const Eigen::Matrix< double, 2, dimspace > &generators)
Constructor.
Definition basis.hpp:970
const VectorRd & normal() const
Return the normal to the face used in the computation of the curl.
Definition basis.hpp:278
BasisType::FunctionValue FunctionValue
Definition basis.hpp:1042
static const bool hasDivergence
Definition basis.hpp:792
BasisType AncestorType
Definition basis.hpp:1400
static const bool hasCurlCurl
Definition basis.hpp:1057
void HessianValue
Definition basis.hpp:784
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.cpp:127
BasisType AncestorType
Definition basis.hpp:1059
GradientValue gradient(size_t i, size_t iqn, const boost::multi_array< GradientValue, 2 > &ancestor_gradient_quad) const
Evaluate the gradient of the i-th function at a quadrature point iqn, knowing all the gradients of an...
Definition basis.hpp:468
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:1213
HessianValue hessian(size_t i, const VectorRd &x) const
Evaluate the hessian of the i-th basis function at point x.
Definition basis.hpp:1237
static const bool hasCurlCurl
Definition basis.hpp:1701
static boost::multi_array< typename detail::basis_evaluation_traits< Family< BasisType >, BasisFunction >::ReturnValue, 2 > compute(const Family< BasisType > &basis, const QuadratureRule &quad)
Evaluate a 'Family' of functions at quadrature nodes (optimised compared the generic basis evaluation...
Definition basis.hpp:2248
HessianValue hessian(size_t i, size_t iqn, const boost::multi_array< HessianValue, 2 > &ancestor_hessian_quad) const
Evaluate the hessian of the i-th function at a quadrature point iqn, knowing all the hessian of ances...
Definition basis.hpp:546
Eigen::Matrix< double, dimspace, dimspace > CurlValue
Definition basis.hpp:1686
ScalarFamilyType AncestorType
Definition basis.hpp:797
static const bool hasCurl
Definition basis.hpp:1054
Eigen::Matrix< double, N, N > FunctionValue
Definition basis.hpp:780
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.cpp:158
ShiftedBasis(const BasisType &basis, const int shift)
Constructor.
Definition basis.hpp:1062
static const bool hasCurl
Definition basis.hpp:1395
static const bool hasCurlCurl
Definition basis.hpp:795
BasisType::CurlValue CurlValue
Definition basis.hpp:393
BasisType::GradientValue ReturnValue
Definition basis.hpp:1869
BasisType::CurlValue CurlValue
Definition basis.hpp:1044
static std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> skew_symmetrise_matrix
Function to skew-symmetrise a matrix (useful together with transform_values_quad)
Definition basis.hpp:1805
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.cpp:22
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.cpp:296
CurlBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1340
static const bool hasHessian
Definition basis.hpp:1397
constexpr const ScalarFamilyType & ancestor() const
Return the ancestor (family that has been tensorized)
Definition basis.hpp:748
static const std::vector< VectorRd > basisRd
Definition basis.hpp:52
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1654
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2101
Cell GeometricSupport
Definition basis.hpp:163
size_t dimension() const
Return the dimension of the basis.
Definition basis.hpp:1182
GradientBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1278
void HessianValue
Definition basis.hpp:1262
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1351
static constexpr const bool hasAncestor
Definition basis.hpp:1621
static const bool hasHessian
Definition basis.hpp:244
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th function at point x.
Definition basis.hpp:481
static const bool hasDivergence
Definition basis.hpp:1162
Eigen::MatrixXd compute_gram_matrix(const boost::multi_array< VectorRd, 2 > &B1, const boost::multi_array< double, 2 > &B2, const QuadratureRule &qr)
Compute the Gram-like matrix given a family of vector-valued and one of scalar-valued functions by te...
Definition basis.cpp:353
Face GeometricSupport
Definition basis.hpp:1690
Cell GeometricSupport
Definition basis.hpp:1468
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1420
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1846
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1694
static const bool hasDivergence
Definition basis.hpp:626
QuadratureRule m_nodes
Nodes for the interpolation.
Definition basis.hpp:2893
static constexpr const TensorRankE tensorRank
Definition basis.hpp:620
MatrixRd HessianValue
Definition basis.hpp:232
static const bool hasCurlCurl
Definition basis.hpp:1398
double DivergenceValue
Definition basis.hpp:953
double HessianValue
Definition basis.hpp:326
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.cpp:122
static const bool hasHessian
Definition basis.hpp:1272
size_t dimension() const
Return the dimension of the family.
Definition basis.hpp:825
const Eigen::MatrixXd symmetricBasis() const
Returns the matrix that can be used, in a Family of this MatrixBasis, to create a basis of the subspa...
Definition basis.hpp:889
VectorRd FunctionValue
Definition basis.hpp:1684
Cell GeometricSupport
Definition basis.hpp:1539
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1916
constexpr const ScalarFamilyType & ancestor() const
Return the ancestor (family used for the tangent)
Definition basis.hpp:1010
static const bool hasCurlCurl
Definition basis.hpp:1548
Eigen::Matrix< double, 2, dimspace > JacobianType
Definition basis.hpp:236
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.cpp:29
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1194
@ CurlCurl
Definition basis.hpp:1747
@ Function
Definition basis.hpp:1742
@ Hessian
Definition basis.hpp:1746
@ Gradient
Definition basis.hpp:1743
@ Divergence
Definition basis.hpp:1745
@ Curl
Definition basis.hpp:1744
@ Matrix
Definition basis.hpp:67
@ Vector
Definition basis.hpp:66
@ Scalar
Definition basis.hpp:65
size_t L
Definition HHO_DiffAdvecReac.hpp:46
std::vector< QuadratureNode > QuadratureRule
Definition quadraturerule.hpp:61
Definition PastixInterface.hpp:16
Definition ddr-magnetostatics.hpp:41
Structure to decompose a set of polynomials on a basis on a face or a cell.
Definition basis.hpp:2788
Compute vectors listing the powers of monomial basis functions (for a cell or face,...
Definition basis.hpp:77
Basis dimensions for various polynomial spaces on edges/faces/elements (when relevant): Pk,...
Definition polynomialspacedimension.hpp:54
Basis evaluation traits. Only specialization of 'BasisFunction' (=Function, Gradient,...
Definition basis.hpp:1837
Evaluate a basis at quadrature nodes. 'BasisFunction' (=Function, Gradient, Curl, Divergence or Hessi...
Definition basis.hpp:2220