TABLE OF CONTENTS


ABINIT/ortho_reim [ Functions ]

[ Top ] [ Functions ]

NAME

 ortho_reim

FUNCTION

 This routine computes the overlap of two wavefunctions (for a given number of bands)
 and orthonormalizes it:
      - Computes the products of two rectangular matrices
         containing the wavefunctions psi and S.psi (where S is the
         overlap (with the PAW terms if necessary)).
      - Does a Cholesky decomposition of this overlap
      - rotates the initial matrix blockvectorx by the triangular matrix to
         have an orthonormal set of wavefunctions

 This version operates on arrays in which the real and the imaginary part
 are packed together (real parts first, them imaginary parts), used when istwfk=2

INPUTS

  blockvectorbx = matrix of dimension (blocksize,vectsize)
                  (e.g. block of overlap*wavefunction)
  blocksize     = dimension of matrices (e.g number of bands)
  spaceComm     = communicator used for MPI parallelization
  vectsize      = dimension of matrices (e.g number of G vector)

OUTPUT

  sqgram        = Choleski decomposition of transpose(blockvector)*blockvectorx

SIDE EFFECTS

  blockvectorx  = on input, matrix of dimension (vectsize,blocksize)
                  (e.g block of wavefunction)
  blockvectorx  = on output, orthonormalized wavefunction.

SOURCE

142 subroutine ortho_reim(blockvectorx,blockvectorbx,blocksize,spaceComm,sqgram,vectsize)
143 
144 !Arguments ------------------------------------
145 !scalars
146  integer,intent(in) :: blocksize,vectsize,spaceComm
147 !arrays
148  real(dp),intent(in) :: blockvectorbx(vectsize,blocksize)
149  real(dp),intent(inout) :: blockvectorx(vectsize,blocksize)
150  real(dp),intent(out) :: sqgram(blocksize,blocksize)
151 
152 !Local variables-------------------------------
153 !scalars
154  integer :: ierr,info
155  character(len=500) :: message
156 
157 ! *********************************************************************
158 
159  call abi_xgemm('t','n',blocksize,blocksize,vectsize,cone,blockvectorx,&
160 &   vectsize,blockvectorbx,vectsize,czero,sqgram,blocksize)
161 
162  call xmpi_sum(sqgram,spaceComm,ierr)
163 
164  !Cholesky factorization of sqgram (ouside upper Triangular of sqgram)
165  call abi_d2zpotrf('u',blocksize,sqgram,blocksize,info) !vz_d
166 
167  if (info /= 0 )  then
168    write(message,'(a,i0)')'dpotrf, info=',info
169    ABI_ERROR(message)
170  end if
171 
172 !Find X  X*sqgram=blockvectorx
173  call abi_xtrsm('r','u','n','n',vectsize,blocksize,one,sqgram,blocksize,blockvectorx,vectsize)
174 
175 end subroutine ortho_reim

ABINIT/zorthonormalize [ Functions ]

[ Top ] [ Functions ]

NAME

 zorthonormalize

FUNCTION

 This routine computes the overlap of two complex wavefunctions (for a given number of bands)
 and orthonormalizes it:
      - Computes the products of two rectangular matrices
         containing the wavefunctions psi and S.psi (where S is the
         overlap (with the PAW terms if necessary)).
      - Does a Cholesky decomposition of this overlap
      - rotates the initial matrix blockvectorx by the triangular matrix to
         have an orthonormal set of wavefunctions

INPUTS

  blockvectorbx = matrix of dimension (blocksize,vectsize)
                  (e.g. block of overlap*wavefunction)
  blocksize     = dimension of matrices (e.g number of bands)
  spaceComm     = communicator used for MPI parallelization
  vectsize      = dimension of matrices (e.g number of G vector)

OUTPUT

  sqgram        = Choleski decomposition of transpose(blockvector)*blockvectorx

SIDE EFFECTS

  blockvectorx  = on input, matrix of dimension (vectsize,blocksize)
                  (e.g block of wavefunction)
  blockvectorx  = on output, orthonormalized wavefunction.

SOURCE

211 subroutine zorthonormalize(blockvectorx,blockvectorbx,blocksize,spaceComm,sqgram,vectsize)
212 
213 !Arguments ------------------------------------
214 !scalars
215  integer,intent(in) :: blocksize,spaceComm,vectsize
216 !arrays
217  complex(dpc),intent(in) :: blockvectorbx(vectsize,blocksize)
218  complex(dpc),intent(inout) :: blockvectorx(vectsize,blocksize)
219  complex(dpc),intent(out) :: sqgram(blocksize,blocksize)
220 
221 !Local variables-------------------------------
222 !scalars
223  integer :: ierr,info
224  character(len=500) :: message
225 
226 ! *********************************************************************
227 
228  call abi_xgemm('c','n',blocksize,blocksize,vectsize,cone,blockvectorx,&
229 & vectsize,blockvectorbx,vectsize,czero,sqgram,blocksize)
230 
231  call xmpi_sum(sqgram,spaceComm,ierr)
232 
233  call abi_xpotrf('u',blocksize,sqgram,blocksize,info)
234 
235  if (info /= 0 )  then
236    write(message,'(a,i0)')'zpotrf, info=',info
237    ABI_ERROR(message)
238  end if
239 
240  call abi_xtrsm('r','u','n','n',vectsize,blocksize,cone,sqgram,blocksize,blockvectorx,vectsize)
241 
242 end subroutine zorthonormalize

m_abi_linalg/abi_xorthonormalize [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

  abi_xorthonormalize

FUNCTION

  abi_xorthonormalize is the generic function for computing the
  overlap of two complex wavefunctions (for a given number of bands)
  and orthonormalizes it:

COPYRIGHT

  Copyright (C) 2001-2024 ABINIT group (LNguyen,FDahm (CS), FBottin, GZ, AR, MT)
  This file is distributed under the terms of the
  GNU General Public License, see ~ABINIT/Infos/copyright
  or http://www.gnu.org/copyleft/gpl.txt .

SOURCE


m_abi_linalg/xorthonormalize [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

  xorthonormalize

FUNCTION

  This routine computes the overlap of two complex wavefunctions (for a given number of bands)
  and orthonormalizes it:
      - Computes the products of two rectangular matrices
         containing the wavefunctions psi and S.psi (where S is the
         overlap (with the PAW terms if necessary)).
      - Does a Cholesky decomposition of this overlap
      - rotates the initial matrix blockvectorx by the triangular matrix to
         have an orthonormal set of wavefunctions

INPUTS

  blockvectorbx = matrix of dimension (blocksize,vectsize)
                  (e.g. block of overlap*wavefunction)
  blocksize     = dimension of matrices (e.g number of bands)
  spaceComm     = communicator used for MPI parallelization
  vectsize      = dimension of matrices (e.g number of G vector)

OUTPUT

  sqgram        = Choleski decomposition of transpose(blockvector)*blockvectorx

SIDE EFFECTS

  blockvectorx  = on input, matrix of dimension (vectsize,blocksize)
                  (e.g block of wavefunction)
  blockvectorx  = on output, orthonormalized wavefunction.

SOURCE

 53 subroutine xorthonormalize(blockvectorx,blockvectorbx,blocksize,spaceComm,sqgram,vectsize,&
 54 &                          x_cplx,timopt,tim_xortho) ! optional arguments
 55 
 56 !Arguments ------------------------------------
 57 !scalars
 58  integer,intent(in) :: blocksize,vectsize,spaceComm,x_cplx
 59  integer, intent(in), optional :: timopt,tim_xortho
 60  !arrays
 61  real(dp),intent(in) :: blockvectorbx(vectsize,blocksize)
 62  real(dp),intent(inout) :: blockvectorx(vectsize,blocksize)
 63  real(dp),intent(out) :: sqgram(x_cplx*blocksize,blocksize)
 64 
 65 !Local variables-------------------------------
 66  real(dp) :: tsec(2)
 67  integer  :: ierr,info
 68  character(len=500) :: message
 69  character, dimension(2) :: cparam
 70 
 71  ! *********************************************************************
 72 
 73  if (present(tim_xortho).and.present(timopt)) then
 74    if(abs(timopt)==3) then
 75      call timab(tim_xortho,1,tsec)
 76    end if
 77  end if
 78 
 79  cparam(1)='t'
 80  cparam(2)='c'
 81 
 82  call abi_xgemm(cparam(x_cplx),'n',blocksize,blocksize,vectsize,cone,blockvectorx,&
 83 &   vectsize,blockvectorbx,vectsize,czero,sqgram,blocksize,x_cplx=x_cplx)
 84 
 85  call xmpi_sum(sqgram,spaceComm,ierr)
 86 
 87  !Cholesky factorization of sqgram (ouside upper Triangular of sqgram)
 88  call abi_xpotrf('u',blocksize,sqgram,blocksize,info,x_cplx=x_cplx)
 89 
 90  if (info /= 0 )  then
 91    write(message,'(a,i0)')'abi_xpotrf, info=',info
 92    ABI_ERROR(message)
 93  end if
 94 
 95  !Find X  X*sqgram=blockvectorx
 96  call abi_xtrsm('r','u','n','n',vectsize,blocksize,cone,sqgram,blocksize,&
 97 &   blockvectorx,vectsize,x_cplx=x_cplx)
 98 
 99  if (present(tim_xortho).and.present(timopt)) then
100    if(abs(timopt)==3) then
101      call timab(tim_xortho,2,tsec)
102    end if
103  end if
104 
105 end subroutine xorthonormalize