Initiate
This subroutine constructs the quadrature points.
Interface 1
INTERFACE Initiate
MODULE SUBROUTINE Initiate(obj, refElem, order, quadratureType, &
alpha, beta, lambda)
TYPE(QuadraturePoint_), INTENT(INOUT) :: obj
!! Total number of xidimension
CLASS(ReferenceElement_), INTENT(IN) :: refElem
!! Reference element
INTEGER(I4B), INTENT(IN) :: order
!! order of integrand
CHARACTER(*), INTENT(IN) :: quadratureType
!! Type of quadrature points
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda
!! Ultraspherical parameter
END SUBROUTINE Initiate
END INTERFACE Initiate
Type of quadrature points:
| Argument | Type | Intent | Optional | Description |
|---|---|---|---|---|
| obj | TYPE(QuadraturePoint_) | INOUT | No | The quadrature point object to be initialized |
| refElem | CLASS(ReferenceElement_) | IN | No | Reference element |
| order | INTEGER(I4B) | IN | No | Order of integrand |
| quadratureType | CHARACTER(*) | IN | No | Type of quadrature points |
| alpha | REAL(DFP) | IN | Yes | Jacobi parameter (required for Jacobi quadrature) |
| beta | REAL(DFP) | IN | Yes | Jacobi parameter (required for Jacobi quadrature) |
| lambda | REAL(DFP) | IN | Yes | Ultraspherical parameter (required for Ultraspherical quadrature) |
In the case of
Jacobialphaandbetashould be given, and in the case ofUltrasphericallambdashould be given. Details aboutquadratureTypecan be found in the QuadraturePoint documentation.
Example
!> author: Vikas Sharma, Ph. D.
! date: 2025-06-17
! summary: This program tests the Initiate method of the QuadraturePoint_ class.
! It runs five tests with different quadrature types and reference line domains:
! GaussLegendre (BIUNIT, UNIT), GaussLegendreLobatto (BIUNIT),
! GaussLegendreRadauLeft (BIUNIT), and GaussLegendreRadauRight (BIUNIT).
! Only Line elements are tested in this example.
PROGRAM main
USE GlobalData
USE QuadraturePoint_Method
USE BaseType
USE Display_Method
USE ReferenceLine_Method
IMPLICIT NONE
TYPE(QuadraturePoint_) :: obj
TYPE(ReferenceLine_) :: refelem
INTEGER(I4B) :: order, nip(1), quadratureType
CHARACTER(:), ALLOCATABLE :: test_name, domainName, quadratureTypeChar
CALL test1
CALL test2
CONTAINS
!----------------------------------------------------------------------------
! test1
!----------------------------------------------------------------------------
! This subroutine tests the Initiate method with GaussLegendre quadrature
! on a reference line in the BIUNIT domain.
! It initializes the QuadraturePoint_ object and displays the results.
SUBROUTINE test1
domainName = "BIUNIT"
refelem = ReferenceLine(nsd=1_I4B, xij=RefCoord_Line(domainName))
order = 4_I4B
quadratureType = GaussLegendre
quadratureTypeChar = "GaussLegendre"
nip(1) = 3
test_name = "order="//tostring(order)//" quadratureType=GaussLegendre"// &
" refelemDomain="//domainName
CALL Initiate(obj=obj, refelem=refelem, order=order, &
quadratureType=quadratureType)
CALL BlankLines()
CALL Display(obj, test_name)
CALL BlankLines()
! In the above interface we are creating quadrature points of
! given accuracy. We can also specify the number of quadrature point
! as show below
CALL Initiate(obj=obj, refelem=refelem, nips=nip, &
quadratureType=quadratureType)
! We can also specify quadratureType as string
CALL Initiate(obj=obj, refelem=refelem, order=order, &
quadratureType=quadratureTypeChar)
CALL Initiate(obj=obj, refelem=refelem, nips=nip, &
quadratureType=quadratureTypeChar)
END SUBROUTINE test1
!----------------------------------------------------------------------------
!
!----------------------------------------------------------------------------
SUBROUTINE test2
domainName = "UNIT"
refelem = ReferenceLine(nsd=1_I4B, xij=RefCoord_Line(domainName))
order = 4_I4B
quadratureType = GaussLegendre
test_name = "order="//tostring(order)//" quadratureType=GaussLegendre"// &
" refelemDomain="//domainName
CALL Initiate(obj=obj, refelem=refelem, order=order, &
quadratureType=quadratureType)
CALL BlankLines()
CALL Display(obj, test_name)
CALL BlankLines()
END SUBROUTINE test2
!----------------------------------------------------------------------------
! test3
!----------------------------------------------------------------------------
! This subroutine tests the Initiate method with GaussLegendre quadrature
! on a reference line in the UNIT domain.
! It initializes the QuadraturePoint_ object and displays the results.
SUBROUTINE test3
domainName = "BIUNIT"
refelem = ReferenceLine(nsd=1_I4B, xij=RefCoord_Line(domainName))
order = 4_I4B
quadratureType = GaussLegendreLobatto
test_name = "order="//tostring(order)//" quadratureType=GaussLegendreLobatto"// &
" refelemDomain="//domainName
CALL Initiate(obj=obj, refelem=refelem, order=order, &
quadratureType=quadratureType)
CALL BlankLines()
CALL Display(obj, test_name)
CALL BlankLines()
END SUBROUTINE test3
!----------------------------------------------------------------------------
! test4
!----------------------------------------------------------------------------
! This subroutine tests the Initiate method with GaussLegendreRadauLeft quadrature
! on a reference line in the BIUNIT domain.
! It initializes the QuadraturePoint_ object and displays the results.
SUBROUTINE test4
domainName = "BIUNIT"
refelem = ReferenceLine(nsd=1_I4B, xij=RefCoord_Line(domainName))
order = 4_I4B
quadratureType = GaussLegendreRadauLeft
test_name = "order="//tostring(order)//" quadratureType=GaussLegendreRadauLeft"// &
" refelemDomain="//domainName
CALL Initiate(obj=obj, refelem=refelem, order=order, &
quadratureType=quadratureType)
CALL BlankLines()
CALL Display(obj, test_name)
CALL BlankLines()
END SUBROUTINE test4
!----------------------------------------------------------------------------
! test5
!----------------------------------------------------------------------------
! This subroutine tests the Initiate method with GaussLegendreRadauRight quadrature
! on a reference line in the BIUNIT domain.
! It initializes the QuadraturePoint_ object and displays the results.
SUBROUTINE test5
domainName = "BIUNIT"
refelem = ReferenceLine(nsd=1_I4B, xij=RefCoord_Line(domainName))
order = 4_I4B
quadratureType = GaussLegendreRadauRight
test_name = "order="//tostring(order)//" quadratureType=GaussLegendreRadauRight"// &
" refelemDomain="//domainName
CALL Initiate(obj=obj, refelem=refelem, order=order, &
quadratureType=quadratureType)
CALL BlankLines()
CALL Display(obj, test_name)
CALL BlankLines()
END SUBROUTINE test5
!----------------------------------------------------------------------------
!
!----------------------------------------------------------------------------
END PROGRAM main
Interface 2
INTERFACE
MODULE SUBROUTINE Initiate(obj, refElem, nips, quadratureType, &
alpha, beta, lambda)
TYPE(QuadraturePoint_), INTENT(INOUT) :: obj
!! Total number of xidimension
CLASS(ReferenceElement_), INTENT(IN) :: refElem
!! Reference element
INTEGER(I4B), INTENT(IN) :: nips(1)
!! order of integrand
CHARACTER(*), INTENT(IN) :: quadratureType
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda
!! Ultraspherical parameter
END SUBROUTINE Initiate
END INTERFACE
This method allows us to initiate the quadrature points by specifying the number of integration points, and quadrature type.
| Argument | Type | Intent | Optional | Description |
|---|---|---|---|---|
| obj | TYPE(QuadraturePoint_) | INOUT | No | The quadrature point object to be initialized |
| refElem | CLASS(ReferenceElement_) | IN | No | Reference element |
| nips | INTEGER(I4B) | IN | No | Number of integration points (array of size 1) |
| quadratureType | CHARACTER(*) | IN | No | Type of quadrature points |
| alpha | REAL(DFP) | IN | Yes | Jacobi parameter (required for Jacobi quadrature) |
| beta | REAL(DFP) | IN | Yes | Jacobi parameter (required for Jacobi quadrature) |
| lambda | REAL(DFP) | IN | Yes | Ultraspherical parameter (required for Ultraspherical quadrature) |
Example
See the example given in the Interface 1
Interface 3
This interface is similar to the Interface 1, but in this interface quadratureType are in integer. You can read more about the quadrature types in the QuadraturePoint documentation.
INTERFACE Initiate
MODULE SUBROUTINE Initiate(obj, refElem, order, quadratureType, alpha, beta, lambda)
TYPE(QuadraturePoint_), INTENT(INOUT) :: obj
!! Total number of xidimension
CLASS(ReferenceElement_), INTENT(IN) :: refElem
!! Reference-element
INTEGER(I4B), INTENT(IN) :: order
!! order of integrand
INTEGER(I4B), INTENT(IN) :: quadratureType
!! Type of quadrature points
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda
!! Ultraspherical parameter
END SUBROUTINE Initiate
END INTERFACE Initiate
| Argument | Type | Intent | Optional | Description |
|---|---|---|---|---|
| obj | TYPE(QuadraturePoint_) | INOUT | No | The quadrature point object to be initialized |
| refElem | CLASS(ReferenceElement_) | IN | No | Reference element |
| order | INTEGER(I4B) | IN | No | Order of integrand |
| quadratureType | INTEGER(I4B) | IN | No | Type of quadrature points (integer enumeration) |
| alpha | REAL(DFP) | IN | Yes | Jacobi parameter (required for Jacobi quadrature) |
| beta | REAL(DFP) | IN | Yes | Jacobi parameter (required for Jacobi quadrature) |
| lambda | REAL(DFP) | IN | Yes | Ultraspherical parameter (required for Ultraspherical quadrature) |
Example
- See the example given in the interface 1
Interface 4
- This inteface is similar to the Interface 3, but in this interface the number of integration points is given.
- This interface is simular to the Interface 2, but in this interface quadratureType are in integer.
INTERFACE Initiate
MODULE SUBROUTINE Initiate(obj, refElem, nips, quadratureType, alpha, beta, lambda)
TYPE(QuadraturePoint_), INTENT(INOUT) :: obj
!! Total number of xidimension
CLASS(ReferenceElement_), INTENT(IN) :: refElem
!! Reference element
INTEGER(I4B), INTENT(IN) :: nips(1)
!! order of integrand
INTEGER(I4B), INTENT(IN) :: quadratureType
!! Type of quadrature points
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda
!! Ultraspherical parameter
END SUBROUTINE Initiate
END INTERFACE Initiate
| Argument | Type | Intent | Optional | Description |
|---|---|---|---|---|
| obj | TYPE(QuadraturePoint_) | INOUT | No | The quadrature point object to be initialized |
| refElem | CLASS(ReferenceElement_) | IN | No | Reference element |
| nips | INTEGER(I4B) | IN | No | Number of integration points (array of size 1) |
| quadratureType | INTEGER(I4B) | IN | No | Type of quadrature points (integer enumeration) |
| alpha | REAL(DFP) | IN | Yes | Jacobi parameter (required for Jacobi quadrature) |
| beta | REAL(DFP) | IN | Yes | Jacobi parameter (required for Jacobi quadrature) |
| lambda | REAL(DFP) | IN | Yes | Ultraspherical parameter (required for Ultraspherical quadrature) |
Interface 5
- This inteface is mainly for cartesian elements. In this interface you can specify order of integrand in each direction.
- The method returns the quadrature points which are obtained by the outer product of the quadrature points in each direction.
INTERFACE Initiate
MODULE SUBROUTINE quad_initiate7(obj, refElem, p, q, r, quadratureType1, quadratureType2, quadratureType3, &
alpha1, beta1, lambda1, alpha2, beta2, lambda2, alpha3, beta3, lambda3)
TYPE(QuadraturePoint_), INTENT(INOUT) :: obj
!! Total number of xidimension
CLASS(ReferenceElement_), INTENT(IN) :: refElem
!! Reference-element
INTEGER(I4B), INTENT(IN) :: p
!! order of integrand in x
INTEGER(I4B), INTENT(IN) :: q
!! order of integrand in y
INTEGER(I4B), INTENT(IN) :: r
!! order of integrand in z direction
INTEGER(I4B), INTENT(IN) :: quadratureType1
!! Type of quadrature points
INTEGER(I4B), INTENT(IN) :: quadratureType2
!! Type of quadrature points
INTEGER(I4B), INTENT(IN) :: quadratureType3
!! Type of quadrature points
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha1, beta1, lambda1
!! Jacobi parameter and Ultraspherical parameters
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha2, beta2, lambda2
!! Jacobi parameter and Ultraspherical parameters
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha3, beta3, lambda3
!! Jacobi parameter and Ultraspherical parameters
END SUBROUTINE quad_initiate7
END INTERFACE Initiate
| Argument | Type | Intent | Optional | Description |
|---|---|---|---|---|
| obj | TYPE(QuadraturePoint_) | INOUT | No | The quadrature point object to be initialized |
| refElem | CLASS(ReferenceElement_) | IN | No | Reference element |
| p | INTEGER(I4B) | IN | No | Order of integrand in x direction |
| q | INTEGER(I4B) | IN | No | Order of integrand in y direction |
| r | INTEGER(I4B) | IN | No | Order of integrand in z direction |
| quadratureType1 | INTEGER(I4B) | IN | No | Type of quadrature points for x direction |
| quadratureType2 | INTEGER(I4B) | IN | No | Type of quadrature points for y direction |
| quadratureType3 | INTEGER(I4B) | IN | No | Type of quadrature points for z direction |
| alpha1, beta1, lambda1 | REAL(DFP) | IN | Yes | Parameters for x direction quadrature |
| alpha2, beta2, lambda2 | REAL(DFP) | IN | Yes | Parameters for y direction quadrature |
| alpha3, beta3, lambda3 | REAL(DFP) | IN | Yes | Parameters for z direction quadrature |
Interface 6
This interface is similar to the Interface 5, but in this interface the number of integration points is given.
INTERFACE Initiate
MODULE SUBROUTINE quad_initiate8(obj, refElem, nipsx, nipsy, nipsz, quadratureType1, quadratureType2, quadratureType3, &
alpha1, beta1, lambda1, alpha2, beta2, lambda2, alpha3, beta3, lambda3)
TYPE(QuadraturePoint_), INTENT(INOUT) :: obj
!! Total number of xidimension
CLASS(ReferenceElement_), INTENT(IN) :: refElem
!! Reference element
INTEGER(I4B), INTENT(IN) :: nipsx(1)
!! number of integration points in x direction
INTEGER(I4B), INTENT(IN) :: nipsy(1)
!! number of integration points in y direction
INTEGER(I4B), INTENT(IN) :: nipsz(1)
!! number of integration points in z direction
INTEGER(I4B), INTENT(IN) :: quadratureType1
!! Type of quadrature points
INTEGER(I4B), INTENT(IN) :: quadratureType2
!! Type of quadrature points
INTEGER(I4B), INTENT(IN) :: quadratureType3
!! Type of quadrature points
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha1, beta1, lambda1
!! Jacobi parameter and Ultraspherical parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha2, beta2, lambda2
!! Jacobi parameter and Ultraspherical parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha3, beta3, lambda3
!! Jacobi parameter and Ultraspherical parameter
END SUBROUTINE quad_initiate8
END INTERFACE Initiate
| Argument | Type | Intent | Optional | Description |
|---|---|---|---|---|
| obj | TYPE(QuadraturePoint_) | INOUT | No | The quadrature point object to be initialized |
| refElem | CLASS(ReferenceElement_) | IN | No | Reference element |
| nipsx | INTEGER(I4B) | IN | No | Number of integration points in x direction (array of size 1) |
| nipsy | INTEGER(I4B) | IN | No | Number of integration points in y direction (array of size 1) |
| nipsz | INTEGER(I4B) | IN | No | Number of integration points in z direction (array of size 1) |
| quadratureType1 | INTEGER(I4B) | IN | No | Type of quadrature points for x direction |
| quadratureType2 | INTEGER(I4B) | IN | No | Type of quadrature points for y direction |
| quadratureType3 | INTEGER(I4B) | IN | No | Type of quadrature points for z direction |
| alpha1, beta1, lambda1 | REAL(DFP) | IN | Yes | Parameters for x direction quadrature |
| alpha2, beta2, lambda2 | REAL(DFP) | IN | Yes | Parameters for y direction quadrature |
| alpha3, beta3, lambda3 | REAL(DFP) | IN | Yes | Parameters for z direction quadrature |
Interface 7
This routine is for developers only, normal user should not worry about this interface.
MODULE PURE SUBROUTINE Initiate(obj, Points)
CLASS(QuadraturePoint_), INTENT(INOUT) :: obj
REAL(DFP), INTENT(IN) :: Points(:, :)
!! Points contains the quadrature points and weights
!! Points( :, ipoint ) contains quadrature points and weights of ipoint
!! quadrature point. The last row contains the weight. The rest of the
!! rows contains the coordinates of quadrature.
END SUBROUTINE Initiate
| Argument | Type | Intent | Optional | Description |
|---|---|---|---|---|
| obj | CLASS(QuadraturePoint_) | INOUT | No | The quadrature point object to be initialized |
| Points | REAL(DFP) | IN | No | Matrix containing quadrature points and weights. Last row contains weights, other rows contain coordinates. |
This subroutine converts the quadrature points into QuadraturePoint_ data type.
This routine is for developers only, normal user should not worry about this interface.
PointsPointscontains the quadrature points and weights- for example,
Points( :, ipoint )contains quadrature points and weights ofipointquadrature point. - The last row always contains the weight.
- The rest of the rows contains the coordinates of quadrature.
Interface 8
This routine is for developers only, normal user should not worry about this interface.
MODULE PURE SUBROUTINE Initiate(obj, tXi, tPoints)
CLASS(QuadraturePoint_), INTENT(INOUT) :: obj
INTEGER(I4B), INTENT(IN) :: tXi
!! Total number of xidimension
!! For line tXi=1
!! For 2D element tXi=2
!! For 3D element tXi=3
INTEGER(I4B), INTENT(IN) :: tPoints
!! Total number quadrature points
END SUBROUTINE Initiate
| Argument | Type | Intent | Optional | Description |
|---|---|---|---|---|
| obj | CLASS(QuadraturePoint_) | INOUT | No | The quadrature point object to be initialized |
| tXi | INTEGER(I4B) | IN | No | Total number of dimensions (1 for line, 2 for 2D element, 3 for 3D element) |
| tPoints | INTEGER(I4B) | IN | No | Total number of quadrature points |
The following interface allocates the memory for storing the quadrature points.
QuadraturePoint (Constructor function)
MODULE PURE FUNCTION QuadraturePoint(Points) RESULT(obj)
TYPE(QuadraturePoint_) :: obj
REAL(DFP), INTENT(IN) :: Points(:, :)
END FUNCTION QuadraturePoint
This function converts Points to an instance of QuadraturePoint.
QuadraturePoint_Pointer
MODULE PURE FUNCTION QuadraturePoint_Pointer(Points) RESULT(obj)
CLASS(QuadraturePoint_), POINTER :: obj
REAL(DFP), INTENT(IN) :: Points(:, :)
END FUNCTION QuadraturePoint_Pointer