Discontinuous Galerkin Library
#include "dg/algorithm.h"

\( f_i = f(\vec x_i) \) More...

Collaboration diagram for evaluate:

Functions

template<class UnaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate (UnaryOp f, const RealGrid1d< real_type > &g)
 Evaluate a 1d function on grid coordinates. More...
 
template<class BinaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate (const BinaryOp &f, const aRealTopology2d< real_type > &g)
 Evaluate a 2d function on grid coordinates. More...
 
template<class TernaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate (const TernaryOp &f, const aRealTopology3d< real_type > &g)
 Evaluate a 3d function on grid coordinates. More...
 
template<class real_type >
thrust::host_vector< real_type > dg::integrate (const thrust::host_vector< real_type > &in, const RealGrid1d< real_type > &g, dg::direction dir=dg::forward)
 Indefinite integral of a function on a grid. More...
 
template<class UnaryOp , class real_type >
thrust::host_vector< real_type > dg::integrate (UnaryOp f, const RealGrid1d< real_type > &g, dg::direction dir=dg::forward)
 Indefinite integral of a function on a grid. More...
 
template<class UnaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate (UnaryOp f, const RealGridX1d< real_type > &g)
 Evaluate a 1d function on grid coordinates. More...
 
template<class BinaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate (const BinaryOp &f, const aRealTopologyX2d< real_type > &g)
 Evaluate a 2d function on grid coordinates. More...
 
template<class TernaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate (const TernaryOp &f, const aRealTopologyX3d< real_type > &g)
 Evaluate a 3d function on grid coordinates. More...
 
template<class BinaryOp , class real_type >
MPI_Vector< thrust::host_vector< real_type > > dg::evaluate (const BinaryOp &f, const aRealMPITopology2d< real_type > &g)
 Evaluate a 2d function on mpi distributed grid coordinates. More...
 
template<class TernaryOp , class real_type >
MPI_Vector< thrust::host_vector< real_type > > dg::evaluate (const TernaryOp &f, const aRealMPITopology3d< real_type > &g)
 Evaluate a 3d function on mpi distributed grid coordinates. More...
 

Detailed Description

\( f_i = f(\vec x_i) \)

The function discretisation routines compute the dG discretisation of analytic functions on a given grid. In 1D the discretisation simply consists of n function values per grid cell ( where n is the number of Legendre coefficients used; currently \( 1 <= n <= 20\) ) evaluated at the Gaussian abscissas in the respective cell. In 2D and 3D we simply use the product space.

Note
We choose x to be the contiguous direction. E.g. in 2D the first element of the resulting vector lies in the grid corner \( (x_0,y_0)\) and the last in \((x_1, y_1)\) .

Function Documentation

◆ evaluate() [1/8]

template<class BinaryOp , class real_type >
MPI_Vector< thrust::host_vector< real_type > > dg::evaluate ( const BinaryOp &  f,
const aRealMPITopology2d< real_type > &  g 
)

Evaluate a 2d function on mpi distributed grid coordinates.

Evaluate is equivalent to the following:

  1. generate the list of grid coordinates \( x_i\), \( y_i\) representing the local part of the given computational space discretization (the grid) that the calling process owns
  2. evaluate the given function or functor at these coordinates and store the result in the output vector \( v_i = f(x_i, y_i)\) for all \( i\)

This code snippet demonstrates how to discretize and compute the norm of a function on a distributed memory system

//define the function to integrate
double function(double x, double y, double amp){
return amp*exp(x)*exp(y);
}
// create a grid of the domain [0,2]x[0,2] with 20 cells in x and y and 3 polynomial coefficients
//... use MPI_Cart_create to create 2d Cartesian communicator
dg::MPIGrid2d g2d( 0, 2, 0, 2, 3, 20, 20, comm2d);
// create the Gaussian weights (volume form) for the integration
const dg::MHVec w2d = dg::create::weights( g2d);
// our function needs to depend only on x and y so let's fix the amplitude to 2
using namespace std::placeholders; //for _1, _2
auto functor = std::bind( function, _1, _2, 2.);
// discretize the function on the grid
const dg::MHVec vec = dg::evaluate( functor, g2d);
// multiply and sum the results
double norm = dg::blas2::dot( vec, w2d, vec);
// norm is now: (e^4-1)^2
get_value_type< MatrixType > dot(const ContainerType1 &x, const MatrixType &m, const ContainerType2 &y)
; Binary reproducible general dot product
Definition: blas2.h:85
thrust::host_vector< real_type > evaluate(UnaryOp f, const RealGrid1d< real_type > &g)
Evaluate a 1d function on grid coordinates.
Definition: evaluation.h:67
MPI_Vector< thrust::host_vector< real_type > > weights(const aRealMPITopology2d< real_type > &g)
Nodal weight coefficients.
Definition: mpi_weights.h:22
mpi Vector class
Definition: mpi_vector.h:32
The simplest implementation of aRealMPITopology2d.
Definition: mpi_grid.h:691
Template Parameters
BinaryOpA class or function type with a member/signature equivalent to
  • real_type operator()(real_type, real_type) const
Parameters
fThe function to evaluate: f = f(x,y)
gThe 2d grid on which to evaluate f
Returns
The output vector v as an MPI host Vector
Note
Use the elementary function \( f(x,y) = x \) (dg::cooX2d) to generate the list of grid coordinates in x direction (or analogous in y, dg::cooY2d)
See also
Introduction to dg methods

◆ evaluate() [2/8]

template<class BinaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate ( const BinaryOp &  f,
const aRealTopology2d< real_type > &  g 
)

Evaluate a 2d function on grid coordinates.

Evaluate is equivalent to the following:

  1. generate the list of grid coordinates \( x_i\), \( y_i\) representing the given computational space discretization (the grid)
  2. evaluate the given function or functor at these coordinates and store the result in the output vector \( v_i = f(x_i, y_i)\) for all \( i \)

This code snippet demonstrates how to discretize and compute the norm of a function on a shared memory system

// define the function to integrate
double function(double x, double y, double amp){
return amp*exp(x)*exp(y);
}
// create a grid of the domain [0,2]x[0,2] with 20 cells in x and y and 3 polynomial coefficients
dg::Grid2d g2d( 0, 2, 0, 2, 3, 20, 20);
// create the Gaussian weights (volume form) for the integration
const dg::HVec w2d = dg::create::weights( g2d);
// our function needs to depend only on x and y so let's fix the amplitude to 2
using namespace std::placeholders; //for _1, _2
auto functor = std::bind( function, _1, _2, 2.);
// discretize the function on the grid
const dg::HVec vec = dg::evaluate( functor, g2d);
// now compute the scalar product (the L2 norm)
double norm = dg::blas2::dot( vec, w2d, vec);
// norm is now: (e^4-1)^2
thrust::host_vector< double > HVec
Host Vector.
Definition: typedefs.h:19
Template Parameters
BinaryOpA class or function type with a member/signature equivalent to
  • real_type operator()(real_type, real_type) const
Parameters
fThe function to evaluate: \( f = f(x,y)\), see A large collection for a host of predefined functors to evaluate
gThe 2d grid on which to evaluate f
Returns
The output vector v as a host vector
Note
Use the elementary function \( f(x,y) = x \) (dg::cooX2d) to generate the list of grid coordinates in x direction (or analogous in y, dg::cooY2d)
See also
Introduction to dg methods
dg::pullback if you want to evaluate a function in physical space

◆ evaluate() [3/8]

template<class BinaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate ( const BinaryOp &  f,
const aRealTopologyX2d< real_type > &  g 
)

Evaluate a 2d function on grid coordinates.

Evaluate is equivalent to the following:

  1. generate the list of grid coordinates \( x_i\), \( y_i\) representing the given computational space discretization (the grid)
  2. evaluate the given function or functor at these coordinates and store the result in the output vector \( v_i = f(x_i, y_i)\) for all \( i \)

This code snippet demonstrates how to discretize and compute the norm of a function on a shared memory system

// define the function to integrate
double function(double x, double y, double amp){
return amp*exp(x)*exp(y);
}
// create a grid of the domain [0,2]x[0,2] with 20 cells in x and y and 3 polynomial coefficients
dg::Grid2d g2d( 0, 2, 0, 2, 3, 20, 20);
// create the Gaussian weights (volume form) for the integration
const dg::HVec w2d = dg::create::weights( g2d);
// our function needs to depend only on x and y so let's fix the amplitude to 2
using namespace std::placeholders; //for _1, _2
auto functor = std::bind( function, _1, _2, 2.);
// discretize the function on the grid
const dg::HVec vec = dg::evaluate( functor, g2d);
// now compute the scalar product (the L2 norm)
double norm = dg::blas2::dot( vec, w2d, vec);
// norm is now: (e^4-1)^2
Template Parameters
BinaryOpA class or function type with a member/signature equivalent to
  • real_type operator()(real_type, real_type) const
Parameters
fThe function to evaluate: \( f = f(x,y)\), see A large collection for a host of predefined functors to evaluate
gThe 2d grid on which to evaluate f
Returns
The output vector v as a host vector
Note
Use the elementary function \( f(x,y) = x \) (dg::cooX2d) to generate the list of grid coordinates in x direction (or analogous in y, dg::cooY2d)
See also
Introduction to dg methods
dg::pullback if you want to evaluate a function in physical space

◆ evaluate() [4/8]

template<class TernaryOp , class real_type >
MPI_Vector< thrust::host_vector< real_type > > dg::evaluate ( const TernaryOp &  f,
const aRealMPITopology3d< real_type > &  g 
)

Evaluate a 3d function on mpi distributed grid coordinates.

Evaluate is equivalent to the following:

  1. generate the list of grid coordinates \( x_i\), \( y_i\), \( z_i \) representing the local part of the given computational space discretization (the grid) that the calling process owns
  2. evaluate the given function or functor at these coordinates and store the result in the output vector \( v_i = f(x_i, y_i, z_i)\) for all \( i\)

This code snippet demonstrates how to discretize and compute the norm of a function

// define the function to integrate
double function(double x, double y, double z, double amp){
return amp*exp(x)*exp(y)*exp(z);
}
// create a grid of the domain [0,2]x[0,2]x[0,2] with 20 cells in x, y and z and 3 polynomial coefficients and x and y and 1 in z
//... use MPI_Cart_create to create 3d Cartesian communicator
dg::MPIGrid3d g3d( 0, 2, 0, 2, 0, 2, 3, 20, 20, 20, comm3d);
// create the Gaussian weights (volume form) for the integration
const dg::MHVec w3d = dg::create::weights( g3d);
// our function needs to depend only on x, y and z so let's fix the amplitude to 2
using namespace std::placeholders; //for _1, _2, _3
auto functor = std::bind( function, _1, _2, _3, 2.);
// discretize the function on the grid
const dg::MHVec vec = dg::evaluate( functor, g3d);
// now compute the scalar product (the L2 norm)
double norm = dg::blas2::dot(vec, w3d, vec);
// norm is now: (exp(4)-exp(0))^3/2
The simplest implementation of aRealMPITopology3d.
Definition: mpi_grid.h:727
Template Parameters
TernaryOpA class or function type with a member/signature equivalent to
  • real_type operator()(real_type, real_type, real_type) const
Parameters
fThe function to evaluate: f = f(x,y,z)
gThe 3d grid on which to evaluate f
Returns
The output vector v as an MPI host Vector
Note
Use the elementary function \( f(x,y,z) = x \) (dg::cooX3d) to generate the list of grid coordinates in x direction (or analogous in y, dg::cooY3d or z, dg::cooZ3d)
See also
Introduction to dg methods

◆ evaluate() [5/8]

template<class TernaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate ( const TernaryOp &  f,
const aRealTopology3d< real_type > &  g 
)

Evaluate a 3d function on grid coordinates.

Evaluate is equivalent to the following:

  1. generate the list of grid coordinates \( x_i\), \( y_i\), \( z_i \) representing the given computational space discretization (the grid)
  2. evaluate the given function or functor at these coordinates and store the result in the output vector \( v_i = f(x_i, y_i, z_i)\) for all \( i\)

This code snippet demonstrates how to discretize and compute the norm of a function

// define the function to integrate
double function(double x, double y, double z, double amp){
return amp*exp(x)*exp(y)*exp(z);
}
// create a grid of the domain [0,2]x[0,2]x[0,2] with 20 cells in x, y and z and 3 polynomial coefficients and x and y and 1 in z
dg::Grid3d g3d( 0, 2, 0, 2, 0, 2, 3, 20, 20, 20);
// create the Gaussian weights (volume form) for the integration
const dg::HVec w3d = dg::create::weights( g3d);
// our function needs to depend only on x, y and z so let's fix the amplitude to 2
using namespace std::placeholders; //for _1, _2, _3
auto functor = std::bind( function, _1, _2, _3, 2.);
// discretize the function on the grid
const dg::HVec vec = dg::evaluate( functor, g3d);
// now compute the scalar product (the L2 norm)
double norm = dg::blas2::dot(vec, w3d, vec);
// norm is now: (exp(4)-exp(0))^3/2
The simplest implementation of aRealTopology3d.
Definition: grid.h:844
Template Parameters
TernaryOpA class or function type with a member/signature equivalent to
  • real_type operator()(real_type, real_type, real_type) const
Parameters
fThe function to evaluate: \( f = f(x,y,z) \), see A large collection for a host of predefined functors to evaluate
gThe 3d grid on which to evaluate f
Returns
The output vector v as a host vector
Note
Use the elementary function \( f(x,y,z) = x \) (dg::cooX3d) to generate the list of grid coordinates in x direction (or analogous in y, dg::cooY3d or z, dg::cooZ3d)
See also
Introduction to dg methods
dg::pullback if you want to evaluate a function in physical space

◆ evaluate() [6/8]

template<class TernaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate ( const TernaryOp &  f,
const aRealTopologyX3d< real_type > &  g 
)

Evaluate a 3d function on grid coordinates.

Evaluate is equivalent to the following:

  1. generate the list of grid coordinates \( x_i\), \( y_i\), \( z_i \) representing the given computational space discretization (the grid)
  2. evaluate the given function or functor at these coordinates and store the result in the output vector \( v_i = f(x_i, y_i, z_i)\) for all \( i\)

This code snippet demonstrates how to discretize and compute the norm of a function

// define the function to integrate
double function(double x, double y, double z, double amp){
return amp*exp(x)*exp(y)*exp(z);
}
// create a grid of the domain [0,2]x[0,2]x[0,2] with 20 cells in x, y and z and 3 polynomial coefficients and x and y and 1 in z
dg::Grid3d g3d( 0, 2, 0, 2, 0, 2, 3, 20, 20, 20);
// create the Gaussian weights (volume form) for the integration
const dg::HVec w3d = dg::create::weights( g3d);
// our function needs to depend only on x, y and z so let's fix the amplitude to 2
using namespace std::placeholders; //for _1, _2, _3
auto functor = std::bind( function, _1, _2, _3, 2.);
// discretize the function on the grid
const dg::HVec vec = dg::evaluate( functor, g3d);
// now compute the scalar product (the L2 norm)
double norm = dg::blas2::dot(vec, w3d, vec);
// norm is now: (exp(4)-exp(0))^3/2
Template Parameters
TernaryOpA class or function type with a member/signature equivalent to
  • real_type operator()(real_type, real_type, real_type) const
Parameters
fThe function to evaluate: \( f = f(x,y,z) \), see A large collection for a host of predefined functors to evaluate
gThe 3d grid on which to evaluate f
Returns
The output vector v as a host vector
Note
Use the elementary function \( f(x,y,z) = x \) (dg::cooX3d) to generate the list of grid coordinates in x direction (or analogous in y, dg::cooY3d or z, dg::cooZ3d)
See also
Introduction to dg methods
dg::pullback if you want to evaluate a function in physical space

◆ evaluate() [7/8]

template<class UnaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate ( UnaryOp  f,
const RealGrid1d< real_type > &  g 
)

Evaluate a 1d function on grid coordinates.

Evaluate is equivalent to the following:

  1. generate a list of grid coordinates \( x_i\) representing the given computational space discretization (the grid)
  2. evaluate the given function or functor at these coordinates and store the result in the output vector \( v_i = f(x_i)\) for all i

This code snippet demonstrates how to discretize and integrate a function

// create a one-dimensional grid on the domain [0,2] with 3 polynomial coefficients and 20 cells
dg::Grid1d g1d( 0, 2, 3, 20);
// create the Gaussian weights (volume form) for the integration
const dg::HVec w1d = dg::create::weights( g1d);
// discretize the exponential function on the grid
const dg::HVec vec = dg::evaluate( exp, g1d);
// now compute the scalar product (the integral)
double integral = dg::blas1::dot( w1d, vec);
// integral is now: e^2-1
get_value_type< ContainerType1 > dot(const ContainerType1 &x, const ContainerType2 &y)
Binary reproducible Euclidean dot product between two vectors
Definition: blas1.h:87
Template Parameters
UnaryOpModel of Unary Function real_type f(real_type)
Parameters
fThe function to evaluate, see A large collection for a host of predefined functors to evaluate
gThe grid that defines the computational space on which to evaluate f
Returns
The output vector v as a host vector
Note
Use the elementary function \( f(x) = x \) (dg::cooX1d() ) to generate the list of grid coordinates
See also
Introduction to dg methods
dg::pullback if you want to evaluate a function in physical space

◆ evaluate() [8/8]

template<class UnaryOp , class real_type >
thrust::host_vector< real_type > dg::evaluate ( UnaryOp  f,
const RealGridX1d< real_type > &  g 
)

Evaluate a 1d function on grid coordinates.

Evaluate is equivalent to the following:

  1. generate a list of grid coordinates \( x_i\) representing the given computational space discretization (the grid)
  2. evaluate the given function or functor at these coordinates and store the result in the output vector \( v_i = f(x_i)\) for all i

This code snippet demonstrates how to discretize and integrate a function

// create a one-dimensional grid on the domain [0,2] with 3 polynomial coefficients and 20 cells
dg::Grid1d g1d( 0, 2, 3, 20);
// create the Gaussian weights (volume form) for the integration
const dg::HVec w1d = dg::create::weights( g1d);
// discretize the exponential function on the grid
const dg::HVec vec = dg::evaluate( exp, g1d);
// now compute the scalar product (the integral)
double integral = dg::blas1::dot( w1d, vec);
// integral is now: e^2-1
Template Parameters
UnaryOpModel of Unary Function real_type f(real_type)
Parameters
fThe function to evaluate, see A large collection for a host of predefined functors to evaluate
gThe grid that defines the computational space on which to evaluate f
Returns
The output vector v as a host vector
Note
Use the elementary function \( f(x) = x \) (dg::cooX1d() ) to generate the list of grid coordinates
See also
Introduction to dg methods
dg::pullback if you want to evaluate a function in physical space

◆ integrate() [1/2]

template<class real_type >
thrust::host_vector< real_type > dg::integrate ( const thrust::host_vector< real_type > &  in,
const RealGrid1d< real_type > &  g,
dg::direction  dir = dg::forward 
)

Indefinite integral of a function on a grid.

\[ F_h(x) = \int_a^x f_h(x') dx' \]

This function computes the indefinite integral of a given input

Parameters
inHost vector discretized on g
gThe grid
dirIf dg::backward then the integral starts at the right boundary (i.e. goes in the reverse direction)

\[ F_h(x) = \int_b^x f_h(x') dx' = \int_a^x f_h(x') dx' - \int_a^b f_h(x') dx' \]

Returns
integral of in on the grid g
See also
Introduction to dg methods

◆ integrate() [2/2]

template<class UnaryOp , class real_type >
thrust::host_vector< real_type > dg::integrate ( UnaryOp  f,
const RealGrid1d< real_type > &  g,
dg::direction  dir = dg::forward 
)

Indefinite integral of a function on a grid.

\[ F_h(x) = \int_a^x f_h(x') dx' \]

This function first evaluates f on the given grid and then computes and returns its indefinite integral

Parameters
fThe function to evaluate and then integrate
gThe grid
dirIf dg::backward then the integral starts at the right boundary (i.e. goes in the reverse direction)

\[ F_h(x) = \int_b^x f_h(x') dx' = \int_a^x f_h(x') dx' - \int_a^b f_h(x') dx' \]

Returns
integral of f on the grid g
See also
Introduction to dg methods