wann_optional.f 4.26 KB
Newer Older
1 2 3 4 5 6
!--------------------------------------------------------------------------------
! Copyright (c) 2016 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
! This file is part of FLEUR and available as free software under the conditions
! of the MIT license as expressed in the LICENSE file in more detail.
!--------------------------------------------------------------------------------

7 8 9 10
      module m_wann_optional
      use m_juDFT
      contains
      subroutine wann_optional(
11
     >               input,atoms,sym,cell,oneD,noco,wann,
12
     >               l_ms,l_sgwf,l_socgwf,
13
     >               aux_latt_const,param_file,l_dim)
14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30
c**************************************************
c     Make preparations for the calculation of 
c     Wannier functions.
c     Frank Freimuth
c**************************************************
      use m_types
      use m_wann_read_inp
      use m_wann_projgen
      use m_wann_kpointgen
      use m_wann_w90kpointgen
      use m_wann_kptsreduc
      use m_wann_kptsreduc2
      use m_wann_wan90prep
      use m_wann_dipole3
      use m_wann_dipole

      implicit none
31

32 33 34 35 36 37 38
      TYPE(t_input),     INTENT(IN)    :: input
      TYPE(t_atoms),     INTENT(IN)    :: atoms
      TYPE(t_sym),       INTENT(IN)    :: sym
      TYPE(t_cell),      INTENT(IN)    :: cell
      TYPE(t_oneD),      INTENT(IN)    :: oneD
      TYPE(t_noco),      INTENT(IN)    :: noco
      TYPE(t_wann),      INTENT(INOUT) :: wann
39 40

      real,intent(in)     :: aux_latt_const
41
      character(len=20),intent(in) :: param_file
42 43 44 45

      integer       :: num_wann(2)
      logical       :: l_nocosoc

46 47 48
      logical, intent(in) :: l_ms,l_sgwf,l_socgwf
      logical,intent(in) :: l_dim(3)

49
      l_nocosoc=noco%l_noco.or.noco%l_soc
50 51 52 53 54 55 56 57

c-----read the input file to determine what to do
      call wann_read_inp(
     >         .true.,wann)

c-----generate projection-definition-file
      if(wann%l_projgen) then
         call wann_projgen(
58
     >      atoms%ntype,atoms%neq,atoms%nat,atoms%zatom,l_nocosoc,wann)
59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74
         wann%l_stopopt=.true.
      endif

c-----generate k-point-files
      if(wann%l_kpointgen) then
             call wann_kpointgen()
             wann%l_stopopt=.true.
      endif
      if(wann%l_w90kpointgen) then
             call wann_w90kpointgen()
             wann%l_stopopt=.true.
      endif

c-----find Wannier-irreducible part of BZ
      if(wann%l_kptsreduc)then
         call wann_kptsreduc(
75 76
     >            sym%nop,sym%mrot,cell%bmat,sym%tau,input%film,
     >            oneD%odi%d1,l_nocosoc)
77 78 79 80 81 82 83
         wann%l_stopopt=.true.
      endif

c-----find Wannier-irreducible part of BZ
      if(wann%l_kptsreduc2)then
         call wann_kptsreduc2(
     >            wann%mhp,
84 85
     >            sym%nop,sym%mrot,cell%bmat,sym%tau,input%film,
     >            oneD%odi%d1,l_nocosoc)
86 87 88 89 90 91
         wann%l_stopopt=.true.
      endif

c-----generate WF1.win and bkpts
      if(wann%l_prepwan90)then
         call wann_wan90prep(
92 93 94 95
     >            input%jspins,cell%amat,cell%bmat,
     >            atoms%nat,atoms%taual,atoms%zatom,atoms%ntype,
     >            atoms%ntype,atoms%neq,wann%l_bzsym,input%film,
     >            oneD%odi%d1,l_ms,l_sgwf,l_socgwf,aux_latt_const,
96
     >            param_file,l_dim)
97 98 99 100 101 102 103 104
      endif

c-----calculate polarization, if not wannierize
c-----if wannierize, then calculate polarization later (after wannierize)
      if(wann%l_dipole3.and..not.wann%l_wannierize)then
         num_wann(1)=wann%band_max(1)-wann%band_min(1)+1
         num_wann(2)=wann%band_max(2)-wann%band_min(2)+1
         call wann_dipole3(
105 106 107
     >            input%jspins,cell%omtil,atoms%nat,atoms%pos,
     >            cell%amat,cell%bmat,atoms%taual,num_wann,
     >            atoms%ntype,atoms%neq,atoms%zatom,l_nocosoc)
108 109 110 111 112 113 114
         wann%l_stopopt=.true.
      endif

c-----calculate polarization, if not wannierize
c-----if wannierize, then calculate polarization later (after wannierize)
      if(wann%l_dipole.and..not.wann%l_wannierize)then
         call wann_dipole(
115 116
     >            input%jspins,cell%omtil,atoms%nat,atoms%pos,
     >            cell%amat,atoms%ntype,atoms%neq,atoms%zatom)
117 118 119
         wann%l_stopopt=.true.
      endif

120
      IF(wann%l_stopopt)  CALL juDFT_end("wann_optional done",1) ! The 1 is temporarily. Should be mpi%irank.
121 122 123

      end subroutine wann_optional
      end module m_wann_optional