TABLE OF CONTENTS


ABINIT/calc_fc [ Functions ]

[ Top ] [ Functions ]

NAME

 calc_fc

FUNCTION

 calculation and output of Fermi-contact term at each atomic site

COPYRIGHT

 Copyright (C) 2009-2018 ABINIT group (JWZ,MT)
 This file is distributed under the terms of the
 GNU General Public License, see ~ABINIT/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

  mpi_atmtab(:)=--optional-- indexes of the atoms treated by current proc
  comm_atom=--optional-- MPI communicator over atoms
  my_natom=number of atoms treated by current processor
  natom=number of atoms in cell.
  nspden=number of spin density components
  ntypat=number of atom types
  pawrad(ntypat) <type(pawrad_type)>=paw radial mesh and related data
  pawrhoij(natom) <type(pawrhoij_type)>= paw rhoij occupancies and related data
  pawtab(ntypat) <type(pawtab_type)>=paw tabulated starting data
  typat(natom)=type (integer) for each atom
  usepaw=1 if PAW is activated

OUTPUT

  (only writing, printing)

SIDE EFFECTS

NOTES

PARENTS

      outscfcv

CHILDREN

      free_my_atmtab,get_my_atmtab,make_fc_paw,wrtout

SOURCE

 44 #if defined HAVE_CONFIG_H
 45 #include "config.h"
 46 #endif
 47 
 48 #include "abi_common.h"
 49 
 50  subroutine calc_fc(my_natom,natom,nspden,ntypat,pawrad,pawrhoij,pawtab,typat,usepaw,&
 51  &                  mpi_atmtab,comm_atom) ! optional arguments (parallelism)
 52 
 53  use defs_basis
 54  use m_profiling_abi
 55  use m_errors
 56  use m_xmpi, only : xmpi_comm_self
 57 
 58  use m_paral_atom, only : get_my_atmtab, free_my_atmtab
 59 
 60  use m_pawrad, only : pawrad_type
 61  use m_pawtab, only : pawtab_type
 62  use m_pawrhoij, only : pawrhoij_type
 63 
 64 !This section has been created automatically by the script Abilint (TD).
 65 !Do not modify the following lines by hand.
 66 #undef ABI_FUNC
 67 #define ABI_FUNC 'calc_fc'
 68  use interfaces_14_hidewrite
 69  use interfaces_65_paw
 70 !End of the abilint section
 71 
 72  implicit none
 73 
 74 !Arguments ------------------------------------
 75 !scalars 
 76  integer,intent(in) :: my_natom,natom,nspden,ntypat,usepaw
 77  integer,optional,intent(in) :: comm_atom
 78 !arrays 
 79  integer,intent(in) :: typat(natom)
 80  integer,optional,target,intent(in) :: mpi_atmtab(:)
 81  type(pawrad_type),intent(in) :: pawrad(ntypat)
 82  type(pawrhoij_type),intent(in) :: pawrhoij(my_natom)
 83  type(pawtab_type),intent(in) :: pawtab(ntypat)
 84 
 85 !Local variables-------------------------------
 86 !scalars
 87  integer :: iatom,my_comm_atom
 88  logical :: my_atmtab_allocated,paral_atom
 89  character(len=500) :: message
 90 !arrays 
 91  integer,pointer :: my_atmtab(:)
 92  real(dp),allocatable :: fc(:,:)
 93 
 94 !***********************************************************************
 95 
 96 
 97 !Compatibility tests
 98  if (usepaw /= 1) then
 99    message = ' usepaw /= 1 but Fermi-contact calculation requires PAW '
100    MSG_ERROR(message)
101  end if
102 
103 !Set up parallelism over atoms
104  paral_atom=(present(comm_atom).and.(my_natom/=natom))
105  nullify(my_atmtab);if (present(mpi_atmtab)) my_atmtab => mpi_atmtab
106  my_comm_atom=xmpi_comm_self;if (present(comm_atom)) my_comm_atom=comm_atom
107  call get_my_atmtab(my_comm_atom,my_atmtab,my_atmtab_allocated,paral_atom,natom,my_natom_ref=my_natom)
108 
109 !Initialization
110  ABI_ALLOCATE(fc,(nspden,natom))
111 
112 !Computation
113  if (paral_atom) then
114    call make_fc_paw(fc,my_natom,natom,nspden,ntypat,pawrhoij,pawrad,pawtab,&
115 &   comm_atom=my_comm_atom,mpi_atmtab=my_atmtab)
116  else
117    call make_fc_paw(fc,my_natom,natom,nspden,ntypat,pawrhoij,pawrad,pawtab)
118  end if
119 
120 !Printing
121  write(message,'(a,a,a)' ) ch10,' Fermi-contact Term Calculation ',ch10
122  call wrtout(ab_out,message,'COLL')
123 
124  do iatom = 1, natom
125    if (nspden == 2) then
126      write(message,'(a,i3,a,i3,a,f12.4)') ' Atom ',iatom,', typat ',typat(iatom),': FC total = ',&
127 &     fc(1,iatom)+fc(2,iatom)
128      call wrtout(ab_out,message,'COLL')
129      write(message,'(a,i3,a,i3,a,f12.4)') ' Atom ',iatom,', typat ',typat(iatom),': FC up - down = ',&
130 &     fc(1,iatom)-fc(2,iatom)
131      call wrtout(ab_out,message,'COLL')
132    else
133      write(message,'(a,i3,a,i3,a,f12.4)') ' Atom ',iatom,', typat ',typat(iatom),': FC = ',&
134 &     fc(1,iatom)
135      call wrtout(ab_out,message,'COLL')
136    end if
137  end do
138 
139  write(message,'(3a)')ch10,ch10,ch10
140  call wrtout(ab_out,message,'COLL')
141 
142 !Memory deallocation
143  ABI_DEALLOCATE(fc)
144 
145 !Destroy atom table used for parallelism
146  call free_my_atmtab(my_atmtab,my_atmtab_allocated)
147 
148 end subroutine calc_fc