TABLE OF CONTENTS


ABINIT/init_e_field_vars [ Functions ]

[ Top ] [ Functions ]

NAME

 init_e_field_vars

FUNCTION

 Initialization of variables and data structures used in polarization
 calculations

COPYRIGHT

 Copyright (C) 2004-2018 ABINIT group 
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  dtset <type(dataset_type)> = all input variables in this dataset
  gmet(3,3) = reciprocal space metric tensor in bohr**-2
  gprimd(3,3) = primitive translations in recip space
  kg(3,mpw*mkmem) = reduced (integer) coordinates of G vecs in basis sphere
  mpi_enreg=information about MPI parallelization
  npwarr(nkpt) = number of planewaves in basis and boundary at this k point
  occ(mband*nkpt*nsppol) = occup number for each band at each k point
  pawang <type(pawang_type)>=paw angular mesh and related data
  pawrad(ntypat) <type(pawrad_type)>=paw radial mesh and related data
  pawtab(ntypat) <type(pawtab_type)>=paw tabulated starting data
  psps <type(pseudopotential_type)>=variables related to pseudopotentials
  rprimd(3,3) = dimensional primitive vectors
  symrec(3,3,nsym) = symmetries in reciprocal space in terms of
    reciprocal space primitive translations
  xred(3,natom) = location of atoms in reduced units

OUTPUT

  dtefield <type(efield_type)> :: initialized polarization variables
  pwind(pwind_alloc,2,3) = array used to compute the overlap matrix smat
                         between k-points k and k +- dk where dk is
                         parallel to the direction idir
  pwind_alloc = first dimension of pwind and pwnsfac
  pwnsfac(2,pwind_alloc) = phase factors for non-symmorphic translations

SIDE EFFECTS

 TO DO

NOTES

PARENTS

      gstate

CHILDREN

      initberry

SOURCE

 56 #if defined HAVE_CONFIG_H
 57 #include "config.h"
 58 #endif
 59 
 60 #include "abi_common.h"
 61 
 62 subroutine init_e_field_vars(dtefield,dtset,gmet,gprimd,kg,&
 63 &              mpi_enreg,npwarr,occ,pawang,pawrad,pawtab,psps,&
 64 &              pwind,pwind_alloc,pwnsfac,rprimd,symrec,xred)
 65 
 66  use defs_basis
 67  use defs_datatypes
 68  use defs_abitypes
 69  use m_efield
 70 
 71  use m_pawang, only : pawang_type
 72  use m_pawrad, only : pawrad_type
 73  use m_pawtab, only : pawtab_type
 74 
 75 !This section has been created automatically by the script Abilint (TD).
 76 !Do not modify the following lines by hand.
 77 #undef ABI_FUNC
 78 #define ABI_FUNC 'init_e_field_vars'
 79  use interfaces_67_common, except_this_one => init_e_field_vars
 80 !End of the abilint section
 81 
 82  implicit none
 83 
 84 !Arguments ------------------------------------
 85 !scalars
 86  integer,intent(out) :: pwind_alloc
 87  type(MPI_type),intent(inout) :: mpi_enreg
 88  type(dataset_type),intent(inout) :: dtset
 89  type(efield_type),intent(inout) :: dtefield !vz_i needs efield2
 90  type(pawang_type),intent(in) :: pawang
 91  type(pseudopotential_type),intent(in) :: psps
 92 !arrays
 93  integer,intent(in) :: kg(3,dtset%mpw*dtset%mkmem),npwarr(dtset%nkpt)
 94  integer,intent(in) :: symrec(3,3,dtset%nsym)
 95  integer,pointer :: pwind(:,:,:)
 96  real(dp),intent(in) :: gmet(3,3),gprimd(3,3),occ(dtset%mband*dtset%nkpt*dtset%nsppol)
 97  real(dp),intent(in) :: rprimd(3,3),xred(3,dtset%natom)
 98  real(dp),pointer :: pwnsfac(:,:)
 99  type(pawrad_type),intent(in) :: pawrad(dtset%ntypat*psps%usepaw)
100  type(pawtab_type),intent(in) :: pawtab(dtset%ntypat*psps%usepaw)
101 
102 !Local variables-------------------------------
103  logical :: initfield
104 !scalars
105 
106 ! *************************************************************************
107 
108  initfield = .false.
109 
110 !initialization
111  dtefield%has_qijb = 0
112  dtefield%has_epawf3 = 0
113  dtefield%has_epaws3 = 0
114  dtefield%has_expibi = 0
115  dtefield%has_rij = 0
116  dtefield%usecprj = 0
117  dtefield%berryopt = 0
118 
119  if ((dtset%berryopt < 0).or.(dtset%berryopt == 4) .or. (dtset%berryopt == 6) .or.(dtset%berryopt == 7) .or. &
120 & (dtset%berryopt == 14) .or.(dtset%berryopt == 16) .or.(dtset%berryopt == 17)) then 
121    nullify(pwind,pwnsfac)
122    call initberry(dtefield,dtset,gmet,gprimd,kg,&
123 &   dtset%mband,dtset%mkmem,mpi_enreg,dtset%mpw,&
124 &   dtset%natom,dtset%nkpt,npwarr,dtset%nsppol,&
125 &   dtset%nsym,dtset%ntypat,occ,pawang,pawrad,pawtab,&
126 &   psps,pwind,pwind_alloc,pwnsfac,rprimd,symrec,&
127 &   dtset%typat,psps%usepaw,xred)
128    initfield = .true.
129  end if
130 
131  if (.not. initfield .and. dtset%orbmag == 0) then
132     ! initorbmag.F90 also allocates pwind and pwnsfac
133    pwind_alloc = 1
134    ABI_ALLOCATE(pwind,(pwind_alloc,2,3))
135    ABI_ALLOCATE(pwnsfac,(2,pwind_alloc))
136    pwind(:,:,:)=0
137    pwnsfac(:,:)=zero
138  end if
139 
140 end subroutine init_e_field_vars