PFASST++
pfasst::quadrature Namespace Reference

Functionality related to computing quadrature nodes and weights. More...

Classes

class  ClenshawCurtis
 Quadrature handler for Clenshaw-Curtis quadrature. More...
 
class  GaussLegendre
 Quadrature handler for Gauss-Legendre quadrature. More...
 
class  GaussLobatto
 Quadrature handler for Gauss-Lobatto quadrature. More...
 
class  GaussRadau
 Quadrature handler for Gauss-Radau quadrature. More...
 
class  IQuadrature
 Interface for quadrature handlers. More...
 
class  Polynomial
 Representation of a polynomial including differentiation, integration and root finding. More...
 
class  Uniform
 Quadrature handler for uniform distributed nodes. More...
 

Enumerations

enum  QuadratureType : int {
  QuadratureType::GaussLegendre = 0, QuadratureType::GaussLobatto = 1, QuadratureType::GaussRadau = 2, QuadratureType::ClenshawCurtis = 3,
  QuadratureType::Uniform = 4, QuadratureType::UNDEFINED = -1
}
 Quadrature type descriptors. More...
 

Functions

template<typename scalar >
static Polynomial< scalar > build_polynomial (const size_t node, const vector< scalar > &nodes)
 
template<typename scalar >
static Matrix< scalar > compute_q_matrix (const vector< scalar > &from, const vector< scalar > &to)
 Compute quadrature matrix \( Q \) between two sets of nodes. More...
 
template<typename scalar >
static Matrix< scalar > compute_q_matrix (const vector< scalar > &nodes)
 Compute quadrature matrix \( Q \) for one set of nodes. More...
 
template<typename scalar >
static Matrix< scalar > compute_q_matrix (const Matrix< scalar > &s_mat)
 Compute quadrature matrix \( Q \) from a given node-to-node quadrature matrix \( S \). More...
 
template<typename scalar >
static Matrix< scalar > compute_s_matrix (const Matrix< scalar > &q_mat)
 Compute node-to-node quadrature matrix \( S \) from a given quadrature matrix \( Q \). More...
 
template<typename scalar >
static Matrix< scalar > compute_s_matrix (const vector< scalar > &from, const vector< scalar > &to)
 Compute node-to-node quadrature matrix \( S \) from two given sets of nodes. More...
 
template<typename scalar >
static vector< scalar > compute_q_vec (const vector< scalar > &nodes)
 Compute vector \( q \) for integration from \( 0 \) to \( 1 \) for given set of nodes. More...
 
template<typename precision = pfasst::time_precision>
shared_ptr< IQuadrature< precision > > quadrature_factory (const size_t nnodes, const QuadratureType qtype)
 Instantiates quadrature handler for given number of nodes and type descriptor. More...
 
template<typename precision = pfasst::time_precision>
vector< precision > compute_nodes (size_t nnodes, QuadratureType qtype)
 Compute quadrature nodes for given quadrature type descriptor. More...
 
template<typename precision = time_precision>
Matrix< precision > compute_interp (vector< precision > dst, vector< precision > src)
 

Detailed Description

Functionality related to computing quadrature nodes and weights.

Note
Please note, that all quadrature nodes are in the range \( [0, 1] \).

Enumeration Type Documentation

Quadrature type descriptors.

Since
v0.3.0
Enumerator
GaussLegendre 

Gauss-Legendre quadrature

GaussLobatto 

Gauss-Lobatto quadrature

GaussRadau 

Gauss-Radau quadrature

ClenshawCurtis 

Clenshaw-Curtis quadrature

Uniform 

Uniform quadrature.

UNDEFINED 

Definition at line 32 of file interface.hpp.

Function Documentation

template<typename scalar >
static Polynomial<scalar> pfasst::quadrature::build_polynomial ( const size_t  node,
const vector< scalar > &  nodes 
)
static

Definition at line 43 of file interface.hpp.

Referenced by compute_q_matrix(), and compute_q_vec().

+ Here is the caller graph for this function:

template<typename precision = time_precision>
Matrix<precision> pfasst::quadrature::compute_interp ( vector< precision >  dst,
vector< precision >  src 
)
Todo:
write documentation
Template Parameters
precisionnumerical type of the interpolation (e.g. double)

Definition at line 95 of file quadrature.hpp.

template<typename precision = pfasst::time_precision>
vector<precision> pfasst::quadrature::compute_nodes ( size_t  nnodes,
QuadratureType  qtype 
)

Compute quadrature nodes for given quadrature type descriptor.

Template Parameters
precisionnumerical type of the nodes (e.g. double)
Parameters
[in]nnodesnumber of quadrature nodes to compute
[in]qtypetype descriptor of the quadrature nodes
Returns
std::vector of quadrature nodes of given type
See also
pfasst::quadrature::QuadratureType for valid types
pfasst::quadrature::quadrature_factory for further details

Definition at line 84 of file quadrature.hpp.

template<typename scalar >
static Matrix<scalar> pfasst::quadrature::compute_q_matrix ( const vector< scalar > &  from,
const vector< scalar > &  to 
)
static

Compute quadrature matrix \( Q \) between two sets of nodes.

Computing the quadrature matrix \( Q \) for polynomial-based integration from one set of quadrature nodes (from) to another set of quadrature nodes (to).

Template Parameters
scalarprecision of quadrature (i.e. double)
Parameters
[in]fromfirst set of quadrature nodes
[in]tosecond set of quadrature nodes
Returns
quadrature matrix \( Q \) with to.size() rows and from.size() colums
Precondition
For correctness of the algorithm it is assumed, that both sets of nodes are in the range \( [0, 1] \).
Since
v0.3.0

Definition at line 80 of file interface.hpp.

References build_polynomial(), pfasst::quadrature::Polynomial< CoeffT >::evaluate(), and pfasst::quadrature::Polynomial< CoeffT >::integrate().

Referenced by pfasst::quadrature::IQuadrature< precision >::compute_weights().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

template<typename scalar >
static Matrix<scalar> pfasst::quadrature::compute_q_matrix ( const vector< scalar > &  nodes)
static

Compute quadrature matrix \( Q \) for one set of nodes.

Template Parameters
scalarprecision of quadrature (i.e. double)
Parameters
[in]nodesquadrature nodes to compute \( Q \) matrix for
Since
v0.3.0

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 113 of file interface.hpp.

template<typename scalar >
static Matrix<scalar> pfasst::quadrature::compute_q_matrix ( const Matrix< scalar > &  s_mat)
static

Compute quadrature matrix \( Q \) from a given node-to-node quadrature matrix \( S \).

Template Parameters
scalarprecision of quadrature (i.e. double)
Parameters
[in]s_mat\( S \) matrix to compute \( Q \) from
See also
pfasst::quadrature::compute_s_matrix
Since
v0.3.0

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 131 of file interface.hpp.

Referenced by compute_s_matrix().

+ Here is the caller graph for this function:

template<typename scalar >
static vector<scalar> pfasst::quadrature::compute_q_vec ( const vector< scalar > &  nodes)
static

Compute vector \( q \) for integration from \( 0 \) to \( 1 \) for given set of nodes.

This equals to the last row of the quadrature matrix \( Q \) for the given set of nodes.

Template Parameters
scalarprecision of quadrature (i.e. double)
Parameters
[in]nodesquadrature nodes to compute \( Q \) matrix for
Precondition
For correctness of the algorithm it is assumed, that the nodes are in the range \( [0, 1] \).
Since
v0.3.0

Definition at line 200 of file interface.hpp.

References build_polynomial(), pfasst::quadrature::Polynomial< CoeffT >::evaluate(), and pfasst::quadrature::Polynomial< CoeffT >::integrate().

Referenced by pfasst::quadrature::IQuadrature< precision >::compute_weights().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

template<typename scalar >
static Matrix<scalar> pfasst::quadrature::compute_s_matrix ( const Matrix< scalar > &  q_mat)
static

Compute node-to-node quadrature matrix \( S \) from a given quadrature matrix \( Q \).

The \( S \) matrix provides a node-to-node quadrature where the \( i \)-th row of \( S \) represents a quadrature from the \( i-1 \)-th node to the \( i \)-th node.

The procedure is simply subtracting the \( i-1 \)-th row of \( Q \) from the \( i \)-th row of \( Q \).

Template Parameters
scalarprecision of quadrature (i.e. double)
Parameters
[in]q_mat\( Q \) matrix to compute \( S \) of
Returns
\( S \) matrix
Since
v0.3.0

Definition at line 158 of file interface.hpp.

Referenced by pfasst::quadrature::IQuadrature< precision >::compute_weights().

+ Here is the caller graph for this function:

template<typename scalar >
static Matrix<scalar> pfasst::quadrature::compute_s_matrix ( const vector< scalar > &  from,
const vector< scalar > &  to 
)
static

Compute node-to-node quadrature matrix \( S \) from two given sets of nodes.

Template Parameters
scalarprecision of quadrature (i.e. double)
Parameters
[in]fromfirst set of quadrature nodes
[in]tosecond set of quadrature nodes
Since
v0.3.0

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 181 of file interface.hpp.

References compute_q_matrix().

+ Here is the call graph for this function:

template<typename precision = pfasst::time_precision>
shared_ptr<IQuadrature<precision> > pfasst::quadrature::quadrature_factory ( const size_t  nnodes,
const QuadratureType  qtype 
)

Instantiates quadrature handler for given number of nodes and type descriptor.

Template Parameters
precisionnumerical type of the nodes (e.g. double)
Parameters
[in]nnodesnumber of quadrature nodes
[in]qtypetype descriptor of the quadrature
Returns
instance of pfasst::quadrature::IQuadrature of specified type with desired number of nodes
Exceptions
pfasst::ValueErrorif qtype is not a valid quadrature type descriptor

Definition at line 53 of file quadrature.hpp.

Referenced by pfasst::examples::advection_diffusion::run_mpi_pfasst(), pfasst::examples::scalar::run_scalar_sdc(), pfasst::examples::advection_diffusion::run_serial_mlsdc(), pfasst::examples::advection_diffusion::run_vanilla_sdc(), and pfasst::examples::vdp::run_vdp_sdc().

+ Here is the caller graph for this function: