QuadraturePoint
This routine returns the quadrature points on Quadrangle.
Calling example:
program main
ans = QuadraturePoint_Quadrangle(order, quadType, &
& xij, alpha, beta, lambda)
end program main
Interface
- ܀ Interface
- ️܀ GaussLegendre
- GaussLegendreLobatto
- GaussLegendreRadauLeft
- ↢
INTERFACE QuadraturePoint_Quadrangle
MODULE FUNCTION QuadraturePoint_Quadrangle1(order, quadType, &
& xij, alpha, beta, lambda) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: order
!! order of integrand in x and y direction
INTEGER(I4B), INTENT(IN) :: quadType
!! Quadrature point type
REAL(DFP), OPTIONAL, INTENT(IN) :: xij(:, :)
!! four vertices of quadrangle in xij format
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(:, :)
!! interpolation points in xij format
END FUNCTION QuadraturePoint_Quadrangle1
END INTERFACE QuadraturePoint_Quadrangle
xij
xij
contains nodal coordinates of quadrangle in xij format.- The number of rows in xij can be 2 or 3 (for 2D or 3D)
- The number of columns in xij is 4
- If xij is absent then biunit quadrangle is assumed.
quadType
quadType
is interpolation point type, it can take following valuesGaussLegendre
GaussLegendreLobatto
GaussLegendreRadauLeft
GaussLegendreRadauRight
GaussChebyshev
GaussChebyshevLobatto
GaussChebyshevRadauLeft
GaussChebyshevRadauRight
GaussJacobi
GaussJacobiLobatto
GaussJacobiRadauLeft
GaussJacobiRadauRight
GaussUltraspherical
GaussUltrasphericalLobatto
GaussUltrasphericalRadauLeft
GaussUltrasphericalRadauRight
alpha, beta, lambda
alpha
andbeta
are parameters for Jacobi quadrature pointslambda
is the parameter for Ultraspherical quadrature points
program main
use easifembase
implicit none
integer( i4b ) :: i1, i2, order, quadType
real( dfp ), allocatable :: x(:,:), xij(:, :)
order=2
xij = RefQuadrangleCoord("BIUNIT")
quadType=GaussLegendre
x =QuadraturePoint_Quadrangle( &
& order=order, &
& quadType=quadType)
call display( Mdencode(x), "xij (order="//tostring(order)//")=" )
end program main
See results
xij (order=2)=
-0.57735 | -0.57735 | 0.57735 | 0.57735 | |
-0.57735 | 0.57735 | -0.57735 | 0.57735 | |
w | 1 | 1 | 1 | 1 |
program main
use easifembase
implicit none
integer( i4b ) :: i1, i2, order, quadType
real( dfp ), allocatable :: x(:,:), xij(:, :)
order=2
xij = RefQuadrangleCoord("BIUNIT")
quadType=GaussLegendreLobatto
x =QuadraturePoint_Quadrangle( &
& order=order, &
& quadType=quadType)
call display( Mdencode(x), "xij (order="//tostring(order)//")=" // char_lf // char_lf )
end program main
See results
xij (order=2)=
-1 | -1 | -1 | -2.66578E-17 | -2.66578E-17 | -2.66578E-17 | 1 | 1 | 1 |
-1 | -2.66578E-17 | 1 | -1 | -2.66578E-17 | 1 | -1 | -2.66578E-17 | 1 |
0.11111 | 0.44444 | 0.11111 | 0.44444 | 1.7778 | 0.44444 | 0.11111 | 0.44444 | 0.11111 |
etails>
program main
use easifembase
implicit none
integer( i4b ) :: i1, i2, order, quadType
real( dfp ), allocatable :: x(:,:), xij(:, :)
order=2
xij = RefQuadrangleCoord("BIUNIT")
quadType=GaussLegendreRadauLeft
x =QuadraturePoint_Quadrangle( &
& order=order, &
& quadType=quadType)
call display( Mdencode(x), "xij (order="//tostring(order)//")=" // char_lf // char_lf )
end program main
See results
xij (order=2)=
-1 | -1 | -1 | -2.66578E-17 | -2.66578E-17 | -2.66578E-17 | 1 | 1 | 1 |
-1 | -2.66578E-17 | 1 | -1 | -2.66578E-17 | 1 | -1 | -2.66578E-17 | 1 |
0.11111 | 0.44444 | 0.11111 | 0.44444 | 1.7778 | 0.44444 | 0.11111 | 0.44444 | 0.11111 |
etails>
Interface 2
Interface
INTERFACE QuadraturePoint_Quadrangle
MODULE FUNCTION QuadraturePoint_Quadrangle2( &
& p, q, quadType1, quadType2, xij, alpha1, beta1, &
& lambda1, alpha2, beta2, lambda2) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: p
!! order of integrand in x direction
INTEGER(I4B), INTENT(IN) :: q
!! order of integrand in y direction
INTEGER(I4B), INTENT(IN) :: quadType1
!! quadrature point type in x direction
INTEGER(I4B), INTENT(IN) :: quadType2
!! quadrature point type in y direction
REAL(DFP), OPTIONAL, INTENT(IN) :: xij(:, :)
!! four vertices of quadrangle in xij format
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha1
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta1
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda1
!! Ultraspherical parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha2
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta2
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda2
!! Ultraspherical parameter
REAL(DFP), ALLOCATABLE :: ans(:, :)
!! interpolation points in xij format
END FUNCTION QuadraturePoint_Quadrangle2
END INTERFACE QuadraturePoint_Quadrangle
Interface 3
INTERFACE QuadraturePoint_Quadrangle
MODULE FUNCTION QuadraturePoint_Quadrangle3(nips, quadType, &
& xij, alpha, beta, lambda) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: nips(1)
!! number of integration points in x and y direction
INTEGER(I4B), INTENT(IN) :: quadType
!! quadrature point type
REAL(DFP), OPTIONAL, INTENT(IN) :: xij(:, :)
!! four vertices of quadrangle in xij format
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(:, :)
!! interpolation points in xij format
END FUNCTION QuadraturePoint_Quadrangle3
END INTERFACE QuadraturePoint_Quadrangle
Interface 4
INTERFACE QuadraturePoint_Quadrangle
MODULE FUNCTION QuadraturePoint_Quadrangle4( &
& nipsx, nipsy, quadType1, quadType2, xij, alpha1, beta1, &
& lambda1, alpha2, beta2, lambda2) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: nipsx(1)
!! order of integrand in x direction
INTEGER(I4B), INTENT(IN) :: nipsy(1)
!! order of integrand in y direction
INTEGER(I4B), INTENT(IN) :: quadType1
!! quadrature point type in x direction
INTEGER(I4B), INTENT(IN) :: quadType2
!! quadrature point type in y direction
REAL(DFP), OPTIONAL, INTENT(IN) :: xij(:, :)
!! four vertices of quadrangle in xij format
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha1
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta1
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda1
!! Ultraspherical parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha2
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta2
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda2
!! Ultraspherical parameter
REAL(DFP), ALLOCATABLE :: ans(:, :)
!! interpolation points in xij format
END FUNCTION QuadraturePoint_Quadrangle4
END INTERFACE QuadraturePoint_Quadrangle