FE-Project
|
Module common / Polynominal. More...
Functions/Subroutines | |
real(rp) function, dimension(size(x), nord+1), public | polynominal_genlagrangepoly (nord, x_lgl, x) |
A function to obtain the Lagrange basis functions related to the Gauss-Legendre-Lobatto (GLL) points. | |
real(rp) function, dimension(nord+1, nord+1), public | polynominal_gendlagrangepoly_lglpt (nord, x_lgl) |
A function to obtain the differential values of Lagrange basis functions at the GLL points. | |
subroutine, public | polynominal_genlegendrepoly_sub (nord, x, p) |
A function to obtain the values of Legendre polynomials which are evaluated at arbitrary points. | |
real(rp) function, dimension(size(x), nord+1), public | polynominal_genlegendrepoly (nord, x) |
A function to obtain the values of Legendre polynomials which are evaluated at arbitrary points. | |
real(rp) function, dimension(size(x), nord+1), public | polynominal_gendlegendrepoly (nord, x, p) |
A function to obtain differential values of Legendre polynomials which are evaluated at arbitrary points. | |
real(rp) function, dimension(nord+1), public | polynominal_gengausslobattopt (nord) |
A function to calculate the Legendre-Gauss-Lobatto (LGL) points. | |
real(rp) function, dimension(nord+1), public | polynominal_gengausslobattoptintweight (nord) |
A function to calculate the Gauss-Lobbato weights. | |
real(rp) function, dimension(nord), public | polynominal_gengausslegendrept (nord) |
A function to calculate the Gauss-Legendre (GL) points. | |
real(rp) function, dimension(nord), public | polynominal_gengausslegendreptintweight (nord) |
A function to calculate the Gauss-Legendre (GL) weights. |
Module common / Polynominal.
real(rp) function, dimension(size(x), nord+1), public scale_polynominal::polynominal_genlagrangepoly | ( | integer, intent(in) | nord, |
real(rp), dimension(nord+1), intent(in) | x_lgl, | ||
real(rp), dimension(:), intent(in) | x ) |
A function to obtain the Lagrange basis functions related to the Gauss-Legendre-Lobatto (GLL) points.
Nord | Order of Lagrange polynomial |
x_lgl | Positions of GLL points |
x | Positions where the Lagrange basis functions are evaluated |
Definition at line 67 of file scale_polynominal.F90.
References polynominal_gendlegendrepoly(), and polynominal_genlegendrepoly().
Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), scale_element_quadrilateral::quadrilateralelement_init(), and scale_element_siacfilter::siac_filter_init().
real(rp) function, dimension(nord+1, nord+1), public scale_polynominal::polynominal_gendlagrangepoly_lglpt | ( | integer, intent(in) | nord, |
real(rp), dimension(nord+1), intent(in) | x_lgl ) |
A function to obtain the differential values of Lagrange basis functions at the GLL points.
Nord | Order of Lagrange polynomial |
x_lgl | Positions of GLL points |
Definition at line 126 of file scale_polynominal.F90.
References polynominal_genlegendrepoly().
Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_init().
subroutine, public scale_polynominal::polynominal_genlegendrepoly_sub | ( | integer, intent(in) | nord, |
real(rp), dimension(:), intent(in) | x, | ||
real(rp), dimension(size(x), nord+1), intent(out) | p ) |
A function to obtain the values of Legendre polynomials which are evaluated at arbitrary points.
Nord | Order of Lagrange polynomial |
x | Positions where the Legendre polynomials are evaluated |
P | Values of the Legendre polynomials at x |
Definition at line 169 of file scale_polynominal.F90.
Referenced by scale_file_common_meshfield::file_common_meshfield_put_field1d_cartesbuf(), scale_file_common_meshfield::file_common_meshfield_put_field2d_cartesbuf(), scale_file_common_meshfield::file_common_meshfield_put_field3d_cartesbuf(), and polynominal_genlegendrepoly().
real(rp) function, dimension(size(x), nord+1), public scale_polynominal::polynominal_genlegendrepoly | ( | integer, intent(in) | nord, |
real(rp), dimension(:), intent(in) | x ) |
A function to obtain the values of Legendre polynomials which are evaluated at arbitrary points.
Nord | Order of Lagrange polynomial |
x | Positions where the Legendre polynomials are evaluated |
P | Values of the Legendre polynomials at x |
Definition at line 200 of file scale_polynominal.F90.
References polynominal_genlegendrepoly_sub().
Referenced by scale_file_common_meshfield::file_common_meshfield_put_field2d_cubedsphere_cartesbuf(), scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), polynominal_gendlagrangepoly_lglpt(), polynominal_gengausslegendreptintweight(), polynominal_gengausslobattoptintweight(), polynominal_genlagrangepoly(), scale_element_quadrilateral::quadrilateralelement_init(), and scale_element_siacfilter::siac_filter_init().
real(rp) function, dimension(size(x), nord+1), public scale_polynominal::polynominal_gendlegendrepoly | ( | integer, intent(in) | nord, |
real(rp), dimension(:), intent(in) | x, | ||
real(rp), dimension(:,:), intent(in) | p ) |
A function to obtain differential values of Legendre polynomials which are evaluated at arbitrary points.
Nord | Order of Lagrange polynomial |
x | Positions where the Legendre polynomials are evaluated |
P | Values of the Legendre polynomials at x |
Definition at line 219 of file scale_polynominal.F90.
Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), polynominal_gengausslegendreptintweight(), polynominal_genlagrangepoly(), and scale_element_quadrilateral::quadrilateralelement_init().
real(rp) function, dimension(nord+1), public scale_polynominal::polynominal_gengausslobattopt | ( | integer, intent(in) | nord | ) |
A function to calculate the Legendre-Gauss-Lobatto (LGL) points.
Nord | Order of Lagrange polynomial |
pts | Position of the LGL points |
Definition at line 250 of file scale_polynominal.F90.
Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), polynominal_gengausslobattoptintweight(), and scale_element_quadrilateral::quadrilateralelement_init().
real(rp) function, dimension(nord+1), public scale_polynominal::polynominal_gengausslobattoptintweight | ( | integer, intent(in) | nord | ) |
A function to calculate the Gauss-Lobbato weights.
Nord | Order of Lagrange polynomial |
int_weight_lgl | Gauss-Lobbato weights |
Definition at line 271 of file scale_polynominal.F90.
References polynominal_gengausslobattopt(), and polynominal_genlegendrepoly().
Referenced by scale_atm_dyn_dgm_trcadvect3d_heve::atm_dyn_dgm_trcadvect3d_heve_init(), scale_atm_phy_mp_dgm_common::atm_phy_mp_dgm_common_gen_intweight(), scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_init().
real(rp) function, dimension(nord), public scale_polynominal::polynominal_gengausslegendrept | ( | integer, intent(in) | nord | ) |
A function to calculate the Gauss-Legendre (GL) points.
Nord | Order of the Legendre polynomial |
pts | Position of the GL points |
Definition at line 294 of file scale_polynominal.F90.
Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), polynominal_gengausslegendreptintweight(), and scale_element_quadrilateral::quadrilateralelement_init().
real(rp) function, dimension(nord), public scale_polynominal::polynominal_gengausslegendreptintweight | ( | integer, intent(in) | nord | ) |
A function to calculate the Gauss-Legendre (GL) weights.
Nord | Order of the Legendre polynomial |
int_weight_gl | Gauss-Legendre weights |
Definition at line 310 of file scale_polynominal.F90.
References polynominal_gendlegendrepoly(), polynominal_gengausslegendrept(), and polynominal_genlegendrepoly().
Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_init().