ply_l2p_header_module Module

Parameters for the plain L2 projection method to transform between Legendre modes and nodal representation.

This method utilizes the L2 projection from Legendre to Lagrange polynomials or the other way around. A numerical Gauss-Legendre Quadrature is used to compute the integral over the product of both functions. The Lagrange polynomials can be defined on any nodeset, see also ply_nodeset_module.

Available options for the nodes to project onto are:

  • 'gauss-legendre' these are the Gauss-Legendre integration points that are also used for the numerical integration (this is the default).
  • 'chebyshev' these are the nodes from the Chebyshev integration.

The set of nodes to use is configured by the nodes_kind option, and if nodes_kind = 'chebyshev' it is also possible to make use of Lobatto points to include the interval boundaries in the nodal representation. This is achieved by setting lobattoPoints = true, by default this is false.

The configuration table for a projection with L2P may, for example, look as follows:

  projection = {
    kind = 'l2p',
    factor = 1.5,
    nodes_kind = 'chebyshev',
    lobattoPoints = true
  }

The example illustrates the three possible settings for the L2P transformation method:

  • factor - Oversampling factor to avoid aliasing.
  • nodes_kind - Selection of set of nodes to use in the nodal representation.
  • lobattoPoints - Whether to include interval bounds, only available for Chebyshev nodes.

Uses

Used by

  • module~~ply_l2p_header_module~~UsedByGraph module~ply_l2p_header_module ply_l2p_header_module module~ply_l2p_module ply_l2p_module module~ply_l2p_module->module~ply_l2p_header_module module~ply_prj_header_module ply_prj_header_module module~ply_prj_header_module->module~ply_l2p_header_module program~ply_l2p_test ply_l2p_test program~ply_l2p_test->module~ply_l2p_header_module program~ply_l2p_test->module~ply_l2p_module module~ply_dynarray_project_module ply_dynarray_project_module module~ply_dynarray_project_module->module~ply_prj_header_module module~ply_poly_project_module ply_poly_project_module module~ply_poly_project_module->module~ply_l2p_module module~ply_poly_project_module->module~ply_prj_header_module program~ply_project_2d_fpt_lobattopoints_test ply_project_2d_fpt_lobattoPoints_test program~ply_project_2d_fpt_lobattopoints_test->module~ply_prj_header_module program~ply_project_2d_fpt_test ply_project_2d_fpt_test program~ply_project_2d_fpt_test->module~ply_prj_header_module program~ply_project_fpt_lobattopoints_test ply_project_fpt_lobattoPoints_test program~ply_project_fpt_lobattopoints_test->module~ply_prj_header_module program~ply_project_fpt_test ply_project_fpt_test program~ply_project_fpt_test->module~ply_prj_header_module program~test_fxtd_n2m2n test_fxtd_n2m2n program~test_fxtd_n2m2n->module~ply_prj_header_module

Interfaces

public interface assignment(=)

public interface operator(==)

  • private pure function isEqual(left, right) result(equality)

    This function provides the test for equality of two projections.

    Two l2p header are considered to be equal, if their node_header, and the factor are equal.

    Arguments

    Type IntentOptional Attributes Name
    type(ply_l2p_header_type), intent(in) :: left

    projection to compare

    type(ply_l2p_header_type), intent(in) :: right

    projection to compare against

    Return Value logical

    is equal??

public interface operator(/=)

  • private pure function isUnequal(left, right) result(unequality)

    This function provides the test for unequality of two projections.

    Two l2p header are considered to be unequal, if their node_header, or the factor are not equal.

    Arguments

    Type IntentOptional Attributes Name
    type(ply_l2p_header_type), intent(in) :: left

    projection to compare

    type(ply_l2p_header_type), intent(in) :: right

    projection to compare against

    Return Value logical

    is unequal??

public interface operator(<)

  • private pure function isSmaller(left, right) result(small)

    This function provides a < comparison of two projections.

    Sorting of l2p header is given by node_header and by the factor.

    Arguments

    Type IntentOptional Attributes Name
    type(ply_l2p_header_type), intent(in) :: left

    projection to compare

    type(ply_l2p_header_type), intent(in) :: right

    projection to compare against

    Return Value logical

    is smaller??

public interface operator(<=)

  • private pure function isSmallerOrEqual(left, right) result(small)

    This function provides a <= comparison of two projections.

    Sorting of l2p header is given by node_header, l2p_blocksize and last by factor.

    Arguments

    Type IntentOptional Attributes Name
    type(ply_l2p_header_type), intent(in) :: left

    projection to compare

    type(ply_l2p_header_type), intent(in) :: right

    projection to compare against

    Return Value logical

    is smaller??

public interface operator(>)

  • private pure function isGreater(left, right) result(great)

    This function provides a > comparison of two projections.

    Sorting of l2p header is given by node_header, l2p_blocksize and last by factor.

    Arguments

    Type IntentOptional Attributes Name
    type(ply_l2p_header_type), intent(in) :: left

    projection to compare

    type(ply_l2p_header_type), intent(in) :: right

    projection to compare against

    Return Value logical

    is greater??

public interface operator(>=)

  • private pure function isGreaterOrEqual(left, right) result(great)

    This function provides a >= comparison of two projections.

    Sorting of l2p header is given by node_header, l2p_blocksize and last by factor.

    Arguments

    Type IntentOptional Attributes Name
    type(ply_l2p_header_type), intent(in) :: left

    projection to compare

    type(ply_l2p_header_type), intent(in) :: right

    projection to compare against

    Return Value logical

    is greater??


Derived Types

type, public ::  ply_l2p_header_type

l2p projection header type, consisting of the node header which give information about the type and number of points for the projection

Components

Type Visibility Attributes Name Initial
type(ply_nodes_header_type), public :: nodes_header
real(kind=rk), public :: factor

Functions

private pure function isEqual(left, right) result(equality)

This function provides the test for equality of two projections.

Read more…

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(in) :: left

projection to compare

type(ply_l2p_header_type), intent(in) :: right

projection to compare against

Return Value logical

is equal??

private pure function isUnequal(left, right) result(unequality)

This function provides the test for unequality of two projections.

Read more…

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(in) :: left

projection to compare

type(ply_l2p_header_type), intent(in) :: right

projection to compare against

Return Value logical

is unequal??

private pure function isSmaller(left, right) result(small)

This function provides a < comparison of two projections.

Read more…

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(in) :: left

projection to compare

type(ply_l2p_header_type), intent(in) :: right

projection to compare against

Return Value logical

is smaller??

private pure function isSmallerOrEqual(left, right) result(small)

This function provides a <= comparison of two projections.

Read more…

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(in) :: left

projection to compare

type(ply_l2p_header_type), intent(in) :: right

projection to compare against

Return Value logical

is smaller??

private pure function isGreater(left, right) result(great)

This function provides a > comparison of two projections.

Read more…

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(in) :: left

projection to compare

type(ply_l2p_header_type), intent(in) :: right

projection to compare against

Return Value logical

is greater??

private pure function isGreaterOrEqual(left, right) result(great)

This function provides a >= comparison of two projections.

Read more…

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(in) :: left

projection to compare

type(ply_l2p_header_type), intent(in) :: right

projection to compare against

Return Value logical

is greater??


Subroutines

public subroutine ply_l2p_header_load(me, conf, thandle)

Load settings to describe a projection method from a Lua table.

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(out) :: me
type(flu_State), intent(inout) :: conf
integer, intent(in) :: thandle

public subroutine ply_l2p_header_define(me, factor, nodes_kind, lobattoPoints)

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(out) :: me

L2P header to define.

integer, intent(in), optional :: factor

Oversampling factor to use in the projection, defaults to 1.

character(len=*), optional :: nodes_kind

Set of nodes to use in the nodal representation.

Read more…
logical, intent(in), optional :: lobattoPoints

Wether to use Lobatto points (include interval bounds) when using the chebyshev nodes, defaults to .false..

public subroutine ply_l2p_header_out(me, conf)

Write L2P settings into a Lua table.

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(in) :: me
type(aot_out_type), intent(inout) :: conf

public subroutine ply_l2p_header_display(me)

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(in) :: me

private pure subroutine Copy_l2p_header(left, right)

Arguments

Type IntentOptional Attributes Name
type(ply_l2p_header_type), intent(out) :: left

fpt to copy to

type(ply_l2p_header_type), intent(in) :: right

fpt to copy from