InterpolationPoint
This routine returns the interpolation points on triangle.
Interface
- ܀ Interface
- ️܀ Equidistance
- BlythPozChebyshev
- BlythPozLegendre
- ↢
INTERFACE
MODULE FUNCTION InterpolationPoint_Triangle(order, ipType, &
& layout, xij) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: order
!! order
INTEGER(I4B), INTENT(IN) :: ipType
!! interpolation point type
REAL(DFP), OPTIONAL, INTENT(IN) :: xij(:, :)
!! Coord of domain in xij format
CHARACTER(*), INTENT(IN) :: layout
!! local node numbering layout, always VEFC
REAL(DFP), ALLOCATABLE :: ans(:, :)
!! xij coordinates
END FUNCTION InterpolationPoint_Triangle
END INTERFACE
xij
xij
contains nodal coordinates of triangle inxij
format.SIZE(xij,1) = nsd
, andSIZE(xij,2)=3
- If
xij
is absent then unit triangle is assumed
ipType
-
ipType
is interpolation point type, it can take following values -
Equidistance
, uniformly/evenly distributed points -
BlythPozChebyshev
-
BlythPozLegendre
-
GaussLegendreLobatto
, which is same asIsaacLegendre
-
GaussChebyshevLobatto
, which is same asIsaacChebyshev
-
GaussJacobiLobatto
-
GaussUltrasphericalLobatto
-
IsaacChebyshev
-
IsaacLegendre
-
ChenBabuska
TODO -
Hesthaven
TODO -
Feket
TODO
layout
layout
specifies the arrangement of points. The nodes are always returned in “VEFC” format (vertex, edge, face, cell). 1:3 are vertex points, then edge, and then internal nodes. The internal nodes also follow the same convention. Please read Gmsh manual on this topic.
program main
use easifembase
implicit none
integer( i4b ) :: i1, i2, order, ipType
real( dfp ), allocatable :: x(:,:), xij(:, :)
CHARACTER( 20 ) :: layout
order=5
xij = RefTriangleCoord("UNIT")
ipType=Equidistance
layout="VEFC"
x =InterpolationPoint_Triangle( order=order, &
& ipType=ipType, layout=layout, xij=xij)
call display( Mdencode(TRANSPOSE(x)), "xij (order="//tostring(order)//")=" )
call blanklines(nol=2)
end program main
See results
xij (order=5)=
0 | 0 |
1 | 0 |
0 | 1 |
0.2 | 0 |
0.4 | 0 |
0.6 | 0 |
0.8 | 0 |
0.8 | 0.2 |
0.6 | 0.4 |
0.4 | 0.6 |
0.2 | 0.8 |
0 | 0.8 |
0 | 0.6 |
0 | 0.4 |
0 | 0.2 |
0.2 | 0.2 |
0.6 | 0.2 |
0.2 | 0.6 |
0.4 | 0.2 |
0.4 | 0.4 |
0.2 | 0.4 |
program main
use easifembase
implicit none
integer( i4b ) :: i1, i2, order, ipType
real( dfp ), allocatable :: x(:,:), xij(:, :)
CHARACTER( 20 ) :: layout
order=5
xij = RefTriangleCoord("UNIT")
ipType=BlythPozChebyshev
layout="VEFC"
x =InterpolationPoint_Triangle( order=order, &
& ipType=ipType, layout=layout, xij=xij)
call display( Mdencode(TRANSPOSE(x)), "xij (order="//tostring(order)//")=" )
call blanklines(nol=2)
end program main
See results
xij (order=5)=
0 | 0 |
1 | 0 |
0 | 1 |
9.54915E-02 | -2.02431E-14 |
0.34549 | -3.27886E-14 |
0.65451 | -3.28071E-14 |
0.90451 | -2.02431E-14 |
0.90451 | 9.54915E-02 |
0.65451 | 0.34549 |
0.34549 | 0.65451 |
9.54915E-02 | 0.90451 |
-2.02431E-14 | 0.90451 |
-3.28071E-14 | 0.65451 |
-3.27886E-14 | 0.34549 |
-2.02431E-14 | 9.54915E-02 |
0.14699 | 0.14699 |
0.70601 | 0.14699 |
0.14699 | 0.70601 |
0.41667 | 0.16667 |
0.41667 | 0.41667 |
0.16667 | 0.41667 |
program main
use easifembase
implicit none
integer( i4b ) :: i1, i2, order, ipType
real( dfp ), allocatable :: x(:,:), xij(:, :)
CHARACTER( 20 ) :: layout
order=5
xij = RefTriangleCoord("UNIT")
ipType=BlythPozLegendre
layout="VEFC"
x =InterpolationPoint_Triangle( order=order, &
& ipType=ipType, layout=layout, xij=xij)
call display( Mdencode(TRANSPOSE(x)), "xij (order="//tostring(order)//")=" )
call blanklines(nol=2)
end program main
See results
xij (order=5)=
-1.85037E-16 | -1.85037E-16 |
1 | -1.85037E-16 |
-1.85037E-16 | 1 |
0.11747 | -1.4803E-16 |
0.35738 | -3.70074E-17 |
0.64262 | -3.70074E-17 |
0.88253 | -1.29526E-16 |
0.88253 | 0.11747 |
0.64262 | 0.35738 |
0.35738 | 0.64262 |
0.11747 | 0.88253 |
-1.29526E-16 | 0.88253 |
-3.70074E-17 | 0.64262 |
-3.70074E-17 | 0.35738 |
-1.4803E-16 | 0.11747 |
0.15829 | 0.15829 |
0.68343 | 0.15829 |
0.15829 | 0.68343 |
0.4133 | 0.17339 |
0.4133 | 0.4133 |
0.17339 | 0.4133 |
- IsaacChebyshev
- IsaacLegendre
- ↢
program main
use easifembase
implicit none
integer( i4b ) :: i1, i2, order, ipType
real( dfp ), allocatable :: x(:,:), xij(:, :)
CHARACTER( 20 ) :: layout
order=5
xij = RefTriangleCoord("UNIT")
ipType=IsaacLegendre
layout="VEFC"
x =InterpolationPoint_Triangle( order=order, &
& ipType=ipType, layout=layout, xij=xij)
call display( Mdencode(TRANSPOSE(x)), "xij (order="//tostring(order)//")=" )
call blanklines(nol=2)
end program main
See results
xij (order=5)=
-8.32667E-17 | -8.32667E-17 |
1 | -8.32667E-17 |
-8.32667E-17 | 1 |
0.11747 | -4.89901E-17 |
0.35738 | -1.58335E-17 |
0.64262 | -1.58335E-17 |
0.88253 | -4.89901E-17 |
0.88253 | 0.11747 |
0.64262 | 0.35738 |
0.35738 | 0.64262 |
0.11747 | 0.88253 |
-4.89901E-17 | 0.88253 |
-1.58335E-17 | 0.64262 |
-1.58335E-17 | 0.35738 |
-4.89901E-17 | 0.11747 |
0.15599 | 0.15599 |
0.68802 | 0.15599 |
0.15599 | 0.68802 |
0.41807 | 0.16387 |
0.41807 | 0.41807 |
0.16387 | 0.41807 |
program main
use easifembase
implicit none
integer( i4b ) :: i1, i2, order, ipType
real( dfp ), allocatable :: x(:,:), xij(:, :)
CHARACTER( 20 ) :: layout
order=5
xij = RefTriangleCoord("UNIT")
ipType=IsaacLegendre
layout="VEFC"
x =InterpolationPoint_Triangle( order=order, &
& ipType=ipType, layout=layout, xij=xij)
call display( Mdencode(TRANSPOSE(x)), "xij (order="//tostring(order)//")=" )
call blanklines(nol=2)
end program main
See results
xij (order=5)=
-8.32667E-17 | -8.32667E-17 |
1 | -8.32667E-17 |
-8.32667E-17 | 1 |
0.11747 | -4.89901E-17 |
0.35738 | -1.58335E-17 |
0.64262 | -1.58335E-17 |
0.88253 | -4.89901E-17 |
0.88253 | 0.11747 |
0.64262 | 0.35738 |
0.35738 | 0.64262 |
0.11747 | 0.88253 |
-4.89901E-17 | 0.88253 |
-1.58335E-17 | 0.64262 |
-1.58335E-17 | 0.35738 |
-4.89901E-17 | 0.11747 |
0.15599 | 0.15599 |
0.68802 | 0.15599 |
0.15599 | 0.68802 |
0.41807 | 0.16387 |
0.41807 | 0.41807 |
0.16387 | 0.41807 |