The Elk Code
findngkmax.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2012 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: findngkmax
8 ! !INTERFACE:
9 pure subroutine findngkmax(nkpt,vkc,nspnfv,vqcss,ngv,vgc,gkmax,ngkmax)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! nkpt : number of k-points (in,integer)
12 ! vkc : k-point vectors in Cartesian coordinates (in,real(3,nkpt))
13 ! nspnfv : number of first-variational spin components: 1 normal case, 2 for
14 ! spin-spiral case (in,integer)
15 ! vqcss : spin-spiral q-vector, not referenced if nspnfv=1 (in,integer)
16 ! ngv : number of G-vectors (in,integer)
17 ! vgc : G-vectors in Cartesian coordinates (in,real(3,ngv))
18 ! gkmax : maximum allowed |G+k| (in,real)
19 ! ngkmax : maximum number of G+k-vectors over all k-points (out,integer)
20 ! !DESCRIPTION:
21 ! Determines the largest number of ${\bf G+k}$-vectors with length less than
22 ! {\tt gkmax} over all the $k$-points. This variable is used for allocating
23 ! arrays.
24 !
25 ! !REVISION HISTORY:
26 ! Created October 2004 (JKD)
27 ! Modified, August 2012 (JKD)
28 ! Removed modmain and added arguments, September 2012 (JKD)
29 !EOP
30 !BOC
31 implicit none
32 ! arguments
33 integer, intent(in) :: nkpt
34 real(8), intent(in) :: vkc(3,nkpt)
35 integer, intent(in) :: nspnfv
36 real(8), intent(in) :: vqcss(3)
37 integer, intent(in) :: ngv
38 real(8), intent(in) :: vgc(3,ngv),gkmax
39 integer, intent(out) :: ngkmax
40 ! local variables
41 integer ispn,ik,n,ig
42 real(8) v1,v2,v3,t0,t1
43 t0=gkmax**2+1.d-6
44 ngkmax=0
45 do ispn=1,nspnfv
46  do ik=1,nkpt
47  if (nspnfv == 2) then
48 ! spin-spiral case
49  if (ispn == 1) then
50  v1=vkc(1,ik)+0.5d0*vqcss(1)
51  v2=vkc(2,ik)+0.5d0*vqcss(2)
52  v3=vkc(3,ik)+0.5d0*vqcss(3)
53  else
54  v1=vkc(1,ik)-0.5d0*vqcss(1)
55  v2=vkc(2,ik)-0.5d0*vqcss(2)
56  v3=vkc(3,ik)-0.5d0*vqcss(3)
57  end if
58  else
59  v1=vkc(1,ik)
60  v2=vkc(2,ik)
61  v3=vkc(3,ik)
62  end if
63  n=0
64  do ig=1,ngv
65  t1=(vgc(1,ig)+v1)**2+(vgc(2,ig)+v2)**2+(vgc(3,ig)+v3)**2
66  if (t1 < t0) n=n+1
67  end do
68  if (n > ngkmax) ngkmax=n
69  end do
70 end do
71 end subroutine
72 !EOC
73 
pure subroutine findngkmax(nkpt, vkc, nspnfv, vqcss, ngv, vgc, gkmax, ngkmax)
Definition: findngkmax.f90:10