Dendro  5.01
Dendro in Greek language means tree. The Dendro library is a large scale (262K cores on ORNL's Titan) distributed memory adaptive octree framework. The main goal of Dendro is to perform large scale multiphysics simulations efficeiently in mordern supercomputers. Dendro consists of efficient parallel data structures and algorithms to perform variational ( finite element) methods and finite difference mthods on 2:1 balanced arbitary adaptive octrees which enables the users to perform simulations raning from black holes (binary black hole mergers) to blood flow in human body, where applications ranging from relativity, astrophysics to biomedical engineering.
Functions
basis Namespace Reference

Functions

void jacobip (double alpha, double beta, unsigned int N, double *x, double *p, unsigned int np)
 
void gradjacobip (double alpha, double beta, int N, double *x, double *dp, unsigned int np)
 
void jacobiglq (double alpha, double beta, int N, double *x, double *w)
 
void jacobigq (double alpha, double beta, int N, double *x, double *w)
 
void lagrange (const double *x0, int N, int at, const double *x, double *px, int m)
 computes the Lagrange polynomials evalueated at x coords. More...
 

Detailed Description

Author
Milinda Fernando School of Computing, University of Utah. Contains implementation of Jacobi polynomials, and Jacobi Gauss-Lobatto quadrature. reference: Mangll implementation.

Function Documentation

◆ gradjacobip()

void basis::gradjacobip ( double  alpha,
double  beta,
int  N,
double *  x,
double *  dp,
unsigned int  np 
)

Evaluate the derivative of the N'th order Jacobi Polynomial of type (alpha, beta ) > -1.

This is C reimplemetation of the function in JacobiP from Nodal Discontinuous Galerkin Methods by Jan S. Hesthaven and Tim Warburton.

Parameters
[in]alphaJacobi polynomial parameter
[in]betaJacobi polynomial parameter
[in]NOrder of the polynomial
[in]npsize of x and p.
[in]xRow vector of locations to evaluate the polynomial.
[out]dpRow vector of the evaluated derivative of the N'th order polynomial
Note
The matrix P is assumed to have the correct size of (1, len(x)). Also then polynomials are normalized to be orthonormal.

◆ jacobiglq()

void basis::jacobiglq ( double  alpha,
double  beta,
int  N,
double *  x,
double *  w 
)

Compute the N'th order Gauss Lobatto quadrature points and weights.

Where (alpha, beta ) > -1 and (alpha + beta <> -1).

Parameters
[in]alphaJacobi polynomial parameter
[in]betaJacobi polynomial parameter
[in]NOrder of the polynomial
[out]xRow vector of Gauss Lobatto node locations
[out]wRow vector of Gauss Lobatto weights
Note
The matrices x and w are assumed to have the correct size of (1, N+1). Also the node location will be between [-1,1].

◆ jacobigq()

void basis::jacobigq ( double  alpha,
double  beta,
int  N,
double *  x,
double *  w 
)

Compute the N'th order Gauss quadrature points and weights.

Where (alpha, beta ) > -1.

This is C reimplemetation of the function in JacobiGQ from Nodal Discontinuous Galerkin Methods by Jan S. Hesthaven and Tim Warburton.

Parameters
[in]alphaJacobi polynomial parameter
[in]betaJacobi polynomial parameter
[in]NOrder of the polynomial
[out]xRow vector of Gauss node locations
[out]wRow vector of Gauss weights
Note
The matrices x and w are assumed to have the correct size of (1, N+1). Also the node locaiton will be between [-1,1].

◆ jacobip()

void basis::jacobip ( double  alpha,
double  beta,
unsigned int  N,
double *  x,
double *  p,
unsigned int  np 
)

Evaluate N'th order Jacobi Polynomial of type (alpha, beta ) > -1.

Also (alpha + beta <> -1).

This is C reimplemetation of the function in JacobiP from Nodal Discontinuous Galerkin Methods by Jan S. Hesthaven and Tim Warburton.

Parameters
[in]alphaJacobi polynomial parameter
[in]betaJacobi polynomial parameter
[in]NOrder of the polynomial
[in]npsize of x and p.
[in]xRow vector of locations to evaluate the polynomial.
[out]pRow vector of the evaluated N'th order polynomial
Note
The matrix P is assumed to have the correct size of (1, len(x)). Also then polynomials are normalized to be orthonormal.

◆ lagrange()

void basis::lagrange ( const double *  x0,
int  N,
int  at,
const double *  x,
double *  px,
int  m 
)

computes the Lagrange polynomials evalueated at x coords.

Parameters
x0nodal locations. N+1 points.
Norder of the Lagrange polynomial
xpoints which Lagrange evaluated at.
pxP(x)
msize of x points, (i.e. similar to px)