TABLE OF CONTENTS


ABINIT/ioniondist [ Functions ]

[ Top ] [ Functions ]

NAME

 ioniondist

FUNCTION

  Compute ion-ion distances

INPUTS

  natom= number of atoms in unit cell
  rprimd(3,3)=dimensional real space primitive translations (bohr)
  xred(3,natom)=dimensionless reduced coordinates of atoms
  inm(natom,natom)=index (m,n) of the atom
  option= 1 output ion-ion distances / 2 output ordering of ion-ion
          distances / 3 output variables in varlist
          according to ion-ion distances * magnetic ordering
          magv magnetic ordering of atoms given als 1 and -1, if not
          given fm is assumed
  varlist=List of variables
  magv(natom)= magnetic ordering of atoms
  atp=atom on which the perturbation was done

OUTPUT

SIDE EFFECTS

PARENTS

      pawuj_utils,shellstruct

CHILDREN

      prmat,wrtout

SOURCE

 35 #if defined HAVE_CONFIG_H
 36 #include "config.h"
 37 #endif
 38 
 39 #include "abi_common.h"
 40 
 41 
 42 subroutine ioniondist(natom,rprimd,xred,inm,option,varlist,magv,atp,prtvol)
 43 
 44  use defs_basis
 45  use m_profiling_abi
 46 
 47  use m_pptools,    only : prmat
 48 
 49 !This section has been created automatically by the script Abilint (TD).
 50 !Do not modify the following lines by hand.
 51 #undef ABI_FUNC
 52 #define ABI_FUNC 'ioniondist'
 53  use interfaces_14_hidewrite
 54  use interfaces_41_geometry, except_this_one => ioniondist
 55 !End of the abilint section
 56 
 57  implicit none
 58 
 59 !Arguments ------------------------------------
 60 !scalars
 61  integer,intent(in)              :: natom,option
 62  integer,intent(in),optional     :: atp                   !atom on which the perturbation was done
 63 !arrays
 64  real(dp),intent(in)             :: rprimd(3,3)
 65  real(dp),intent(in)             :: xred(3,natom)
 66  real(dp),intent(out)            :: inm(natom,natom)      
 67  integer,intent(in),optional     :: magv(natom)
 68  real(dp),intent(in),optional    :: varlist(natom)
 69  integer,intent(in),optional     :: prtvol
 70 
 71 !Local variables-------------------------------
 72 !scalars
 73  integer                      :: iatom,jatom,katom,kdum,atpp,prtvoll
 74  !character(len=500)           :: message
 75 !arrays
 76  integer                      :: interq(natom)
 77  real(dp)                     :: hxcart(3,natom),distm(natom,natom)
 78  real(dp)                     :: magvv(natom)
 79 
 80 ! *************************************************************************
 81 
 82  hxcart=matmul(rprimd,xred)
 83  interq=(/(iatom,iatom=1,natom)/)
 84  inm=0
 85 
 86  if (present(magv)) then
 87    magvv=magv
 88  else
 89    magvv=(/ (1, iatom=1,natom)  /)
 90  end if
 91 
 92  if (present(atp)) then
 93    atpp=atp
 94  else
 95    atpp=1
 96  end if
 97 
 98  if (present(prtvol)) then
 99    prtvoll=prtvol
100  else
101    prtvoll=1
102  end if
103 
104  if (option==3.and.(.not.present(varlist))) then
105    call  wrtout(std_out,'ioniondist error: option=3 but no variable list provided for symmetrization','COLL')
106    return
107  end if
108 
109 
110 !DEBUG
111 !write(message, '(a,a)' ) ch10,' ioniondist start '
112 !call wrtout(std_out,message,'COLL')
113 !END DEBUG
114 
115  distm=0
116  katom=atpp-1
117  do iatom=1,natom
118    katom=katom+1
119    if (katom > natom) katom=1
120    distm(iatom,iatom)=0
121    do jatom=iatom,natom
122      distm(iatom,jatom)=dist2(xred(:,katom),xred(:,jatom),rprimd,1)*magvv(katom)*magvv(jatom)
123      distm(jatom,iatom)=distm(iatom,jatom)
124    end do
125  end do
126 
127  if (prtvoll>=3) then
128    call  wrtout(std_out,'ioniondist: ionic distances:','COLL')
129    call prmat(distm,natom,natom,natom,std_out)
130  end if
131 
132  distm=anint(distm*10000_dp)/10000_dp           ! rounding needed else distm(iatom,jatom)/= distm(1,kdum) sometimes fails
133 
134  do iatom=1,natom
135    if (option==1) then
136      inm(iatom,:)=distm(iatom,:)
137    else
138      do jatom=iatom,natom
139        kdum=1
140        do while ( (kdum <= natom) .and. (distm(iatom,jatom)/= distm(1,kdum)) )
141          kdum=kdum+1
142        end do
143        if (option==2) then
144          inm(iatom,jatom)=interq(kdum)
145        else if (option==3) then
146          inm(iatom,jatom)=varlist(kdum)
147        end if
148        inm(jatom,iatom)=inm(iatom,jatom)
149      end do
150    end if
151  end do
152 
153  if (prtvoll==2) then
154    call  wrtout(std_out,'ioniondist: symmetrized matrix:','COLL')
155    call prmat(distm,1,natom,natom,std_out)
156  else if (prtvoll>=3) then
157    call  wrtout(std_out,'ioniondist: symmetrized matrix:','COLL')
158    call prmat(distm,natom,natom,natom,std_out)
159  end if
160 
161 end subroutine ioniondist