4
1
mirror of https://github.com/pfloos/quack synced 2024-07-04 18:36:03 +02:00
quack/src/eDFT/quadrature_grid.f90

80 lines
1.7 KiB
Fortran
Raw Normal View History

2019-03-13 11:07:31 +01:00
subroutine quadrature_grid(nRad,nAng,nGrid,root,weight)
! Build roots and weights of quadrature grid
implicit none
include 'parameters.h'
! Input variables
2020-03-23 17:26:12 +01:00
integer,intent(in) :: nRad
integer,intent(in) :: nAng
integer,intent(in) :: nGrid
2019-03-13 11:07:31 +01:00
! Local variables
integer :: i,j,k
double precision :: scale
double precision,allocatable :: Radius(:)
double precision,allocatable :: RadWeight(:)
double precision,allocatable :: XYZ(:,:)
double precision,allocatable :: XYZWeight(:)
! Output variables
2020-03-23 17:26:12 +01:00
double precision,intent(out) :: root(ncart,nGrid)
2019-03-13 11:07:31 +01:00
double precision,intent(out) :: weight(nGrid)
! Memory allocation
2020-03-23 17:26:12 +01:00
allocate(Radius(nRad),RadWeight(nRad),XYZ(ncart,nAng),XYZWeight(nAng))
2019-03-13 11:07:31 +01:00
! Findthe radial grid
scale = 1d0
call EulMac(Radius,RadWeight,nRad,scale)
write(*,20)
write(*,30)
write(*,20)
do i = 1,nRad
write(*,40) i,Radius(i),RadWeight(i)
end do
write(*,20)
write(*,*)
! Find the angular grid
call Lebdev(XYZ,XYZWeight,nAng)
write(*,20)
write(*,50)
write(*,20)
do j = 1,nAng
2020-03-23 17:26:12 +01:00
write(*,60) j,(XYZ(k,j),k=1,ncart), XYZWeight(j)
2019-03-13 11:07:31 +01:00
end do
write(*,20)
! Form the roots and weights
k = 0
do i=1,nRad
do j=1,nAng
k = k + 1
root(:,k) = Radius(i)*XYZ(:,j)
weight(k) = RadWeight(i)*XYZWeight(j)
enddo
enddo
2020-03-24 12:28:56 +01:00
! Print grids
2019-03-13 11:07:31 +01:00
20 format(T2,58('-'))
30 format(T20,'Radial Quadrature',/, &
T6,'I',T26,'Radius',T50,'Weight')
40 format(T3,I4,T18,F17.10,T35,F25.10)
50 format(T20,'Angular Quadrature',/, &
T6,'I',T19,'X',T29,'Y',T39,'Z',T54,'Weight')
60 format(T3,I4,T13,3F10.5,T50,F10.5)
end subroutine quadrature_grid