ply_split_legendre_matrix Function

public pure function ply_split_legendre_matrix(nModes) result(split_matrix)

Compute the transformation matrix for a projection to the left and right half-interval of Legendre polynomials for the given maximal number of modes.

Note: The transformation matrices to each subinterval are triangular, and the diagonal entries are the same. To save memory both matrices are stored in a single 2 dimensional array of size (nModes, nModes).

This matrix only needs to be computed once for a sufficiently high order, as submatices out of it can by used to perform the transformation for any lower polynomial degree.

The upper triangular matrix is created for the right subinterval, while the lower triangular matrix is used to store the rotated version for the left subinterval. For the right interval we interpret the first index as row index and the second as column. For the left interval this is reverted and we interpret the first index as columns of the matrix.

Note

Why is this a function? The reasoning for making this a function is that we need to return exactly one thing (the split matrix). It is then quite natural to refer to this by ply_split_legendre_matrix. A subroutine on the other hand usually describes something that should be done. Thus the name for a subroutine would then be ply_split_legendre_compute_matrix (describing the action performed by the subroutine). When using OpenMP it sometimes is better to use subroutines, even though it would be more natural to use a function. However, here we do not expect this to be the case, as this is expected to be called only once.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nModes

The maximal number of modes to compute the transformation for.

The resulting matrix v will be max_modes x max_modes large and can be used for the transformation of all polynomials with up to this many modes.

Return Value real(kind=rk), (nModes,nModes)


Calls

proc~~ply_split_legendre_matrix~~CallsGraph proc~ply_split_legendre_matrix ply_split_legendre_matrix proc~alpha alpha proc~ply_split_legendre_matrix->proc~alpha proc~alpha_beta alpha_beta proc~ply_split_legendre_matrix->proc~alpha_beta proc~alpha_frac alpha_frac proc~ply_split_legendre_matrix->proc~alpha_frac proc~beta beta proc~ply_split_legendre_matrix->proc~beta

Called by

proc~~ply_split_legendre_matrix~~CalledByGraph proc~ply_split_legendre_matrix ply_split_legendre_matrix proc~ply_split_element_init ply_split_element_init proc~ply_split_element_init->proc~ply_split_legendre_matrix proc~ply_split_legendre_test ply_split_legendre_test proc~ply_split_legendre_test->proc~ply_split_legendre_matrix program~ply_split_legendre_test_prog ply_split_legendre_test_prog program~ply_split_legendre_test_prog->proc~ply_split_legendre_matrix program~ply_split_legendre_test_prog->proc~ply_split_legendre_test proc~ply_sample_adaptive ply_sample_adaptive proc~ply_sample_adaptive->proc~ply_split_element_init proc~ply_split_element_1d_test ply_split_element_1D_test proc~ply_split_element_1d_test->proc~ply_split_element_init proc~ply_split_element_2d_test ply_split_element_2D_test proc~ply_split_element_2d_test->proc~ply_split_element_init proc~ply_split_element_3d_test ply_split_element_3D_test proc~ply_split_element_3d_test->proc~ply_split_element_init proc~ply_split_element_test ply_split_element_test proc~ply_split_element_test->proc~ply_split_element_init proc~ply_split_element_test->proc~ply_split_element_1d_test proc~ply_split_element_test->proc~ply_split_element_2d_test proc~ply_split_element_test->proc~ply_split_element_3d_test proc~ply_sample_data ply_sample_data proc~ply_sample_data->proc~ply_sample_adaptive program~ply_split_element_test_prog ply_split_element_test_prog program~ply_split_element_test_prog->proc~ply_split_element_test proc~ply_sampled_track_output ply_sampled_track_output proc~ply_sampled_track_output->proc~ply_sample_data