https://mooseframework.inl.gov
Public Member Functions | Private Attributes | List of all members
PolynomialQuadrature::Legendre Class Reference

Uniform distributions use Legendre polynomials. More...

#include <PolynomialQuadrature.h>

Inheritance diagram for PolynomialQuadrature::Legendre:
[legend]

Public Member Functions

 Legendre (const Real lower_bound, const Real upper_bound)
 
virtual void store (std::ostream &stream, void *context) const override
 
virtual void store (nlohmann::json &json) const override
 
virtual Real compute (const unsigned int order, const Real x, const bool normalize=true) const override
 Legendre polynomial using static function then scales by <P_n^2> = 1 / (2n+1) More...
 
virtual Real computeDerivative (const unsigned int order, const Real x, const unsigned int m=1) const override
 Compute derivative of Legendre polynomial Recursive algorithm: d^m/dx^m P_n = (n + m - 1)d^(m-1)/dx^(m-1) P_{n-1} + xd^m/dx^m P_{n-1}. More...
 
Real computeDerivativeRef (const unsigned int order, const Real x, const unsigned int m=1) const
 
virtual Real innerProduct (const unsigned int order) const override
 
virtual void gaussQuadrature (const unsigned int order, std::vector< Real > &points, std::vector< Real > &weights) const override
 Gauss-Legendre quadrature: sum(weights) = 2. More...
 
virtual void clenshawQuadrature (const unsigned int order, std::vector< Real > &points, std::vector< Real > &weights) const override
 Just normal Clenshaw Curtis with shifted points. More...
 
Real productIntegral (const std::vector< unsigned int > order) const
 

Private Attributes

const Real _lower_bound
 
const Real _upper_bound
 

Detailed Description

Uniform distributions use Legendre polynomials.

Definition at line 58 of file PolynomialQuadrature.h.

Constructor & Destructor Documentation

◆ Legendre()

Legendre::Legendre ( const Real  lower_bound,
const Real  upper_bound 
)

Definition at line 117 of file PolynomialQuadrature.C.

Member Function Documentation

◆ clenshawQuadrature()

void Legendre::clenshawQuadrature ( const unsigned int  order,
std::vector< Real > &  points,
std::vector< Real > &  weights 
) const
overridevirtual

Just normal Clenshaw Curtis with shifted points.

Reimplemented from PolynomialQuadrature::Polynomial.

Definition at line 188 of file PolynomialQuadrature.C.

191 {
192  clenshaw_curtis(order, points, weights);
193  Real dx = (_upper_bound - _lower_bound) / 2.0;
194  Real xav = (_upper_bound + _lower_bound) / 2.0;
195  for (unsigned int i = 0; i < points.size(); ++i)
196  points[i] = points[i] * dx + xav;
197 }
void clenshaw_curtis(const unsigned int order, std::vector< Real > &points, std::vector< Real > &weights)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ compute()

Real Legendre::compute ( const unsigned int  order,
const Real  x,
const bool  normalize = true 
) const
overridevirtual

Legendre polynomial using static function then scales by <P_n^2> = 1 / (2n+1)

Reimplemented from PolynomialQuadrature::Polynomial.

Definition at line 140 of file PolynomialQuadrature.C.

141 {
142  if ((x < _lower_bound) || (x > _upper_bound))
143  ::mooseError("The requested polynomial point is outside of distribution bounds.");
144 
145  Real val = legendre(order, x, _lower_bound, _upper_bound);
146  if (normalize)
147  val /= innerProduct(order);
148 
149  return val;
150 }
const std::vector< double > x
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real legendre(const unsigned int order, const Real x, const Real lower_bound=-1.0, const Real upper_bound=1.0)
Legendre polynomial of specified order.
virtual Real innerProduct(const unsigned int order) const override

◆ computeDerivative()

Real Legendre::computeDerivative ( const unsigned int  order,
const Real  x,
const unsigned int  m = 1 
) const
overridevirtual

Compute derivative of Legendre polynomial Recursive algorithm: d^m/dx^m P_n = (n + m - 1)d^(m-1)/dx^(m-1) P_{n-1} + xd^m/dx^m P_{n-1}.

Reimplemented from PolynomialQuadrature::Polynomial.

Definition at line 153 of file PolynomialQuadrature.C.

154 {
155  if ((x < _lower_bound) || (x > _upper_bound))
156  ::mooseError("The requested polynomial point is outside of distribution bounds.");
157 
158  Real xref = 2.0 / (_upper_bound - _lower_bound) * (x - (_upper_bound + _lower_bound) / 2.0);
159  Real Jac = pow(2.0 / (_upper_bound - _lower_bound), m);
160  return Jac * computeDerivativeRef(order, xref, m);
161 }
Real computeDerivativeRef(const unsigned int order, const Real x, const unsigned int m=1) const
const std::vector< double > x
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ computeDerivativeRef()

Real Legendre::computeDerivativeRef ( const unsigned int  order,
const Real  x,
const unsigned int  m = 1 
) const

Definition at line 164 of file PolynomialQuadrature.C.

Referenced by computeDerivative().

167 {
168  if (m == 0)
169  return legendre(order, xref);
170  else if (m > order)
171  return 0.0;
172  else if (order == 1)
173  return 1.0;
174  else
175  return static_cast<Real>(order + m - 1) * computeDerivativeRef(order - 1, xref, m - 1) +
176  xref * computeDerivativeRef(order - 1, xref, m);
177 }
Real computeDerivativeRef(const unsigned int order, const Real x, const unsigned int m=1) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real legendre(const unsigned int order, const Real x, const Real lower_bound=-1.0, const Real upper_bound=1.0)
Legendre polynomial of specified order.

◆ gaussQuadrature()

void Legendre::gaussQuadrature ( const unsigned int  order,
std::vector< Real > &  points,
std::vector< Real > &  weights 
) const
overridevirtual

Gauss-Legendre quadrature: sum(weights) = 2.

Implements PolynomialQuadrature::Polynomial.

Definition at line 180 of file PolynomialQuadrature.C.

183 {
184  gauss_legendre(order, points, weights, _lower_bound, _upper_bound);
185 }
void gauss_legendre(const unsigned int order, std::vector< Real > &points, std::vector< Real > &weights, const Real lower_bound, const Real upper_bound)
Generalized formula for any polynomial order.

◆ innerProduct()

Real Legendre::innerProduct ( const unsigned int  order) const
overridevirtual

Implements PolynomialQuadrature::Polynomial.

Definition at line 200 of file PolynomialQuadrature.C.

Referenced by compute().

201 {
202  return 1. / (2. * static_cast<Real>(order) + 1.);
203 }
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ productIntegral()

Real PolynomialQuadrature::Polynomial::productIntegral ( const std::vector< unsigned int order) const
inherited

Definition at line 89 of file PolynomialQuadrature.C.

90 {
91  const unsigned int nprod = order.size();
92 
93  if (nprod == 1)
94  return (order[0] == 0 ? 1.0 : 0.0);
95  else if (nprod == 2)
96  return (order[0] == order[1] ? innerProduct(order[0]) : 0.0);
97 
98  unsigned int poly_order = std::accumulate(order.begin(), order.end(), 0);
99  unsigned int quad_order = (poly_order % 2 == 0 ? poly_order : poly_order - 1) / 2;
100 
101  std::vector<Real> xq;
102  std::vector<Real> wq;
103  gaussQuadrature(quad_order, xq, wq);
104 
105  Real val = 0.0;
106  for (unsigned int q = 0; q < xq.size(); ++q)
107  {
108  Real prod = wq[q];
109  for (unsigned int i = 0; i < nprod; ++i)
110  prod *= compute(order[i], xq[q], false);
111  val += prod;
112  }
113 
114  return val / std::accumulate(wq.begin(), wq.end(), 0.0);
115 }
virtual void gaussQuadrature(const unsigned int order, std::vector< Real > &points, std::vector< Real > &weights) const =0
virtual Real compute(const unsigned int order, const Real x, const bool normalize=true) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real innerProduct(const unsigned int order) const =0

◆ store() [1/2]

void Legendre::store ( std::ostream &  stream,
void context 
) const
overridevirtual

Reimplemented from PolynomialQuadrature::Polynomial.

Definition at line 123 of file PolynomialQuadrature.C.

124 {
125  std::string type = "Legendre";
126  dataStore(stream, type, context);
127  dataStore(stream, _lower_bound, context);
128  dataStore(stream, _upper_bound, context);
129 }
void dataStore(std::ostream &stream, std::unique_ptr< const PolynomialQuadrature::Polynomial > &ptr, void *context)

◆ store() [2/2]

void Legendre::store ( nlohmann::json &  json) const
overridevirtual

Reimplemented from PolynomialQuadrature::Polynomial.

Definition at line 132 of file PolynomialQuadrature.C.

133 {
134  json["type"] = "Legendre";
135  json["lower_bound"] = _lower_bound;
136  json["upper_bound"] = _upper_bound;
137 }

Member Data Documentation

◆ _lower_bound

const Real PolynomialQuadrature::Legendre::_lower_bound
private

◆ _upper_bound

const Real PolynomialQuadrature::Legendre::_upper_bound
private

The documentation for this class was generated from the following files: