TABLE OF CONTENTS


ABINIT/lincom_cgcprj [ Functions ]

[ Top ] [ Functions ]

NAME

 lincom_cgcprj

FUNCTION

 For one k point and spinpol, compute a set (size nband_out) of linear combinations of nband_in wavefunctions,
 that are known in the cg+cprj representation :
 cgout_n(:,:) <--- Sum_m [ cg_m(:,:) . alpha_mn ]
 cprjout_n(:,:) <--- Sum_m [ cprj_m(:,:) . alpha_mn ]
 If nband_out is smaller or equal to nband_in, the result might be in-place (output in cg instead of cgout, and in cprj instead of cprjout).
 Otherwise, it is contained in the optional cgout+cprjout pair.
 In the present status, the cg and cgout relates to all the k points and spins, and rely on the icg index,
 while it is assumed that cprj and cprjout refer to the specific k point and spin.
 This is not coherent.
 THIS MIGHT BE CHANGED IN THE FUTURE !
 This implementation is NOT band-parallelized
 Also, it is far of being optimal at the level of linear algebra, and involves extra copying
 that are detrimental for performance...

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

  alpha_mn(2,nband_in,nband_out)=complex matrix of coefficients of the linear combinations to be computed
  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 (for input, and possibly for in-place output)
  inplace= if 0, output in cgout and cprjout ; if 1, output in cg and cprj
  mcg=second dimension of cg array (mpw*nspinor*mband*mkmem*nsppol)
  mcprj=second dimension of cprj array 
  natom=number of atoms
  nband_in=number of bands, size of the input set of wavefunctions
  nband_out=number of bands, size of the output set of wavefunctions (should be equal to nband_in if inplace==1)
  npw=number of planewaves in basis at this k point
  nspinor=number of spinor components
  usepaw=1 if PAW is activated
  [icgout= shift in cgout array to locate current k-point and spinpol (for output)]
  [mcgout=second dimension of cgout array (mpw*nspinor*mband*mkmem*nsppol)]
  [mcprjout=second dimension of cprjout array] 

OUTPUT

  [cgout(2,mcgout)= plane wave wavefunction coefficients for the set of output wavefunctions]
  [cprjout(natom,mcprjout) <type(pawcprj_type)>= projected output wave functions <Proj_i|Cnk> with NL projectors]

SIDE EFFECTS

  (this quantities are input, and possibly updated output when inplace==1)
  cg(2,mcg)= plane wave wavefunction coefficients for the set of input wavefunctions (all k points and spinpol)
  cprj(natom,mcprj) <type(pawcprj_type)>= projected input wave functions <Proj_i|Cnk> with NL projectors 

PARENTS

      cgcprj_cholesky,wf_mixing

CHILDREN

      pawcprj_alloc,pawcprj_free,pawcprj_lincom,zgemm

SOURCE

 65 #if defined HAVE_CONFIG_H
 66 #include "config.h"
 67 #endif
 68 
 69 #include "abi_common.h"
 70 
 71  subroutine lincom_cgcprj(alpha_mn,cg,cprj,dimcprj,&
 72 & icg,inplace,mcg,mcprj,natom,nband_in,nband_out,npw,nspinor,usepaw, & 
 73 & cgout,cprjout,icgout) ! optional args
 74 
 75  use defs_basis
 76  use m_errors 
 77  use m_profiling_abi
 78  use m_pawcprj, only : pawcprj_type, pawcprj_alloc, pawcprj_lincom, pawcprj_free
 79 
 80 !This section has been created automatically by the script Abilint (TD).
 81 !Do not modify the following lines by hand.
 82 #undef ABI_FUNC
 83 #define ABI_FUNC 'lincom_cgcprj'
 84 !End of the abilint section
 85 
 86  implicit none
 87 
 88 !Arguments ------------------------------------
 89 !scalars
 90  integer, intent(in) :: icg,inplace,mcg,mcprj
 91  integer, intent(in) :: natom,nband_in,nband_out,npw,nspinor,usepaw
 92  integer, intent(in),optional :: icgout
 93 !arrays
 94  integer, intent(in) :: dimcprj(natom)
 95  real(dp), intent(inout) :: cg(2,mcg)
 96  real(dp), intent(in) :: alpha_mn(2,nband_in,nband_out)
 97  real(dp), intent(out),optional :: cgout(:,:)
 98  type(pawcprj_type),intent(inout) :: cprj(natom,mcprj)
 99  type(pawcprj_type),intent(out),optional :: cprjout(:,:)
100 
101 !Local variables-------------------------------
102 !scalars
103  integer :: iband_in,iband_out,ii
104 !arrays
105  real(dp),allocatable :: al(:,:),cgout_(:,:)
106  type(pawcprj_type),allocatable :: cprjout_(:,:)
107 
108 ! *************************************************************************
109 
110 !DEBUG
111 !write(std_out,*)' lincom_cgcprj : enter '
112 !write(std_out,*)' lincom_cgcprj : npw, nspinor=',npw,nspinor
113 !write(std_out,*)' lincom_cgcprj : icgout=',icgout
114 !ENDDEBUG
115 
116  if(inplace==0)then
117    if(.not.present(cgout))then
118      MSG_ERROR(' inplace==0 while .not.present(cgout) is not permitted ')
119    end if
120    if(usepaw==1) then
121      if(.not.present(cprjout))then 
122        MSG_ERROR(' inplace==0 and usepaw==1 while .not.present(cprjout) is not permitted ')
123      end if
124    end if
125  end if
126 
127 !Take care of the plane wave part
128  ABI_ALLOCATE(cgout_,(2,npw*nspinor*nband_out))
129 
130  call zgemm('N','N',npw*nspinor,nband_out,nband_in,dcmplx(1._dp), &
131 & cg(:,icg+1:icg+npw*nspinor*nband_in),npw*nspinor, &
132 & alpha_mn,nband_in,dcmplx(0._dp),cgout_,npw*nspinor)
133 
134  if(inplace==1)then
135    cg(:,icg+1:icg+npw*nspinor*nband_out)=cgout_
136  else
137    cgout(:,icgout+1:icgout+npw*nspinor*nband_out)=cgout_
138  end if
139  ABI_DEALLOCATE(cgout_)
140 
141 !Take care of the cprj part
142  if(usepaw==1) then
143 
144    ABI_DATATYPE_ALLOCATE(cprjout_,(natom,nspinor*nband_out))
145    call pawcprj_alloc(cprjout_,cprj(1,1)%ncpgr,dimcprj)
146    ABI_ALLOCATE(al,(2,nband_in))
147    do iband_out=1,nband_out
148      ii=(iband_out-1)*nspinor
149      do iband_in=1,nband_in
150        al(1,iband_in)=alpha_mn(1,iband_in,iband_out)
151        al(2,iband_in)=alpha_mn(2,iband_in,iband_out)
152      end do
153      call pawcprj_lincom(al,cprj,cprjout_(:,ii+1:ii+nspinor),nband_in)
154    end do
155    ABI_DEALLOCATE(al)
156 
157    if(inplace==1)then
158      cprj=cprjout_
159    else
160      cprjout=cprjout_
161    end if
162    call pawcprj_free(cprjout_)
163    ABI_DATATYPE_DEALLOCATE(cprjout_)
164 
165  end if
166 
167 end subroutine lincom_cgcprj