TABLE OF CONTENTS


ABINIT/shellstruct [ Functions ]

[ Top ] [ Functions ]

NAME

  shellstruct 

 FUNCTION calculates shell structure (multiplicities, radii) of an atomic configuration 

INPUTS

  natom=number of atoms in unit cell
  xred=reduced coordinates of atoms
  rprimd=unit cell vectors
  magv = magnetic ordering of atoms given as 1 and -1, if not given fm is assumed
  atp = atom on which the perturbation was done

OUTPUT

  sdisv(nat)= distance of each shell to central atom (only the first nsh entries are relevant)
  nsh    = number of shells
  mult(nat) = number of atoms on shell      (only the first nsh entries are relevant)

SIDE EFFECTS

PARENTS

      pawuj_det

CHILDREN

      ioniondist,prmat,sort_dp,sort_int,wrtout

SOURCE

 30 #if defined HAVE_CONFIG_H
 31 #include "config.h"
 32 #endif
 33 
 34 #include "abi_common.h"
 35 
 36 
 37 subroutine shellstruct(xred,rprimd,natom,magv,distv,smult,sdisv,nsh,atp,prtvol) 
 38 
 39  use defs_basis
 40  use m_profiling_abi
 41  use m_sort
 42 
 43  use m_pptools,       only : prmat
 44 
 45 !This section has been created automatically by the script Abilint (TD).
 46 !Do not modify the following lines by hand.
 47 #undef ABI_FUNC
 48 #define ABI_FUNC 'shellstruct'
 49  use interfaces_14_hidewrite
 50  use interfaces_41_geometry, except_this_one => shellstruct
 51 !End of the abilint section
 52 
 53  implicit none
 54 
 55 !Arguments ------------------------------------
 56 !scalars
 57  integer,intent(in)              :: natom
 58  integer,intent(in),optional     :: atp  
 59  integer,intent(in),optional     :: prtvol
 60  integer,intent(out)             :: nsh          
 61 !arrays
 62  real(dp),intent(in)             :: rprimd(3,3)
 63  real(dp),intent(in)             :: xred(3,natom)
 64  integer,intent(out)             :: smult(natom) 
 65  integer,intent(in),optional     :: magv(natom)
 66  real(dp),intent(out)            :: sdisv(natom)
 67  real(dp),intent(out)            :: distv(natom)
 68 
 69 !Local variables-------------------------------
 70 !scalars
 71  integer                      :: iatom,atpp,ish,prtvoll
 72  character(len=500)           :: message
 73  real(dp),parameter           :: rndfact=10000_dp
 74 !arrays
 75  integer                      :: iperm(natom),jperm(natom)
 76  real(dp)                     :: distvh(natom,natom)
 77  real(dp)                     :: magvv(natom)
 78 
 79 ! *************************************************************************
 80 
 81  if (present(magv)) then
 82    magvv=magv
 83  else
 84    magvv=(/ (1, iatom=1,natom)  /) 
 85  end if
 86 
 87  if (present(atp)) then
 88    atpp=atp 
 89  else
 90    atpp=1 
 91  end if
 92 
 93  if (present(prtvol)) then
 94    prtvoll=prtvol 
 95  else
 96    prtvoll=1
 97  end if
 98 
 99 !DEBUB
100  write(std_out,*)'shellstruct start'
101 !END DEBUG 
102 
103 !Calculate ionic distances  
104  call ioniondist(natom,rprimd,xred,distvh,1,magv=int(magvv),atp=atpp)
105  distv=distvh(1,:)
106 
107  if (prtvol>2) then
108    write(std_out,'(a)')' shellstruct ionic distances in cell (distv) : '
109    call prmat(distv(1:natom),1,natom,1,std_out)
110  end if
111 
112  iperm=(/ (iatom, iatom=1,natom ) /)
113  jperm=iperm
114  distv=anint(distv*rndfact)/rndfact
115 !Sort distances
116  call sort_dp(natom,distv,iperm,10d-5)
117  call sort_int(natom,iperm,jperm)
118  
119  smult=0
120  sdisv=dot_product(rprimd(1,:),rprimd(1,:))+dot_product(rprimd(2,:),rprimd(2,:))+dot_product(rprimd(3,:),rprimd(3,:))
121 
122  nsh=1
123  smult(1)=1
124  sdisv(1)=distv(1)
125 
126  do iatom=2,natom
127    do ish=1,natom
128      if (distv(iatom)>sdisv(ish)) then
129        cycle  
130      else if (distv(iatom)==sdisv(ish)) then
131        smult(ish)=smult(ish)+1
132        exit    
133      else if (distv(iatom)<sdisv(ish)) then
134        smult(ish+1:natom)=smult(ish:natom-1)
135        sdisv(ish+1:natom)=sdisv(ish:natom-1)
136        smult(ish)=1
137        sdisv(ish)=distv(iatom)
138        nsh=nsh+1
139        exit
140      end if
141    end do  
142  end do
143 
144  distv=(/ ( distv(jperm(iatom)),iatom=1,natom ) /)
145  
146  if (prtvoll>2) then
147    write(message,'(a,i4,a)')' shellstruct found ',nsh,' shells at distances (sdisv) '
148    call wrtout(std_out,message,'COLL')
149    call prmat(sdisv(1:nsh),1,nsh,1,std_out)  
150    write(message,fmt='(a,150i4)')' and multiplicities (smult) ', smult(1:nsh)
151    call wrtout(std_out,message,'COLL')
152  end if
153  
154 !DEBUB
155  write(std_out,*)'shellstruct leave'
156 !END DEBUG
157 
158 end subroutine shellstruct