Block-Structured AMR Software Framework
sdcquadrature_mod Module Reference

Module to create quadrature matrices and accompanying routines for SDC. More...

Data Types

interface  poly_eval
 

Functions/Subroutines

subroutine sdc_quadrature (qtype_in, nnodes, nnodes0, nodes, nflags, qmats)
 Subroutine to create quadrature matrices. More...
 
logical function not_proper (flags, node)
 Function to decide if the restriction of the nodes is pointwise, e.g. coarse nodes are every other fine node. More...
 
subroutine sdc_qnodes (qnodes, flags, qtype, nnodes)
 Subroutine to compute high precision quadrature nodes. More...
 
subroutine sdc_qmats (qmat, smat, dst, src, flags, ndst, nsrc)
 Subroutine to compute the quadrature matrices. More...
 
real(qp) function poly_eval (p, n, x)
 Polynomial manipulation routines. More...
 
complex(qp) function poly_eval_complex (p, n, x)
 Function to evaluate complex polynomial. More...
 
subroutine poly_diff (p, n)
 Subroutine to differentiate polynomial (in place) More...
 
subroutine poly_int (p, n)
 Subroutine to integrate polynomial (in place) More...
 
subroutine poly_legendre (p, n)
 Subroutine to compute Legendre polynomial coefficients using Bonnet's recursion formula. More...
 
subroutine poly_roots (roots, p0, n)
 Subroutine to compute polynomial roots using the Durand-Kerner algorithm. The roots are assumed to be real. More...
 
recursive subroutine qsort (a)
 Subroutine to sort (inplace) using the quick sort algorithm. Adapted from http://www.fortran.com/qsort_c.f95. More...
 
subroutine, private qsort_partition (a, marker)
 

Variables

integer, parameter qp = c_long_double
 
integer, parameter dp = c_double
 
real(qp), parameter eps = 1.0e-23_qp
 
integer, parameter sdc_gauss_lobatto = 1
 Quadrature node types. More...
 
integer, parameter sdc_gauss_radau = 2
 
integer, parameter sdc_clenshaw_curtis = 3
 
integer, parameter sdc_uniform = 4
 
integer, parameter sdc_gauss_legendre = 5
 
integer, parameter sdc_proper_nodes = 2**8
 
integer, parameter sdc_composite_nodes = 2**9
 
integer, parameter sdc_no_left = 2**10
 
integer, parameter amrex_real = c_double
 

Detailed Description

Module to create quadrature matrices and accompanying routines for SDC.

Function/Subroutine Documentation

◆ not_proper()

logical function sdcquadrature_mod::not_proper ( integer(c_int), dimension(:), intent(in)  flags,
integer, intent(in)  node 
)

Function to decide if the restriction of the nodes is pointwise, e.g. coarse nodes are every other fine node.

◆ poly_diff()

subroutine sdcquadrature_mod::poly_diff ( real(qp), dimension(0:n), intent(inout)  p,
integer, intent(in), value  n 
)

Subroutine to differentiate polynomial (in place)

◆ poly_eval()

real(qp) function sdcquadrature_mod::poly_eval ( real(qp), dimension(0:n), intent(in)  p,
integer, intent(in), value  n,
real(qp), intent(in)  x 
)

Polynomial manipulation routines.

A polynomial p

p(x) = a_n x^n + ... + a_2 x^2 + a_1 x + a_0

is stored as a Fortran array p(0:n) according to

p = [ a_0, a_1, ..., a_n ]. Function to evaluate real polynomial

◆ poly_eval_complex()

complex(qp) function sdcquadrature_mod::poly_eval_complex ( real(qp), dimension(0:n), intent(in)  p,
integer, intent(in), value  n,
complex(qp), intent(in)  x 
)

Function to evaluate complex polynomial.

◆ poly_int()

subroutine sdcquadrature_mod::poly_int ( real(qp), dimension(0:n), intent(inout)  p,
integer, intent(in), value  n 
)

Subroutine to integrate polynomial (in place)

◆ poly_legendre()

subroutine sdcquadrature_mod::poly_legendre ( real(qp), dimension(0:n), intent(out)  p,
integer, intent(in), value  n 
)

Subroutine to compute Legendre polynomial coefficients using Bonnet's recursion formula.

◆ poly_roots()

subroutine sdcquadrature_mod::poly_roots ( real(qp), dimension(n), intent(out)  roots,
real(qp), dimension(0:n), intent(in)  p0,
integer, intent(in), value  n 
)

Subroutine to compute polynomial roots using the Durand-Kerner algorithm. The roots are assumed to be real.

◆ qsort()

recursive subroutine sdcquadrature_mod::qsort ( real(qp), dimension(:), intent(inout)  a)

Subroutine to sort (inplace) using the quick sort algorithm. Adapted from http://www.fortran.com/qsort_c.f95.

◆ qsort_partition()

subroutine, private sdcquadrature_mod::qsort_partition ( real(qp), dimension(:), intent(inout)  a,
integer, intent(out)  marker 
)
private

◆ sdc_qmats()

subroutine sdcquadrature_mod::sdc_qmats ( real(amrex_real), dimension(ndst-1, nsrc), intent(out)  qmat,
real(amrex_real), dimension(ndst-1, nsrc), intent(out)  smat,
real(c_long_double), dimension(ndst), intent(in)  dst,
real(c_long_double), dimension(nsrc), intent(in)  src,
integer(c_int), dimension(nsrc), intent(in)  flags,
integer(c_int), intent(in), value  ndst,
integer(c_int), intent(in), value  nsrc 
)

Subroutine to compute the quadrature matrices.

Parameters
[in]ndstNumber of destination points
[in]nsrcNumber of source points
[in]dstDestination points
[in]srcSource points
[out]qmatO to dst quadrature weights
[out]smatdst(m) to dst(m+1) quadrature weights

◆ sdc_qnodes()

subroutine sdcquadrature_mod::sdc_qnodes ( real(c_long_double), dimension(nnodes), intent(out)  qnodes,
integer(c_int), dimension(nnodes), intent(out)  flags,
integer(c_int), intent(in), value  qtype,
integer(c_int), intent(in), value  nnodes 
)

Subroutine to compute high precision quadrature nodes.

Parameters
[in]nnodesNumber of nodes
[in]qtypeType of nodes (see pf_dtype)
[out]qnodesThe computed quadrature nodes

◆ sdc_quadrature()

subroutine sdcquadrature_mod::sdc_quadrature ( integer, intent(in)  qtype_in,
integer, intent(in)  nnodes,
integer, intent(in)  nnodes0,
real(amrex_real), dimension(nnodes), intent(inout)  nodes,
integer, dimension(nnodes), intent(out)  nflags,
real(amrex_real), dimension(nnodes,nnodes-1,4), intent(inout)  qmats 
)

Subroutine to create quadrature matrices.

Variable Documentation

◆ amrex_real

integer, parameter sdcquadrature_mod::amrex_real = c_double

◆ dp

integer, parameter sdcquadrature_mod::dp = c_double

◆ eps

real(qp), parameter sdcquadrature_mod::eps = 1.0e-23_qp

◆ qp

integer, parameter sdcquadrature_mod::qp = c_long_double

◆ sdc_clenshaw_curtis

integer, parameter sdcquadrature_mod::sdc_clenshaw_curtis = 3

◆ sdc_composite_nodes

integer, parameter sdcquadrature_mod::sdc_composite_nodes = 2**9

◆ sdc_gauss_legendre

integer, parameter sdcquadrature_mod::sdc_gauss_legendre = 5

◆ sdc_gauss_lobatto

integer, parameter sdcquadrature_mod::sdc_gauss_lobatto = 1

Quadrature node types.

◆ sdc_gauss_radau

integer, parameter sdcquadrature_mod::sdc_gauss_radau = 2

◆ sdc_no_left

integer, parameter sdcquadrature_mod::sdc_no_left = 2**10

◆ sdc_proper_nodes

integer, parameter sdcquadrature_mod::sdc_proper_nodes = 2**8

◆ sdc_uniform

integer, parameter sdcquadrature_mod::sdc_uniform = 4