#include <GaussLobatto.hpp>

Public Member Functions | |
| const Vector< real_t > & | vertices () const |
| size_t | numberOfPoints () const |
| real_t | operator() (const size_t &i) const |
| real_t | weight (const size_t &i) const |
| ~GaussLobatto () | |
Private Member Functions | |
| void | __bisection (const Interval &interval, real_t &root) |
| void | __newton (const Interval &interval, real_t &root) |
| GaussLobatto (const size_t &d) | |
| GaussLobatto (const GaussLobatto &) | |
Private Attributes | |
| const size_t | __degree |
| Vector< real_t > | __vertices |
| Vector< real_t > | __weights |
Friends | |
| class | GaussLobattoManager |
Definition at line 36 of file GaussLobatto.hpp.
| GaussLobatto::GaussLobatto | ( | const size_t & | d | ) | [private] |
Constructor
| d | the degree |
Definition at line 68 of file GaussLobatto.cpp.
References __degree, __newton(), __vertices, __weights, and LegendreBasis::getValue().
00069 : __degree(degree), 00070 __vertices(degree+1), 00071 __weights(degree+1) 00072 { 00073 if (degree==0) { 00074 __vertices[0] = 0.; 00075 __weights[0] = 2; 00076 } else { 00077 __vertices[0] = -1; 00078 __vertices[degree] = 1; 00079 00080 if (__degree!=1) { 00081 if (__degree==2) { 00082 Interval intervalRef(-1.,1.); 00083 00084 // Use Bisection method 00085 //__bisection(intervalRef,__vertices[1]); 00086 00087 // Use Newton method 00088 __newton(intervalRef,__vertices[1]); 00089 00090 } else { 00091 Interval intervalRef(-1.,1.); 00092 GaussLobatto gaussN_1(__degree-1); 00093 00094 for (size_t i=1; i<__degree; ++i) { 00095 Interval interval( gaussN_1(i-1),gaussN_1(i)); 00096 00097 // Use Bisection method 00098 //__bisection(interval,__vertices[i]); 00099 00100 // Use Newton method 00101 __newton(interval,__vertices[i]); 00102 } 00103 } 00104 } 00105 00106 const real_t w0= 2./(__degree*(__degree+1)); 00107 __weights[0]= w0; 00108 __weights[__degree]=w0; 00109 00110 LegendreBasis Ln(__degree); 00111 for (size_t i=1;i<__degree; ++i) { 00112 __weights[i]= w0 / (Ln.getValue(__vertices[i],__degree) * Ln.getValue(__vertices[i],__degree)); 00113 } 00114 } 00115 }

| GaussLobatto::GaussLobatto | ( | const GaussLobatto & | ) | [explicit, private] |
| GaussLobatto::~GaussLobatto | ( | ) | [inline] |
| void GaussLobatto::__bisection | ( | const Interval & | interval, | |
| real_t & | root | |||
| ) | [private] |
Bisection to compute recursively qudrature vertices
| interval | interval of research | |
| root | root of the polynom |
Definition at line 26 of file GaussLobatto.cpp.
References __degree, Interval::a(), Interval::b(), and LegendreBasis::getDerivativeValue().
00028 { 00029 LegendreBasis Ln(__degree); 00030 const real_t epsilon = 1e-10; 00031 00032 real_t a = interval.a(); 00033 real_t b = interval.b(); 00034 00035 do { 00036 root= 0.5*(a+b); 00037 00038 if (Ln.getDerivativeValue(root,__degree)== 0) return; 00039 if (Ln.getDerivativeValue(root,__degree)*Ln.getDerivativeValue(a,__degree) <0) { 00040 b=root; 00041 } else { 00042 a=root; 00043 } 00044 } while (std:: abs(b-a) > epsilon ); 00045 }

| void GaussLobatto::__newton | ( | const Interval & | interval, | |
| real_t & | root | |||
| ) | [private] |
Definition at line 48 of file GaussLobatto.cpp.
References __degree, Interval::a(), Interval::b(), LegendreBasis::getDerivativeValue(), and LegendreBasis::getSecondDerivativeValue().
Referenced by GaussLobatto().
00050 { 00051 LegendreBasis Ln(__degree); 00052 00053 const real_t epsilon = 1e-15; 00054 real_t a = interval.a(); 00055 real_t b = interval.b(); 00056 root= 0.5*(a+b); 00057 real_t xn= root; 00058 real_t delta = 0; 00059 do { 00060 xn = root; 00061 if (Ln.getDerivativeValue(xn,__degree)== 0) return; 00062 delta = (Ln.getDerivativeValue(xn,__degree) / Ln.getSecondDerivativeValue(xn,__degree)); 00063 root -= delta; 00064 } while (std:: abs(delta) > epsilon ); 00065 }

| const Vector<real_t>& GaussLobatto::vertices | ( | ) | const [inline] |
Definition at line 56 of file GaussLobatto.hpp.
References __vertices.
Referenced by SpectralFEMPreconditioner::Internal::computes(), and SpectralFEMPreconditioner::Internal::computesTransposed().
00057 { 00058 return __vertices; 00059 }
| size_t GaussLobatto::numberOfPoints | ( | ) | const [inline] |
Read-only access to the number of quadrature nodes
Definition at line 66 of file GaussLobatto.hpp.
References __vertices, and Vector< T >::size().
Referenced by SpectralFEMPreconditioner::Internal::__lagrangeToLegendre(), SpectralFEMPreconditioner::Internal::__lagrangeToLegendreTransposed(), RealExpressionIntegrate::execute(), and SpectralFunction::operator=().
00067 { 00068 return __vertices.size(); 00069 }

| real_t GaussLobatto::operator() | ( | const size_t & | i | ) | const [inline] |
Read-only access to the i-th quadrature point
| i | the number of the quadrature point |
Definition at line 78 of file GaussLobatto.hpp.
References __vertices.
00079 { 00080 return __vertices[i]; 00081 }
| real_t GaussLobatto::weight | ( | const size_t & | i | ) | const [inline] |
Read-only access to the weight of the i-th quadrature point
| i | the number of the quadrature point |
Definition at line 90 of file GaussLobatto.hpp.
References __weights.
Referenced by RealExpressionIntegrate::execute(), and SpectralFunction::operator=().
00091 { 00092 return __weights[i]; 00093 }
friend class GaussLobattoManager [friend] |
Definition at line 96 of file GaussLobatto.hpp.
const size_t GaussLobatto::__degree [private] |
degree
Definition at line 39 of file GaussLobatto.hpp.
Referenced by __bisection(), __newton(), and GaussLobatto().
Vector<real_t> GaussLobatto::__vertices [private] |
nodes for the quadrature
Definition at line 41 of file GaussLobatto.hpp.
Referenced by GaussLobatto(), numberOfPoints(), operator()(), and vertices().
Vector<real_t> GaussLobatto::__weights [private] |
weights for the quadrature
Definition at line 42 of file GaussLobatto.hpp.
Referenced by GaussLobatto(), and weight().
1.5.6