pwdft  0.1
PW-DFT code for education
gvec Module Reference

Data Types

type  fft
 

Functions/Subroutines

subroutine setup_gvec (xfft, g2cut)
 
subroutine base2357 (nfft)
 

Variables

type(fft), save g_rh
 
type(fft), save g_wf
 

Function/Subroutine Documentation

◆ setup_gvec()

subroutine gvec::setup_gvec ( type(fft), intent(out)  xfft,
real(8), intent(in)  g2cut 
)

Definition at line 43 of file gvec.F90.

43  !
44  use constant, only : pi
45  use atm_spec, only : avec, bvec
46  !
47  real(8),intent(in) :: g2cut
48  type(fft),intent(out) :: xfft
49  !
50  integer :: ii, i1, i2, i3, nmax3(3), nmax
51  real(8) :: norm, gv0(3), g2
52  integer,allocatable :: mill0(:,:)
53  !
54  do ii = 1, 3
55  norm = sqrt(dot_product(avec(1:3,ii), avec(1:3,ii)))
56  nmax3(ii) = ceiling(sqrt(g2cut)*norm/(2.0d0*pi))
57  end do
58  !
59  nmax = product(nmax3(1:3)*2+1)
60  allocate(mill0(3,nmax))
61  !
62  xfft%npw = 0
63  do i3 = -nmax3(3), nmax3(3)
64  do i2 = -nmax3(2), nmax3(2)
65  do i1 = -nmax3(1), nmax3(1)
66  gv0(1:3) = matmul(bvec(1:3,1:3), dble((/i1,i2,i3/)))
67  g2 = dot_product(gv0,gv0)
68  if(g2 < g2cut) then
69  xfft%npw = xfft%npw + 1
70  mill0(1:3,xfft%npw) = (/i1,i2,i3/)
71  end if
72  end do
73  end do
74  end do
75  !
76  write(*,*) " Numver of PW inside sphere : ", xfft%npw
77  !
78  ! Find FFT grid
79  !
80  do ii = 1, 3
81  xfft%nft(ii) = maxval(mill0(ii,1:xfft%npw))*2 + 1
82  call base2357(xfft%nft(ii))
83  end do
84  !
85  write(*,*) " FFT grid : ", xfft%nft(1:3)
86  xfft%nr = product(xfft%nft(1:3))
87  write(*,*) " Numver of PW in grid : ", xfft%nr
88  allocate(xfft%mill(3,xfft%nr), xfft%map(xfft%npw))
89  !
90  ! Whole grid : Miller index
91  !
92  ii = 0
93  do i3 = 0, xfft%nft(3) - 1
94  do i2 = 0, xfft%nft(2) - 1
95  do i1 = 0, xfft%nft(1) - 1
96  ii = ii + 1
97  xfft%mill(1:3,ii) = modulo((/i1, i2, i3/) + xfft%nft(1:3)/2, xfft%nft(1:3)) &
98  & - xfft%nft(1:3)/2
99  end do
100  end do
101  end do
102  !
103  ! Mapper
104  !
105  do ii = 1, xfft%npw
106  mill0(1:3,ii) = modulo(mill0(1:3,ii), xfft%nft(1:3))
107  xfft%map(ii) = 1 + mill0(1,ii) + xfft%nft(1)*(mill0(2,ii) + xfft%nft(2)*mill0(3,ii))
108  end do
109  !
110  deallocate(mill0)
111  !

References atm_spec::avec, base2357(), atm_spec::bvec, and constant::pi.

Referenced by stdin::read_stdin().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ base2357()

subroutine gvec::base2357 ( integer, intent(inout)  nfft)

Definition at line 115 of file gvec.F90.

115  !
116  integer,intent(inout) :: nfft
117  !
118  integer :: base(4) = (/2,3,5,7/), ibase, iexp, &
119  & start, last, nexp, nfft1, nfft2
120  !
121  start = nfft
122  last = 2 * nfft
123  !
124  do nfft1 = start, last
125  !
126  nfft2 = nfft1
127  !
128  do ibase = 1, 4
129  !
130  nexp = ceiling(log(dble(nfft2))/log(dble(base(ibase))))
131  !
132  do iexp = 1, nexp
133  if(mod(nfft2, base(ibase)) == 0) then
134  nfft2 = nfft2 / base(ibase)
135  else
136  exit
137  end if
138  end do
139  !
140  if(nfft2 == 1) then
141  nfft = nfft1
142  return
143  end if
144  !
145  end do
146  !
147  end do
148  !

Referenced by setup_gvec().

Here is the caller graph for this function:

Variable Documentation

◆ g_rh

◆ g_wf

atm_spec::avec
real(8), dimension(3, 3), save avec
Unit lattice vector [Bohr].
Definition: atm_spec.F90:50
constant
Definition: constant.F90:23
atm_spec::bvec
real(8), dimension(3, 3), save bvec
Unit reciplocal lattice vector [Bohr^-1].
Definition: atm_spec.F90:50
atm_spec
Definition: atm_spec.F90:23