mpi_lanczos_iteration_d Subroutine

subroutine mpi_lanczos_iteration_d(MpiComm, MatVec, iter, vin, vout, alfa, beta)

Arguments

Type IntentOptional Attributes Name
integer :: MpiComm
subroutine MatVec(Nloc, vin, vout)
Arguments
Type IntentOptional Attributes Name
integer :: Nloc
real(kind=8), dimension(Nloc) :: vin
real(kind=8), dimension(Nloc) :: vout
integer :: iter
real(kind=8), intent(inout), dimension(:) :: vin
real(kind=8), intent(inout), dimension(size(vin)) :: vout
real(kind=8), intent(inout) :: alfa
real(kind=8), intent(inout) :: beta

Calls

proc~~mpi_lanczos_iteration_d~2~~CallsGraph proc~mpi_lanczos_iteration_d~2 mpi_lanczos_iteration_d allreduce_mpi allreduce_mpi proc~mpi_lanczos_iteration_d~2->allreduce_mpi get_master_mpi get_master_mpi proc~mpi_lanczos_iteration_d~2->get_master_mpi

Source Code

subroutine mpi_lanczos_iteration_d(MpiComm,MatVec,iter,vin,vout,alfa,beta)
  integer                                    :: MpiComm
  interface
     subroutine MatVec(Nloc,vin,vout)
       integer                 :: Nloc
       real(8),dimension(Nloc) :: vin
       real(8),dimension(Nloc) :: vout
     end subroutine MatVec
  end interface
  real(8),dimension(:),intent(inout)         :: vin !Nloc
  real(8),dimension(size(vin)),intent(inout) :: vout
  real(8),dimension(size(vin))               :: tmp
  real(8),intent(inout)                      :: alfa,beta
  real(8)                                    :: atmp,btmp
  integer                                    :: iter,nloc
  real(8)                                    :: norm,norm_tmp
  !
  logical                                    :: mpi_master
  !
  if(MpiComm==Mpi_Comm_Null)return
  !
  nloc=size(vin)
  !
  mpi_master=get_master_MPI(MpiComm)
  !
  if(iter==1)then
     norm_tmp=dot_product(vin,vin)
     norm = 0d0 ; call AllReduce_MPI(MpiComm,norm_tmp,norm)
     if(mpi_master.AND.norm==0d0)stop "MPI_LANCZOS_ITERATION_D: norm = 0!!"
     vin=vin/sqrt(norm)
  else
     tmp = vin
     vin = vout/beta
     vout= -beta*tmp
  endif
  call MatVec(nloc,vin,tmp)
  vout = vout + tmp
  atmp = dot_product(vin,vout) ; alfa = 0d0; call AllReduce_MPI(MpiComm,atmp,alfa)
  vout = vout - alfa*vin
  btmp = dot_product(vout,vout); beta = 0d0; call AllReduce_MPI(MpiComm,btmp,beta)
  beta = sqrt(beta)
end subroutine mpi_lanczos_iteration_d