TABLE OF CONTENTS


ABINIT/integvol [ Functions ]

[ Top ] [ Functions ]

NAME

 integvol

FUNCTION

 This routine integrates the volume of the Bader atom

COPYRIGHT

 Copyright (C) 2002-2017 ABINIT group (PCasek,FF,XG)
 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

  (see side effects)

OUTPUT

  (see side effects)

SIDE EFFECTS

  This routine works on the data contained in the aimfields and aimprom modules

 WARNING
 This file does not follow the ABINIT coding rules (yet)

PARENTS

      drvaim

CHILDREN

      coeffs_gausslegint

SOURCE

 36 #if defined HAVE_CONFIG_H
 37 #include "config.h"
 38 #endif
 39 
 40 #include "abi_common.h"
 41 
 42 subroutine integvol()
 43 
 44  use defs_basis
 45  use defs_aimfields
 46  use defs_aimprom
 47  use defs_parameters
 48  use m_errors
 49  use m_profiling_abi
 50 
 51  use m_numeric_tools,   only : coeffs_gausslegint
 52 
 53 !This section has been created automatically by the script Abilint (TD).
 54 !Do not modify the following lines by hand.
 55 #undef ABI_FUNC
 56 #define ABI_FUNC 'integvol'
 57 !End of the abilint section
 58 
 59  implicit none
 60 
 61 !Arguments ------------------------------------
 62 
 63 !Local variables ------------------------------
 64 !scalars
 65  integer :: batom,ii,jj,nph,nth
 66  real(dp) :: chgint,ct1,ct2,nsphe,phimax,phimin
 67  real(dp) :: rsmax,rsmin,themax,themin
 68  logical :: gaus,weit
 69 !arrays
 70  real(dp) :: shift(3)
 71  real(dp),allocatable :: rdint(:,:)
 72  real(dp),allocatable :: wgrs(:,:)
 73 
 74 ! *********************************************************************
 75 
 76  tpi=two_pi
 77  gaus=.true.
 78  weit=.true.
 79 
 80 
 81  rewind(unts)
 82  read(unts,*) batom,shift
 83  read(unts,*) nth,themin,themax
 84  read(unts,*) nph,phimin,phimax
 85 
 86  write(std_out,*) 'NTH NPH ',nth,nph
 87 
 88  ABI_ALLOCATE(wgrs,(nth,nph))
 89  ABI_ALLOCATE(rdint,(nth,nph))
 90 
 91  do ii=1,nth
 92    do jj=1,nph
 93      if (weit) then
 94        read(unts,*) th(ii),ph(jj),rs(ii,jj),wgrs(ii,jj)
 95      else
 96        read(unts,*) th(ii),ph(jj),rs(ii,jj)
 97      end if
 98    end do
 99  end do
100  read(unts,*) rsmin,rsmax
101 
102 
103  if (gaus) then
104    ct1=cos(themin)
105    ct2=cos(themax)
106    call coeffs_gausslegint(ct1,ct2,cth,wcth,nth)
107    call coeffs_gausslegint(phimin,phimax,ph,wph,nph)
108  end if
109 
110  do ii=1,nth
111    do jj=1,nph
112      if (.not.weit) then
113        if (gaus) then
114          wgrs(ii,jj)=wcth(ii)*wph(jj)
115        else
116          wgrs(ii,jj)=1._dp
117        end if
118      end if
119    end do
120  end do
121 
122  nsphe=0._dp
123  do ii=1,nth
124    do jj=1,nph
125      nsphe=nsphe+rs(ii,jj)**3/3._dp*wgrs(ii,jj)
126    end do
127  end do
128  if (gaus.or.weit) then
129    nsphe=nsphe*(pi/(themin-themax))*(tpi/(phimax-phimin))
130  else
131    nsphe=nsphe/(nth*nph)*2.0*tpi
132  end if
133  chgint=nsphe
134 
135  write(std_out,*) ':VOLTOT ',batom,chgint
136  write(untout,'("Volume of the Bader atom: ", I6, F16.8)') batom,chgint
137 
138 end subroutine integvol