TABLE OF CONTENTS


ABINIT/chiscwrt [ Functions ]

[ Top ] [ Functions ]

NAME

  chiscwrt 

FUNCTION

  Distributes values of chi_org on chi_sc according to ion-ion distances and multiplicities in shells

INPUTS

  chi_org  = response matrix in original unit cell
  disv_org = distances (multiplied by magntization) in original unit cell
  nat_org = number of atoms in original unit cell
  sdisv_org = radii of shells in original unit cell
  smult_org = multiplicities of shells in original unit cell
  nsh_org   = number of shells in original unit cell
  disv_sc   = distances (multiplied by magntization) in super-cell
  nat_sc  = number of atoms in super-cell
  sdisv_sc = radii of shells in super-cell (was unused, so suppressed)
  smult_sc = multiplicities of shells in super-cell
  nsh_sc = number of shells in super-cell
  opt =

COPYRIGHT

  Copyright (C) 1998-2017 ABINIT group (DJA)
  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 .

OUTPUT

  chi_sc = first line of reponse matrix in super-cell

SIDE EFFECTS

PARENTS

      pawuj_det

CHILDREN

      prmat,wrtout

SOURCE

 43 #if defined HAVE_CONFIG_H
 44 #include "config.h"
 45 #endif
 46 
 47 #include "abi_common.h"
 48 
 49 
 50 subroutine chiscwrt(chi_org,disv_org,nat_org,sdisv_org,smult_org,nsh_org,chi_sc,&
 51 & disv_sc,nat_sc,smult_sc,nsh_sc,opt,prtvol) 
 52 
 53  use defs_basis
 54  use m_profiling_abi
 55 
 56  use m_pptools,    only : prmat
 57 
 58 !This section has been created automatically by the script Abilint (TD).
 59 !Do not modify the following lines by hand.
 60 #undef ABI_FUNC
 61 #define ABI_FUNC 'chiscwrt'
 62  use interfaces_14_hidewrite
 63 !End of the abilint section
 64 
 65  implicit none
 66 
 67 !Arguments ------------------------------------
 68 !scalars
 69  integer,intent(in)              :: nat_org,nsh_org,nsh_sc 
 70  integer,intent(in)              :: nat_sc
 71  integer,intent(in),optional     :: prtvol                  
 72  integer,intent(in),optional     :: opt
 73 !arrays
 74  real(dp),intent(in)             :: chi_org(nat_org),disv_org(nat_org),sdisv_org(nsh_org)
 75  integer,intent(in)              :: smult_org(nsh_org), smult_sc(nsh_sc)
 76  real(dp),intent(in)             :: disv_sc(nat_sc)
 77  real(dp),intent(out)            :: chi_sc(nat_sc)
 78 
 79 !Local variables-------------------------------
 80 !scalars
 81  integer                      :: iatom,jatom,jsh,optt
 82  character(len=500)           :: message
 83 !arrays
 84  real(dp)                     :: chi_orgl(nat_org)
 85 
 86 ! *************************************************************************
 87 
 88  if (present(opt)) then
 89    optt=opt
 90  else
 91    optt=1
 92  end if
 93 
 94 
 95  do iatom=1,nat_org
 96    do jsh=1,nsh_org
 97      if (disv_org(iatom)==sdisv_org(jsh)) then
 98        if (opt==2) then
 99          chi_orgl(iatom)=chi_org(iatom)
100        else if (opt==1) then
101          chi_orgl(iatom)=chi_org(iatom)*smult_org(jsh)/smult_sc(jsh)
102        end if 
103        exit
104      end if
105    end do !iatom
106  end do  !jsh
107  
108  if (prtvol>1) then
109    write(message,fmt='(a,150f10.5)')' chiscwrt: chi at input ',chi_org
110    call wrtout(std_out,message,'COLL')
111    write(message,fmt='(a,150f10.5)')' chiscwrt: chi after division ',chi_orgl
112    call wrtout(std_out,message,'COLL')
113  end if
114 
115  do iatom=1,nat_sc
116    do jatom=1,nat_org
117      if (disv_org(jatom)==disv_sc(iatom)) then
118        chi_sc(iatom)=chi_orgl(jatom)
119        exit
120      else if (jatom==nat_org) then
121        chi_sc(iatom)=0_dp
122      end if
123    end do
124  end do
125 
126  if (prtvol>1) then
127    write(message,'(a)')' chiscwrt, chi_sc '
128    call wrtout(std_out,message,'COLL')
129    call prmat(chi_sc,1,nat_sc,1,std_out)
130  end if
131 
132 end subroutine chiscwrt