TABLE OF CONTENTS


ABINIT/cgcprj_cholesky [ Modules ]

[ Top ] [ Modules ]

NAME

  cgcprj_cholesky

FUNCTION

 Cholesky orthonormalization of the vectors stored in cg+cprj mode.

 This implementation is NOT band-parallelized
 Also, it is far of being optimal at the level of linear algebra

COPYRIGHT

 Copyright (C) 2017-2018 ABINIT group (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

  atindx1(natom)=index table for atoms, inverse of atindx
  dimcprj(natom)=number of lmn components in the <p_{lmn}^i|\psi> for the i-th atom
  icg=shift in cg array to locate current k-point and spinpol 
  ikpt=current k point index
  isppol=current spin polarization index
  istwf=input option parameter that describes the storage of wfs
  mcg=second dimension of cg array (mpw*nspinor*mband*mkmem*nsppol)
  mcprj=second dimension of cprj_k array
  mkmem=number of k points which can fit in memory
  mpi_enreg=information about MPI parallelization
  natom=number of atoms
  nattyp(ntypat)=number of atoms of each type in cell.
  nband=number of bands
  npw=number of planewaves in basis at this k point
  nspinor=number of spinor components
  nsppol=number of spin polarizations
  ntypat=number of types of atoms
  pawtab(dtset%ntypat*dtset%usepaw) <type(pawtab_type)>=paw tabulated starting data
  usepaw=1 if PAW is activated

SIDE EFFECTS

  cg(2,mcg)= plane wave wavefunction coefficients for the set of input wavefunctions (all k points and spinpol)
  cprj_k(natom,mcprj) <type(pawcprj_type)>= projected input wave functions <Proj_i|Cnk> with NL projectors for the specific k point and spinpol

PARENTS

      wf_mixing

CHILDREN

      dotprod_set_cgcprj,lincom_cgcprj,zpotrf,ztrsm

SOURCE

 53 #if defined HAVE_CONFIG_H
 54 #include "config.h"
 55 #endif
 56 
 57 #include "abi_common.h"
 58 
 59  subroutine cgcprj_cholesky(atindx1,cg,cprj_k,dimcprj,icg,ikpt,isppol,istwf,mcg,mcprj,mkmem,&
 60 &  mpi_enreg,natom,nattyp,nband,npw,nspinor,nsppol,ntypat,pawtab,usepaw)
 61 
 62  use defs_basis
 63  use defs_abitypes
 64  use m_pawtab, only : pawtab_type
 65  use m_pawcprj, only : pawcprj_type
 66 
 67 !This section has been created automatically by the script Abilint (TD).
 68 !Do not modify the following lines by hand.
 69 #undef ABI_FUNC
 70 #define ABI_FUNC 'cgcprj_cholesky'
 71  use interfaces_66_wfs, except_this_one => cgcprj_cholesky
 72 !End of the abilint section
 73 
 74  implicit none
 75 
 76 !Arguments ------------------------------------
 77 !scalars
 78  integer,intent(in) :: icg,ikpt,isppol,istwf,mcg,mcprj,mkmem
 79  integer,intent(in) :: natom,nband,npw,nspinor,nsppol,ntypat,usepaw
 80 !arrays
 81  integer, intent(in) :: atindx1(natom),dimcprj(natom),nattyp(ntypat)
 82  real(dp), intent(inout) :: cg(2,mcg) 
 83  type(pawcprj_type),intent(inout) :: cprj_k(natom,mcprj)
 84  type(MPI_type),intent(in) :: mpi_enreg
 85  type(pawtab_type),intent(in) :: pawtab(ntypat*usepaw)
 86 
 87 !Local variables ------------------------------
 88 !scalars
 89  integer :: hermitian,ierr,ii,inplace
 90 !arrays
 91  real(dp), allocatable :: dmn(:,:,:),smn(:,:,:)
 92 
 93 ! *************************************************************************
 94 
 95  ABI_ALLOCATE(smn,(2,nband,nband))
 96  ABI_ALLOCATE(dmn,(2,nband,nband))
 97 
 98  hermitian=1
 99  call dotprod_set_cgcprj(atindx1,cg,cg,cprj_k,cprj_k,dimcprj,hermitian,&
100 & 0,0,icg,icg,ikpt,isppol,istwf,nband,mcg,mcg,mcprj,mcprj,mkmem,&
101 & mpi_enreg,natom,nattyp,nband,nband,npw,nspinor,nsppol,ntypat,pawtab,smn,usepaw)
102 
103 !Cholesky factorization: O = U^H U with U upper triangle matrix.
104  call ZPOTRF('U',nband,smn,nband,ierr)
105 
106 !Solve X U = 1. 
107  dmn=zero
108  do ii=1,nband
109    dmn(1,ii,ii)=one 
110  end do
111  call ZTRSM('Right','Upper','Normal','Normal',nband,nband,cone,smn,nband,dmn,nband)
112 
113  inplace=1
114 !This call does not take into account the fact that X=dmn is an upper triangular matrix...
115 !The number of operations might be divided by two.
116  call lincom_cgcprj(dmn,cg,cprj_k,dimcprj,&
117 & icg,inplace,mcg,mcprj,natom,nband,nband,npw,nspinor,usepaw)
118 
119  ABI_DEALLOCATE(smn)
120  ABI_DEALLOCATE(dmn)
121 
122  end subroutine cgcprj_cholesky