The Elk Code
spiralsc.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2012 S. Sharma, J. K. Dewhurst and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine spiralsc
7 use modmain
8 use modmpi
9 use moddelf
10 implicit none
11 ! local variables
12 integer nq,iq,jq
13 real(8) q
14 ! store original parameters
15 avec0(:,:)=avec(:,:)
16 natoms0(:)=natoms(:)
17 atposl0(:,:,:)=atposl(:,:,:)
18 bfcmt00(:,:,:)=bfcmt0(:,:,:)
19 mommtfix0(:,:,:)=mommtfix(:,:,:)
21 ngridk0(:)=ngridk
22 ! initialise universal variables
23 call init0
24 ! initialise q-point dependent variables
25 call init2
26 ! store original parameters
27 atposc0(:,:,:)=atposc(:,:,:)
28 10 continue
29 call sstask(80,filext)
30 ! if nothing more to do then restore input parameters and return
31 if (iqss == 0) then
32  filext='.OUT'
33  natoms(:)=natoms0(:)
34  avec(:,:)=avec0(:,:)
35  atposl(:,:,:)=atposl0(:,:,:)
36  bfcmt0(:,:,:)=bfcmt00(:,:,:)
37  mommtfix(:,:,:)=mommtfix0(:,:,:)
39  ngridk(:)=ngridk0(:)
40  return
41 end if
42 ! spiral dry run: just generate empty SS files
43 if (task == 352) goto 10
44 if (mp_mpi) then
45  write(*,'("Info(spiralsc): working on ",A)') 'SS'//trim(filext)
46 end if
47 ! determine k-point grid size from radkpt
48 autokpt=.true.
49 ! generate the spin-spiral supercell
50 call genscss
51 ! initialise or read the charge density and potentials from file
52 trdstate=(task == 351)
53 ! run the ground-state calculation
54 call gndstate
55 if (mp_mpi) then
56  write(80,'(I6,T20," : number of unit cells in supercell")') nscss
57  write(80,'(G18.10,T20," : total energy per unit cell")') engytot/dble(nscss)
58  write(80,*)
59  write(80,'("q-point in lattice and Cartesian coordinates :")')
60  write(80,'(3G18.10)') vql(:,iqss)
61  write(80,'(3G18.10)') vqc(:,iqss)
62  q=sqrt(vqc(1,iqss)**2+vqc(2,iqss)**2+vqc(3,iqss)**2)
63  write(80,'(G18.10,T20," : length of q-vector")') q
64  write(80,*)
65  nq=nint(dble(nqptnr)*wqpt(iqss))
66  write(80,'(I6,T20," : number of equivalent q-points")') nq
67  write(80,'("Equivalent q-points in lattice and Cartesian coordinates :")')
68  do iq=1,nqptnr
69  jq=ivqiq(ivq(1,iq),ivq(2,iq),ivq(3,iq))
70  if (jq == iqss) then
71  write(80,'(3G18.10)') vql(:,iq)
72  write(80,'(3G18.10)') vqc(:,iq)
73  write(80,*)
74  end if
75  end do
76  close(80)
77 end if
78 ! delete the eigenvector files
79 call delfiles(evec=.true.)
80 ! synchronise MPI processes
81 call mpi_barrier(mpicom,ierror)
82 goto 10
83 end subroutine
84 
real(8), dimension(3, maxatoms, maxspecies) bfcmt0
Definition: modmain.f90:275
subroutine gndstate
Definition: gndstate.f90:10
character(256) filext
Definition: modmain.f90:1301
real(8), dimension(3, maxatoms, maxspecies) atposl0
Definition: modmain.f90:51
integer task
Definition: modmain.f90:1299
logical mp_mpi
Definition: modmpi.f90:17
logical autokpt
Definition: modmain.f90:444
integer, dimension(:,:), allocatable ivq
Definition: modmain.f90:529
integer iqss
Definition: modmain.f90:297
subroutine spiralsc
Definition: spiralsc.f90:7
subroutine genscss
Definition: genscss.f90:7
integer, dimension(3) ngridk0
Definition: modmain.f90:448
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition: modmain.f90:51
real(8), dimension(:,:), allocatable vqc
Definition: modmain.f90:547
real(8) engytot
Definition: modmain.f90:983
subroutine sstask(fnum, fext)
Definition: sstask.f90:7
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
subroutine init2
Definition: init2.f90:7
integer, dimension(3) ngridk
Definition: modmain.f90:448
real(8), dimension(:,:), allocatable vql
Definition: modmain.f90:545
real(8), dimension(3, maxatoms, maxspecies) mommtfix0
Definition: modmain.f90:259
real(8), dimension(3, maxatoms, maxspecies) mommtfix
Definition: modmain.f90:259
integer, dimension(:,:,:), allocatable ivqiq
Definition: modmain.f90:531
real(8), dimension(3, 3) avec0
Definition: modmain.f90:12
integer nqptnr
Definition: modmain.f90:527
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
subroutine delfiles(evec, devec, eval, occ, pmat, epsi)
Definition: moddelf.f90:25
integer nscss
Definition: modmain.f90:299
Definition: modmpi.f90:6
real(8), dimension(3, maxatoms, maxspecies) bfcmt00
Definition: modmain.f90:275
logical trdstate
Definition: modmain.f90:682
subroutine init0
Definition: init0.f90:10
real(8), dimension(:), allocatable wqpt
Definition: modmain.f90:549
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition: modmain.f90:54
real(8), dimension(3, maxatoms, maxspecies) atposc0
Definition: modmain.f90:54
integer, dimension(maxspecies) natoms0
Definition: modmain.f90:36
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19
logical autokpt0
Definition: modmain.f90:444