QuadraturePoint_Line
This routine returns the quadrature points on the line.
Calling example:
ans = QuadraturePoint_Line1( &
& order, &
& quadType, &
& layout, &
& xij, &
& alpha, &
& beta, &
& lambda)
ansThe last row of ans contains the weights. The first few rows contain the quadrature points. If xij is absent then ans has two rows. If xij is present then ans has SIZE(xij, 1)+1 rows.
orderis order of integrandxijcontains nodal coordinates of line in xij format.- SIZE(xij,1) = nsd, and SIZE(xij,2)=2
- If xij is absent then [-1,1] is used
quadTypeis interpolation point type, it can take following valuesGaussLegendreGaussLegendreLobattoGaussLegendreRadauLeftGaussLegendreRadauRightGaussChebyshevGaussChebyshevLobattoGaussChebyshevRadauLeftGaussChebyshevRadauRightGaussJacobiGaussJacobiLobattoGaussJacobiRadauLeftGaussJacobiRadauRightGaussUltrasphericalGaussUltrasphericalLobattoGaussUltrasphericalRadauLeftGaussUltrasphericalRadauRight
layoutspecifies the arrangement of points. Following options are possible:layout=VEFCvertex, edge, face, cell, in this case first two points are boundary points, remaining (from 3 to n) are internal points in increasing order.layout=INCREASINGpoints are arranged in increasing order
Interface
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE QuadraturePoint_Line
MODULE FUNCTION QuadraturePoint_Line1( &
& order, &
& quadType, &
& layout, &
& xij, &
& alpha, &
& beta, &
& lambda) RESULT(ans)
!!
INTEGER(I4B), INTENT(IN) :: order
!! Order of interpolation
INTEGER(I4B), INTENT(IN) :: quadType
!! Quadrature point type
!! Equidistance,
!! GaussLegendre,
!! GaussLegendreLobatto,
!! GaussChebyshev,
!! GaussChebyshevLobatto,
!! GaussJacobi,
!! GaussJacobiLobatto
CHARACTER(*), INTENT(IN) :: layout
!! "VEFC"
!! "INCREASING"
REAL(DFP), OPTIONAL, INTENT(IN) :: xij(:, :)
!! domain of interpolation
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda
!! Ultraspherical parameter
REAL(DFP), ALLOCATABLE :: ans(:, :)
!! quadrature points
!! If xij is present then the number of rows in ans
!! is same as size(xij,1) + 1.
!! If xij is not present then the number of rows in
!! ans is 2
!! The last row of ans contains the weights
!! The first few rows contains the quadrature points
END FUNCTION QuadraturePoint_Line1
END INTERFACE QuadraturePoint_Line
program main
use easifemBase
implicit none
integer(i4b) :: order
real(dfp) :: x
real(dfp), allocatable :: xij(:,:), coeff(:,:), ans(:, :)
character( len = * ), parameter :: layout="INCREASING"
integer(i4b) :: quadType
order = 4_I4B
quadType = GaussLegendre
ans = QuadraturePoint_Line(order=order, quadType=quadType, layout=layout)
call display(mdencode(ans), "GaussLegendre: " //char_lf//char_lf )
quadType = GaussLegendreLobatto
ans = QuadraturePoint_Line(order=order, quadType=quadType, layout=layout, xij=xij)
call display(mdencode(ans), "GaussLegendreLobatto: " //char_lf//char_lf)
quadType = GaussChebyshev
ans = QuadraturePoint_Line(order=order, quadType=quadType, layout=layout, xij=xij)
call display(mdencode(ans), "GaussChebyshev: "//char_lf//char_lf)
quadType = GaussChebyshevLobatto
ans = QuadraturePoint_Line(order=order, quadType=quadType, layout=layout, xij=xij)
call display(mdencode(ans), "GaussChebyshevLobatto: "//char_lf//char_lf)
end program main
See results
GaussLegendre:
| -0.7746 | 3.71231E-16 | 0.7746 |
| 0.55556 | 0.88889 | 0.55556 |
GaussLegendreLobatto:
| -1 | -0.44721 | 0.44721 | 1 |
| 0.16667 | 0.83333 | 0.83333 | 0.16667 |
GaussChebyshev:
| -0.86603 | 1.03412E-13 | 0.86603 |
| 1.0472 | 1.0472 | 1.0472 |
GaussChebyshevLobatto:
| -1 | -0.5 | 0.5 | 1 |
| 0.5236 | 1.0472 | 1.0472 | 0.5236 |
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE QuadraturePoint_Line
MODULE FUNCTION QuadraturePoint_Line2( &
& order, &
& quadType, &
& xij, &
& layout, &
& alpha, &
& beta, &
& lambda) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: order
!! order of interpolation
INTEGER(I4B), INTENT(IN) :: quadType
!! Quadrature point type
!! Equidistance
!! GaussLegendre
!! GaussLegendreLobatto
!! GaussChebyshev,
!! GaussChebyshevLobatto
!! GaussJacobi
!! GaussJacobiLobatto
REAL(DFP), INTENT(IN) :: xij(2)
!! end points
CHARACTER(*), INTENT(IN) :: layout
!! "VEFC"
!! "INCREASING"
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda
!! Ultraspherical parameter
REAL(DFP), ALLOCATABLE :: ans(:,:)
END FUNCTION QuadraturePoint_Line2
END INTERFACE QuadraturePoint_Line
See above example.