better test automatically

This commit is contained in:
Martin Diehl 2021-07-21 18:06:37 +02:00
parent b1ba64e6af
commit 000007ba59
1 changed files with 20 additions and 1 deletions

View File

@ -39,7 +39,7 @@ contains
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief initializes FEM interpolation data !> @brief initializes FEM interpolation data
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine FEM_quadrature_init subroutine FEM_quadrature_init()
print'(/,a)', ' <<<+- FEM_quadrature init -+>>>'; flush(6) print'(/,a)', ' <<<+- FEM_quadrature init -+>>>'; flush(6)
@ -182,6 +182,8 @@ subroutine FEM_quadrature_init
FEM_quadrature_points (3,5)%p(121:156)= permutationStar211([0.3523052600879940_pReal, 0.0992057202494530_pReal]) FEM_quadrature_points (3,5)%p(121:156)= permutationStar211([0.3523052600879940_pReal, 0.0992057202494530_pReal])
FEM_quadrature_points (3,5)%p(157:168)= permutationStar31([0.1344783347929940_pReal]) FEM_quadrature_points (3,5)%p(157:168)= permutationStar31([0.1344783347929940_pReal])
call selfTest
end subroutine FEM_quadrature_init end subroutine FEM_quadrature_init
@ -369,4 +371,21 @@ pure function permutationStar1111(point) result(qPt)
end function permutationStar1111 end function permutationStar1111
!--------------------------------------------------------------------------------------------------
!> @brief Check correctness of quadrature weights.
!--------------------------------------------------------------------------------------------------
subroutine selfTest
integer :: o, d
do d = lbound(FEM_quadrature_weights,1), ubound(FEM_quadrature_weights,1)
do o = lbound(FEM_quadrature_weights(d,:),1), ubound(FEM_quadrature_weights(d,:),1)
if (dNeq(sum(FEM_quadrature_weights(d,o)%p),1.0_pReal,5e-15_pReal)) error stop 'quadrature weights'
enddo
enddo
end subroutine selfTest
end module FEM_quadrature end module FEM_quadrature