Class quark_ovlp_kern_eigmodes_class
In: QuarkDwOvlpClass/quark_ovlp_kern_eigmodes_class.F90
QuarkDwOvlpClass_v0.1/quark_ovlp_kern_eigmodes_class.F90
comlib lattice_class logfile_class timer_class error_class field_fermion_class field_gauge_class field_5dfermion_class quark_wilson_class solver_class counter_class ks_class f95_lapack quark_ovlp_kern_eigmodes_class dot/f_23.png

Defines low/high mode container for overlap fermion kernel operator

Kernel operator definition (Hw)

 Hw = gamma5 Dw(m0) / ( Dw(m0)a5/2 + 1 ),

 Dw(m0) = 4 - m0 - 1/2 M_{4dhop},

 M_{4dhop} = sum_{mu=1,4}[ (1-gamma_mu)U_mu(n)        delta_{n+mu,m}
                          +(1+gamma_mu)U_mu(n-mu)^dag delta_{n-mu,m} ].

Version

$Id: quark_ovlp_kern_eigmodes_class.F90,v 1.4 2011/01/13 13:30:06 ishikawa Exp $

Methods

Included Modules

comlib lattice_class logfile_class timer_class error_class field_fermion_class field_gauge_class field_5dfermion_class quark_wilson_class solver_class counter_class ks_class f95_lapack

Public Instance methods

SIGN_KERNEL_HIGH_MODES
Constant :
SIGN_KERNEL_HIGH_MODES = 1 :integer, parameter
SIGN_KERNEL_HIGH_MODES
Constant :
SIGN_KERNEL_HIGH_MODES = 1 :integer, parameter
SIGN_KERNEL_LOW_MODES
Constant :
SIGN_KERNEL_LOW_MODES = 0 :integer, parameter
SIGN_KERNEL_LOW_MODES
Constant :
SIGN_KERNEL_LOW_MODES = 0 :integer, parameter
Subroutine :
this :type(quark_ovlp_kern_eigmodes), intent(in)
iout :integer, intent(in)
tol :real(DP), intent(in)
iiter :integer, intent(inout)
Hv :type(field_quark_wg), intent(inout)
v :type(field_quark_wg), intent(in)
u :type(vfield_gluon_wg), intent(in)
: even/odd site ordered gluon field

Multiply inverse Hermitian sign function kernel operaotr

Hv = Hw \ v

where

 Hw = gamma_5 Dw / (Dw a5/2 + 1)

or

 Hw = gamma_5 Dw

 Hw : the kernel for sign function

Hw \ v = (Dw a5/2 + 1 ) / Dw gamma_5 v

Dw = 4 - m0 - 1/2 Mhoping

m0 : negative mass term (m0 >0) a5 : kernel parameter

[Source]

subroutine assign_inv_mult_sign_kernel(this,iout,tol,iiter,Hv,v,u)
!
! Multiply inverse Hermitian sign function kernel operaotr
!
! Hv = Hw \ v 
!
! where
!
!  Hw = gamma_5 Dw / (Dw a5/2 + 1) 
!
! or 
!
!  Hw = gamma_5 Dw
!
!  Hw : the kernel for sign function
!
! Hw \ v = (Dw a5/2 + 1 ) / Dw gamma_5 v
!
! Dw = 4 - m0 - 1/2 Mhoping
!
! m0 : negative mass term (m0 >0)
! a5 : kernel parameter
!
  use solver_class
  use counter_class
  implicit none
  type(quark_ovlp_kern_eigmodes), intent(in) :: this
  integer,  intent(in) :: iout
  integer,  intent(inout) :: iiter
  real(DP), intent(in) :: tol
  type(field_quark_wg),    intent(inout) :: Hv
  type(field_quark_wg),    intent(in)    :: v
  type(vfield_gluon_wg),   intent(in)    :: u   ! even/odd site ordered gluon field

#ifdef _EOPREC_LOWMODE_SOLVER_
  integer, parameter :: EIGEN_NSITE=COL*SPIN*NTH*NZ*NY*NX
#else
  integer, parameter :: EIGEN_NSITE=COL*SPIN*NT*NZ*NY*NX
#endif

  type(timer) :: solver_time
  type(field_quark_wg), allocatable :: w,s,b
  type(quark_wilson) :: mw
  type(cg_alg), allocatable :: solver
  character(CHARLEN) :: str
  integer    :: ieo,ioe, mult_count
  type(counter) :: imult
  real(DP) :: M,etime,rtmp,kappa
  real(DP) :: a5,m0

  call new(solver_time)
  call tic(solver_time)

  a5 = this%a5
  m0 = this%m0
  imult = mult_iter

  allocate(w,s)
  call new(w)
  call new(s)

#ifdef _EOPREC_LOWMODE_SOLVER_
  ieo = 0
  ioe = 1-ieo

  M = m0
  rtmp  = 1.0_DP/(4.0_DP-M)
  kappa = rtmp/2.0d0
  call set_kappa(mw,kappa)

  call assign_mult(s,v,rtmp)   ! s := 1/(4-M) v
  call accum_mult_gamma5(s)    ! s := gamma5 s
  call assign_mult_hop(mw,w%eo(ieo),s%eo(ioe),u)  ! we = Meo so
  call accum_add_mult(s%eo(ieo),w%eo(ieo),kappa)  ! se = se + kappa * we = se + kappa Meo so

  !      v
  ! we = Dee^dag se
  !
  call accum_mult_gamma5(s%eo(ieo))
  call assign_mult_eoprec_wd(mw,w%eo(ieo),s%eo(ieo),u)
  call accum_mult_gamma5(w%eo(ieo))

#else
  !-----------------------------------------
  !  w  = gamma_5 Dw v = Dw^dag gamma5 v
  !-----------------------------------------
  M = m0
  call assign(s,v)
  call assign_mult_wd(M,w,s,u)
  call  accum_mult_gamma5(w)
#endif

  !=============================================
  ! Solve:
  !
  !   (Dw^dag Dw) x = w
  !
  ! for x
  !
  !=============================================
  allocate(solver)
  call new(solver,NSIZE=EIGEN_NSITE,mode=MODE_NORMAL,max_iter=iiter,tol=tol)

#ifdef _EOPREC_LOWMODE_SOLVER_
  call pack(w%eo(ieo),solver%src_vec(:))
#else
  call pack(w,solver%src_vec(:))
#endif

  !=======================
  ! start CG Solver
  !=======================
  do 
    call solve(solver)
    select case(get_status(solver))
    case (OP_NOP)
      cycle
    case (OP_DO_MATVEC)

      !---------------------------------------------
      ! multiply gamma_5 Dw gamma_5 Dw
      !---------------------------------------------
#ifdef _EOPREC_LOWMODE_SOLVER_
      M = m0
      call unpack(solver%src_vec,w%eo(ieo))
      call assign_mult_eoprec_wd(mw,s%eo(ieo),w%eo(ieo),u)   ! s := D w
      call  accum_mult_gamma5(s%eo(ieo))                     ! s := gamma_5 s
      call assign_mult_eoprec_wd(mw,w%eo(ieo),s%eo(ieo),u)   ! w := D s
      call  accum_mult_gamma5(w%eo(ieo))                     ! w := gamma_5 w
      call pack(w%eo(ieo),solver%dst_vec)
#else
      M = m0
      call unpack(solver%src_vec,w)
      call assign_mult_wd(M,s,w,u)    ! s := D w
      call  accum_mult_gamma5(s)      ! s := gamma_5 s
      call assign_mult_wd(M,w,s,u)    ! w := D s
      call  accum_mult_gamma5(w)      ! w := gamma_5 w
      call pack(w,solver%dst_vec)
#endif

    case (OP_PRINT_STATUS)
      if (this%debug) then
      if (nodeid == 0) then
        mult_count = get_count(mult_iter) - get_count(imult)
        write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)') get_current_iteration(solver), get_residual_norm(solver),mult_count
        write(*,'("#",I5," ERR_=",E24.16," MULT_=",I10)') get_current_iteration(solver), get_residual_norm(solver),mult_count
      endif
      endif
      cycle
    case (OP_CONVERGED)
      if (nodeid == 0) then
        mult_count = get_count(mult_iter) - get_count(imult)
        write(iout,'("#",A," iteration converged.")')TRIM(get_name(solver))
        write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), get_residual_norm(solver),mult_count
      endif
      exit      
    case (OP_MAXITER_REACHED)
      if (nodeid == 0) then
        mult_count = get_count(mult_iter) - get_count(imult)
        write(iout,'("#",A," iteration does not converge.")')TRIM(get_name(solver))
        write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), get_residual_norm(solver),mult_count
      endif
      write(str,'("solver did not converge.: ",A,I4)')__FILE__,__LINE__
      call error_stop(str)
    case default
      write(str,'("something is wrong in solver.: ",A,I4)')__FILE__,__LINE__
      call error_stop(str)
    end select
  enddo

#ifdef _EOPREC_LOWMODE_SOLVER_
  call unpack(solver%dst_vec,w%eo(ieo))
  call assign_mult_hop(mw,w%eo(ioe),w%eo(ieo),u)    ! wo = Moe we
  call  accum_mult_add(w%eo(ioe),kappa,s%eo(ioe))   ! wo = kappa * wo + so = kappa Moe we + so
#else
  call unpack(solver%dst_vec,w)
#endif

  if (this%debug) then
  M = m0
  call assign_mult_wd(M,s,w,u)   ! s := Dw w
  call accum_mult_gamma5(s)      ! s := gamma5 s
  allocate(b)
  call new(b)
  call assign_sub(b,v,s)         ! b = v - s
  rtmp = abs2(b)/abs2(v)
  rtmp = sqrt(rtmp)
  deallocate(b)
  if (nodeid == 0) then
      mult_count = get_count(mult_iter) - get_count(imult)
      write(*,*)
      write(*   ,'("%",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), rtmp,mult_count
      write(iout,*)
      write(iout,'("%",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), rtmp,mult_count
  endif
  endif

  !--------------------------------
  ! set solution on Hv
  !--------------------------------
  if (abs(a5) > EPSILON(a5)) then

    !
    ! Hv = (a5/2) (Dw + 2/a5) (gamma_5 Dw) \ v
    !
    M = m0 - 2.0_DP/a5
    call assign_mult_wd(M,Hv,w,u)
    call accum_mult(Hv,a5/2.0_DP)

  else

    !
    ! Hv = w
    !
    call assign(Hv,w)

  endif

  call delete(solver)
  deallocate(solver)
  deallocate(w,s)

  call toc(solver_time)

  !============================================
  ! NOTE: get_elapse contains comlib_sumcast.
  ! This means that get_elapse should not used
  ! in if-(nodeid==0)-then-endif block.
  !============================================
  etime = get_elapse(solver_time)  
  if (nodeid == 0) then
    if (this%debug) write(*   ,'(" ETIME_=",E24.16)')etime

    write(iout,'(" ETIME_=",E24.16)')etime
  endif

  return
end subroutine
Subroutine :
this :type(quark_ovlp_kern_eigmodes), intent(in)
iout :integer, intent(in)
tol :real(DP), intent(in)
iiter :integer, intent(inout)
Hv :type(field_quark_wg), intent(inout)
v :type(field_quark_wg), intent(in)
u :type(vfield_gluon_wg), intent(in)
: even/odd site ordered gluon field

Multiply inverse Hermitian sign function kernel operaotr

Hv = Hw \ v

where

 Hw = gamma_5 Dw / (Dw a5/2 + 1)

or

 Hw = gamma_5 Dw

 Hw : the kernel for sign function

Hw \ v = (Dw a5/2 + 1 ) / Dw gamma_5 v

Dw = 4 - m0 - 1/2 Mhoping

m0 : negative mass term (m0 >0) a5 : kernel parameter

[Source]

subroutine assign_inv_mult_sign_kernel(this,iout,tol,iiter,Hv,v,u)
!
! Multiply inverse Hermitian sign function kernel operaotr
!
! Hv = Hw \ v 
!
! where
!
!  Hw = gamma_5 Dw / (Dw a5/2 + 1) 
!
! or 
!
!  Hw = gamma_5 Dw
!
!  Hw : the kernel for sign function
!
! Hw \ v = (Dw a5/2 + 1 ) / Dw gamma_5 v
!
! Dw = 4 - m0 - 1/2 Mhoping
!
! m0 : negative mass term (m0 >0)
! a5 : kernel parameter
!
  use solver_class
  use counter_class
  implicit none
  type(quark_ovlp_kern_eigmodes), intent(in) :: this
  integer,  intent(in) :: iout
  integer,  intent(inout) :: iiter
  real(DP), intent(in) :: tol
  type(field_quark_wg),    intent(inout) :: Hv
  type(field_quark_wg),    intent(in)    :: v
  type(vfield_gluon_wg),   intent(in)    :: u   ! even/odd site ordered gluon field

#ifdef _EOPREC_LOWMODE_SOLVER_
  integer, parameter :: EIGEN_NSITE=COL*SPIN*NTH*NZ*NY*NX
#else
  integer, parameter :: EIGEN_NSITE=COL*SPIN*NT*NZ*NY*NX
#endif

  type(timer) :: solver_time
  type(field_quark_wg), allocatable :: w,s,b
  type(quark_wilson) :: mw
  type(cg_alg), allocatable :: solver
  character(CHARLEN) :: str
  integer    :: ieo,ioe, mult_count
  type(counter) :: imult
  real(DP) :: M,etime,rtmp,kappa
  real(DP) :: a5,m0

  call new(solver_time)
  call tic(solver_time)

  a5 = this%a5
  m0 = this%m0
  imult = mult_iter

  allocate(w,s)
  call new(w)
  call new(s)

#ifdef _EOPREC_LOWMODE_SOLVER_
  ieo = 0
  ioe = 1-ieo

  M = m0
  rtmp  = 1.0_DP/(4.0_DP-M)
  kappa = rtmp/2.0d0
  call set_kappa(mw,kappa)

  call assign_mult(s,v,rtmp)   ! s := 1/(4-M) v
  call accum_mult_gamma5(s)    ! s := gamma5 s
  call assign_mult_hop(mw,w%eo(ieo),s%eo(ioe),u)  ! we = Meo so
  call accum_add_mult(s%eo(ieo),w%eo(ieo),kappa)  ! se = se + kappa * we = se + kappa Meo so

  !      v
  ! we = Dee^dag se
  !
  call accum_mult_gamma5(s%eo(ieo))
  call assign_mult_eoprec_wd(mw,w%eo(ieo),s%eo(ieo),u)
  call accum_mult_gamma5(w%eo(ieo))

#else
  !-----------------------------------------
  !  w  = gamma_5 Dw v = Dw^dag gamma5 v
  !-----------------------------------------
  M = m0
  call assign(s,v)
  call assign_mult_wd(M,w,s,u)
  call  accum_mult_gamma5(w)
#endif

  !=============================================
  ! Solve:
  !
  !   (Dw^dag Dw) x = w
  !
  ! for x
  !
  !=============================================
  allocate(solver)
  call new(solver,NSIZE=EIGEN_NSITE,mode=MODE_NORMAL,max_iter=iiter,tol=tol)

#ifdef _EOPREC_LOWMODE_SOLVER_
  call pack(w%eo(ieo),solver%src_vec(:))
#else
  call pack(w,solver%src_vec(:))
#endif

  !=======================
  ! start CG Solver
  !=======================
  do 
    call solve(solver)
    select case(get_status(solver))
    case (OP_NOP)
      cycle
    case (OP_DO_MATVEC)

      !---------------------------------------------
      ! multiply gamma_5 Dw gamma_5 Dw
      !---------------------------------------------
#ifdef _EOPREC_LOWMODE_SOLVER_
      M = m0
      call unpack(solver%src_vec,w%eo(ieo))
      call assign_mult_eoprec_wd(mw,s%eo(ieo),w%eo(ieo),u)   ! s := D w
      call  accum_mult_gamma5(s%eo(ieo))                     ! s := gamma_5 s
      call assign_mult_eoprec_wd(mw,w%eo(ieo),s%eo(ieo),u)   ! w := D s
      call  accum_mult_gamma5(w%eo(ieo))                     ! w := gamma_5 w
      call pack(w%eo(ieo),solver%dst_vec)
#else
      M = m0
      call unpack(solver%src_vec,w)
      call assign_mult_wd(M,s,w,u)    ! s := D w
      call  accum_mult_gamma5(s)      ! s := gamma_5 s
      call assign_mult_wd(M,w,s,u)    ! w := D s
      call  accum_mult_gamma5(w)      ! w := gamma_5 w
      call pack(w,solver%dst_vec)
#endif

    case (OP_PRINT_STATUS)
      if (this%debug) then
      if (nodeid == 0) then
        mult_count = get_count(mult_iter) - get_count(imult)
        write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)') get_current_iteration(solver), get_residual_norm(solver),mult_count
        write(*,'("#",I5," ERR_=",E24.16," MULT_=",I10)') get_current_iteration(solver), get_residual_norm(solver),mult_count
      endif
      endif
      cycle
    case (OP_CONVERGED)
      if (nodeid == 0) then
        mult_count = get_count(mult_iter) - get_count(imult)
        write(iout,'("#",A," iteration converged.")')TRIM(get_name(solver))
        write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), get_residual_norm(solver),mult_count
      endif
      exit      
    case (OP_MAXITER_REACHED)
      if (nodeid == 0) then
        mult_count = get_count(mult_iter) - get_count(imult)
        write(iout,'("#",A," iteration does not converge.")')TRIM(get_name(solver))
        write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), get_residual_norm(solver),mult_count
      endif
      write(str,'("solver did not converge.: ",A,I4)')__FILE__,__LINE__
      call error_stop(str)
    case default
      write(str,'("something is wrong in solver.: ",A,I4)')__FILE__,__LINE__
      call error_stop(str)
    end select
  enddo

#ifdef _EOPREC_LOWMODE_SOLVER_
  call unpack(solver%dst_vec,w%eo(ieo))
  call assign_mult_hop(mw,w%eo(ioe),w%eo(ieo),u)    ! wo = Moe we
  call  accum_mult_add(w%eo(ioe),kappa,s%eo(ioe))   ! wo = kappa * wo + so = kappa Moe we + so
#else
  call unpack(solver%dst_vec,w)
#endif

  if (this%debug) then
  M = m0
  call assign_mult_wd(M,s,w,u)   ! s := Dw w
  call accum_mult_gamma5(s)      ! s := gamma5 s
  allocate(b)
  call new(b)
  call assign_sub(b,v,s)         ! b = v - s
  rtmp = abs2(b)/abs2(v)
  rtmp = sqrt(rtmp)
  deallocate(b)
  if (nodeid == 0) then
      mult_count = get_count(mult_iter) - get_count(imult)
      write(*,*)
      write(*   ,'("%",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), rtmp,mult_count
      write(iout,*)
      write(iout,'("%",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), rtmp,mult_count
  endif
  endif

  !--------------------------------
  ! set solution on Hv
  !--------------------------------
  if (abs(a5) > EPSILON(a5)) then

    !
    ! Hv = (a5/2) (Dw + 2/a5) (gamma_5 Dw) \ v
    !
    M = m0 - 2.0_DP/a5
    call assign_mult_wd(M,Hv,w,u)
    call accum_mult(Hv,a5/2.0_DP)

  else

    !
    ! Hv = w
    !
    call assign(Hv,w)

  endif

  call delete(solver)
  deallocate(solver)
  deallocate(w,s)

  call toc(solver_time)

  !============================================
  ! NOTE: get_elapse contains comlib_sumcast.
  ! This means that get_elapse should not used
  ! in if-(nodeid==0)-then-endif block.
  !============================================
  etime = get_elapse(solver_time)  
  if (nodeid == 0) then
    if (this%debug) write(*   ,'(" ETIME_=",E24.16)')etime

    write(iout,'(" ETIME_=",E24.16)')etime
  endif

  return
end subroutine
Subroutine :
this :type(quark_ovlp_kern_eigmodes), intent(in)
iout :integer, intent(in)
: solver logfile id
tol :real(DP), intent(in)
: solver tolerance
iiter :integer, intent(inout)
: solver iteration count
Hv :type(field_quark_wg), intent(inout)
v :type(field_quark_wg), intent(in)
u :type(vfield_gluon_wg), intent(in)

Multiply and assign sign kernel for overlap operator

Dw(-m0) = 4 - m0 - 1/2 Hoping

Hv = gamma_5 Dw / (Dw a5/2 + 1) v = Hw v

Hw : the kernel of sign function for overlap operator

[Source]

subroutine assign_mult_sign_kernel(this,iout,tol,iiter,Hv,v,u)
!
! Multiply and assign sign kernel for overlap operator
!
! Dw(-m0) = 4 - m0 - 1/2 Hoping
!
! Hv = gamma_5 Dw / (Dw a5/2 + 1) v = Hw v
!
! Hw : the kernel of sign function for overlap operator
!
  use solver_class
  use counter_class
  implicit none
  type(quark_ovlp_kern_eigmodes), intent(in)    :: this
  integer,                        intent(in)    :: iout    ! solver logfile id
  integer,                        intent(inout) :: iiter   ! solver iteration count
  real(DP),                       intent(in)    :: tol     ! solver tolerance
  type(field_quark_wg),    intent(inout) :: Hv
  type(field_quark_wg),    intent(in)    :: v
  type(vfield_gluon_wg),   intent(in)    :: u

  type(timer) :: solver_time
  type(field_quark_wg), allocatable :: w,s
  type(quark_wilson) :: mass_param
  type(cg_alg), allocatable :: solver
  character(CHARLEN) :: str
  type(counter) :: imult
  integer  :: mult_count
  real(DP) :: M,etime
  real(DP) :: a5,m0

  call new(solver_time)
  call tic(solver_time)
 
  a5 = this%a5
  m0 = this%m0
  imult = mult_iter

  if (abs(a5) > EPSILON(a5)) then

    allocate(w,s)
    call new(w)
    call new(s)

    !-------------------------------------
    !  w  = (D a5/2 +1) \ v
    !-------------------------------------

    !------------------------------
    !  s = (D a5/2+1)^dag v
    !------------------------------
    call assign_mult_gamma5(s,v)

    M = m0 - 2.0_DP/a5
    call assign_mult_wd(M,w,s,u)
    call accum_mult(w,a5/2.0_DP)

    call accum_mult_gamma5(w)


    !=============================================
    ! Solve:
    !   (D a5/2 + 1)^dag (D a5/2 + 1) x = w
    ! for x
    !=============================================
    allocate(solver)
    call new(solver,NSIZE=NSITE,mode=MODE_NORMAL,max_iter=iiter,tol=tol)
    call pack(w,solver%src_vec(:))

    !=======================
    ! start CG Solver
    !=======================
    do 
      call solve(solver)
      select case(get_status(solver))
      case (OP_NOP)
        cycle
      case (OP_DO_MATVEC)

        !---------------------------------------------
        ! multiply  (D a5/2 + 1)^dag (D a5/2 + 1)
        !---------------------------------------------
        call unpack(solver%src_vec,w)

        M = m0 - 2.0_DP/a5
        call assign_mult_wd(M,s,w,u)   
        call accum_mult(s,a5/2.0_DP)      ! s := (D 2/a5+1) w

        call accum_mult_gamma5(s)         ! s := gamma_5 s
        call assign_mult_wd(M,w,s,u)      ! w := (D + 2/a5) s
        call accum_mult(w,a5/2.0_DP)      ! w := (D a5/2+1) w
        call accum_mult_gamma5(w)         ! w := gamma_5 w

        call pack(w,solver%dst_vec)

      case (OP_PRINT_STATUS)
        if (this%debug) then
        if (nodeid == 0) then
          mult_count = get_count(mult_iter) - get_count(imult)
          write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)') get_current_iteration(solver), get_residual_norm(solver),mult_count
          write(*,'("#",I5," ERR_=",E24.16," MULT_=",I10)') get_current_iteration(solver), get_residual_norm(solver),mult_count
        endif
        endif
        cycle
      case (OP_CONVERGED)
        if (nodeid == 0) then
          mult_count = get_count(mult_iter) - get_count(imult)
          write(iout,'("#",A," iteration converged.")')TRIM(get_name(solver))
          write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), get_residual_norm(solver),mult_count
        endif
        exit      
      case (OP_MAXITER_REACHED)
        if (nodeid == 0) then
          mult_count = get_count(mult_iter) - get_count(imult)
          write(iout,'("#",A," iteration does not converge.")')TRIM(get_name(solver))
          write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), get_residual_norm(solver),mult_count
        endif
        write(str,'("solver did not converge.: ",A,I4)')__FILE__,__LINE__
        call error_stop(str)
      case default
        write(str,'("something is wrong in solver.: ",A,I4)')__FILE__,__LINE__
        call error_stop(str)
      end select

    enddo

    call unpack(solver%dst_vec,w)
    call delete(solver)
    deallocate(solver)

    call toc(solver_time)

    !============================================
    ! NOTE: get_elapse contains comlib_sumcast.
    ! This means that get_elapse should not used
    ! in if-(nodeid==0)-then-endif block.
    !============================================
    etime = get_elapse(solver_time)  
    if (nodeid == 0) then
      write(iout,'(" ETIME_=",E24.16)')etime
    endif

    !
    !  Hv = gamma_5 D w
    !
    M = m0
    call assign_mult_wd(M,Hv,w,u)
    call accum_mult_gamma5(Hv)

    deallocate(w,s)
  else
    allocate(w)
    call new(w)
    !
    !  Hv = gamma_5 D v
    !
    call assign(w,v)
    M = m0
    call assign_mult_wd(M,Hv,w,u)
    call accum_mult_gamma5(Hv)
    deallocate(w)
  endif

  return
end subroutine
Subroutine :
this :type(quark_ovlp_kern_eigmodes), intent(in)
iout :integer, intent(in)
: solver logfile id
tol :real(DP), intent(in)
: solver tolerance
iiter :integer, intent(inout)
: solver iteration count
Hv :type(field_quark_wg), intent(inout)
v :type(field_quark_wg), intent(in)
u :type(vfield_gluon_wg), intent(in)

Multiply and assign sign kernel for overlap operator

Dw(-m0) = 4 - m0 - 1/2 Hoping

Hv = gamma_5 Dw / (Dw a5/2 + 1) v = Hw v

Hw : the kernel of sign function for overlap operator

[Source]

subroutine assign_mult_sign_kernel(this,iout,tol,iiter,Hv,v,u)
!
! Multiply and assign sign kernel for overlap operator
!
! Dw(-m0) = 4 - m0 - 1/2 Hoping
!
! Hv = gamma_5 Dw / (Dw a5/2 + 1) v = Hw v
!
! Hw : the kernel of sign function for overlap operator
!
  use solver_class
  use counter_class
  implicit none
  type(quark_ovlp_kern_eigmodes), intent(in)    :: this
  integer,                        intent(in)    :: iout    ! solver logfile id
  integer,                        intent(inout) :: iiter   ! solver iteration count
  real(DP),                       intent(in)    :: tol     ! solver tolerance
  type(field_quark_wg),    intent(inout) :: Hv
  type(field_quark_wg),    intent(in)    :: v
  type(vfield_gluon_wg),   intent(in)    :: u

  type(timer) :: solver_time
  type(field_quark_wg), allocatable :: w,s
  type(quark_wilson) :: mass_param
  type(cg_alg), allocatable :: solver
  character(CHARLEN) :: str
  type(counter) :: imult
  integer  :: mult_count
  real(DP) :: M,etime
  real(DP) :: a5,m0

  call new(solver_time)
  call tic(solver_time)
 
  a5 = this%a5
  m0 = this%m0
  imult = mult_iter

  if (abs(a5) > EPSILON(a5)) then

    allocate(w,s)
    call new(w)
    call new(s)

    !-------------------------------------
    !  w  = (D a5/2 +1) \ v
    !-------------------------------------

    !------------------------------
    !  s = (D a5/2+1)^dag v
    !------------------------------
    call assign_mult_gamma5(s,v)

    M = m0 - 2.0_DP/a5
    call assign_mult_wd(M,w,s,u)
    call accum_mult(w,a5/2.0_DP)

    call accum_mult_gamma5(w)


    !=============================================
    ! Solve:
    !   (D a5/2 + 1)^dag (D a5/2 + 1) x = w
    ! for x
    !=============================================
    allocate(solver)
    call new(solver,NSIZE=NSITE,mode=MODE_NORMAL,max_iter=iiter,tol=tol)
    call pack(w,solver%src_vec(:))

    !=======================
    ! start CG Solver
    !=======================
    do 
      call solve(solver)
      select case(get_status(solver))
      case (OP_NOP)
        cycle
      case (OP_DO_MATVEC)

        !---------------------------------------------
        ! multiply  (D a5/2 + 1)^dag (D a5/2 + 1)
        !---------------------------------------------
        call unpack(solver%src_vec,w)

        M = m0 - 2.0_DP/a5
        call assign_mult_wd(M,s,w,u)   
        call accum_mult(s,a5/2.0_DP)      ! s := (D 2/a5+1) w

        call accum_mult_gamma5(s)         ! s := gamma_5 s
        call assign_mult_wd(M,w,s,u)      ! w := (D + 2/a5) s
        call accum_mult(w,a5/2.0_DP)      ! w := (D a5/2+1) w
        call accum_mult_gamma5(w)         ! w := gamma_5 w

        call pack(w,solver%dst_vec)

      case (OP_PRINT_STATUS)
        if (this%debug) then
        if (nodeid == 0) then
          mult_count = get_count(mult_iter) - get_count(imult)
          write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)') get_current_iteration(solver), get_residual_norm(solver),mult_count
          write(*,'("#",I5," ERR_=",E24.16," MULT_=",I10)') get_current_iteration(solver), get_residual_norm(solver),mult_count
        endif
        endif
        cycle
      case (OP_CONVERGED)
        if (nodeid == 0) then
          mult_count = get_count(mult_iter) - get_count(imult)
          write(iout,'("#",A," iteration converged.")')TRIM(get_name(solver))
          write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), get_residual_norm(solver),mult_count
        endif
        exit      
      case (OP_MAXITER_REACHED)
        if (nodeid == 0) then
          mult_count = get_count(mult_iter) - get_count(imult)
          write(iout,'("#",A," iteration does not converge.")')TRIM(get_name(solver))
          write(iout,'("#",I5," ERR_=",E24.16," MULT_=",I10)',advance='no') get_current_iteration(solver), get_residual_norm(solver),mult_count
        endif
        write(str,'("solver did not converge.: ",A,I4)')__FILE__,__LINE__
        call error_stop(str)
      case default
        write(str,'("something is wrong in solver.: ",A,I4)')__FILE__,__LINE__
        call error_stop(str)
      end select

    enddo

    call unpack(solver%dst_vec,w)
    call delete(solver)
    deallocate(solver)

    call toc(solver_time)

    !============================================
    ! NOTE: get_elapse contains comlib_sumcast.
    ! This means that get_elapse should not used
    ! in if-(nodeid==0)-then-endif block.
    !============================================
    etime = get_elapse(solver_time)  
    if (nodeid == 0) then
      write(iout,'(" ETIME_=",E24.16)')etime
    endif

    !
    !  Hv = gamma_5 D w
    !
    M = m0
    call assign_mult_wd(M,Hv,w,u)
    call accum_mult_gamma5(Hv)

    deallocate(w,s)
  else
    allocate(w)
    call new(w)
    !
    !  Hv = gamma_5 D v
    !
    call assign(w,v)
    M = m0
    call assign_mult_wd(M,Hv,w,u)
    call accum_mult_gamma5(Hv)
    deallocate(w)
  endif

  return
end subroutine
Subroutine :
mass :real(DP), intent(in)
: supply a positive value for negative mass operator
Dv :type(field_quark_wg), intent(inout)
v :type(field_quark_wg), intent(inout)
u :type(vfield_gluon_wg), intent(in)

Dv = Dw v

Dw : negative mass wilson dirac matrix

Dw = 4 - mass - 1/2 Hop

[Source]

subroutine assign_mult_wd(mass,Dv,v,u)
!
! Dv = Dw v
!
! Dw : negative mass wilson dirac matrix
!
! Dw = 4 - mass - 1/2 Hop
!
  real(DP),                intent(in)    :: mass  ! supply a positive value for negative mass operator
  type(field_quark_wg),    intent(inout) :: Dv
  type(field_quark_wg),    intent(inout) :: v
  type(vfield_gluon_wg),   intent(in)    :: u
  type(field_quark_wg), allocatable :: w
  type(quark_wilson) :: mass_param  ! dummy
  real(DP) :: rtmp

  allocate(w)
  call new(w)

  rtmp = 4.0_DP - mass
  call assign_mult(w,v,rtmp)                !  w := v *(4 - mass)
  call assign_mult_hop(mass_param,Dv,v,u)   ! Dv := Hop v 
  rtmp = -0.5_DP
  call accum_mult_add(Dv,rtmp,w)            ! Dv := Dv * (-0.5) + w = (4 - mass - 0.5 Hop) v

  deallocate(w)

  return
end subroutine
Subroutine :
mass :real(DP), intent(in)
: supply a positive value for negative mass operator
Dv :type(field_quark_wg), intent(inout)
v :type(field_quark_wg), intent(inout)
u :type(vfield_gluon_wg), intent(in)

Dv = Dw v

Dw : negative mass wilson dirac matrix

Dw = 4 - mass - 1/2 Hop

[Source]

subroutine assign_mult_wd(mass,Dv,v,u)
!
! Dv = Dw v
!
! Dw : negative mass wilson dirac matrix
!
! Dw = 4 - mass - 1/2 Hop
!
  real(DP),                intent(in)    :: mass  ! supply a positive value for negative mass operator
  type(field_quark_wg),    intent(inout) :: Dv
  type(field_quark_wg),    intent(inout) :: v
  type(vfield_gluon_wg),   intent(in)    :: u
  type(field_quark_wg), allocatable :: w
  type(quark_wilson) :: mass_param  ! dummy
  real(DP) :: rtmp

  allocate(w)
  call new(w)

  rtmp = 4.0_DP - mass
  call assign_mult(w,v,rtmp)                !  w := v *(4 - mass)
  call assign_mult_hop(mass_param,Dv,v,u)   ! Dv := Hop v 
  rtmp = -0.5_DP
  call accum_mult_add(Dv,rtmp,w)            ! Dv := Dv * (-0.5) + w = (4 - mass - 0.5 Hop) v

  deallocate(w)

  return
end subroutine
Subroutine :
this :type(quark_ovlp_kern_eigmodes), intent(inout)

[Source]

subroutine delete_lowmode(this)
  implicit none
  type(quark_ovlp_kern_eigmodes), intent(inout) :: this
  if (allocated(this%V))  deallocate(this%V)
  if (allocated(this%EV)) deallocate(this%EV)
  this%is_initialized = .FALSE.
  return
end subroutine
Subroutine :
this :type(quark_ovlp_kern_eigmodes), intent(inout)

[Source]

subroutine delete_lowmode(this)
  implicit none
  type(quark_ovlp_kern_eigmodes), intent(inout) :: this
  if (allocated(this%V))  deallocate(this%V)
  if (allocated(this%EV)) deallocate(this%EV)
  this%is_initialized = .FALSE.
  return
end subroutine
Subroutine :
this :type(quark_ovlp_kern_eigmodes), intent(inout)
iout :integer, intent(in)
tol :real(DP), intent(in)
iiter :integer, intent(out)
u :type(vfield_gluon_wg), intent(in)

Compute low/high modes of

Hw = gamma_5 Dw/(Dw*a5/2+1)

[Source]

subroutine get_sign_kernel_eigen_modes(this,iout,tol,iiter,u)
!                    
! Compute low/high modes of 
!                    
! Hw  = gamma_5 Dw/(Dw*a5/2+1)
!
  use ks_class
  use f95_lapack, only : LA_HESV
  implicit none
  type(quark_ovlp_kern_eigmodes), intent(inout) :: this
  integer,  intent(in)  :: iout
  real(DP), intent(in)  :: tol
  integer,  intent(out) :: iiter
  type(vfield_gluon_wg), intent(in)    ::  u

  type(field_quark_wg),  allocatable :: w,v
  type(ks_alg),          allocatable :: eigen_solver
  type(timer) :: solver_time
  complex(DP) :: ctmp
  real(DP) :: etime, rtmp
  real(DP) :: a5,m0
  character(CHARLEN) :: str
  integer, parameter :: INV_MAXITER=10000
  integer :: iter_high
  integer :: NS,nev,i,j,jiter,nnev

  call new(solver_time)
  call tic(solver_time)

  if (.not.this%is_initialized) then
    call error_stop("get_sign_kernel_eigen_modes:type(qurark_ovlp_kern_eigmodes) is not initialized.")
  endif

  a5 = this%a5
  m0 = this%m0

  allocate(w,v)
  call new(w)
  call new(v)

  allocate(eigen_solver)

  iiter = 0
  nev = this%NEV
  write(iout,'("# NEV=",I3)')nev
  call new(eigen_solver,n=NSITE,m=nev*3,nev=nev,maxiter=400,tol=tol,debug=this%debug,mode=ESLV_MODE_HERMITE)

  do
    call solve(eigen_solver)
    select case(get_status(eigen_solver))
    case (ESLV_OP_NOP)
      cycle
    case (ESLV_OP_DO_MATVEC)
      call unpack(eigen_solver%src_vec,w)
      jiter = INV_MAXITER

      select case(this%which_modes)
      case(SIGN_KERNEL_LOW_MODES)

        !
        ! for low-modes, apply  1/Hw
        !
        call assign_inv_mult_sign_kernel(this,iout,tol,jiter,v,w,u)

      case(SIGN_KERNEL_HIGH_MODES)

        !
        ! for high-modes, apply  Hw
        !
        call     assign_mult_sign_kernel(this,iout,tol,jiter,v,w,u)

      end select

      iiter = iiter + jiter
      call pack(v,eigen_solver%dst_vec)
    case (ESLV_OP_CONVERGED)
      exit
    case (ESLV_OP_FAIL)
      write(str,'("Error in eigen solver. ITER:",I5,1X,A,I5)') get_current_iteration(eigen_solver),__FILE__,__LINE__
      call error_stop(TRIM(str))
    end select
  enddo

  nnev = get_num_converged(eigen_solver)

  if (0 == nodeid) write(iout,'("# NNEV=",I3)')nnev

  !--------------------
  ! store eigen modes
  !--------------------
  do i=1,NEV
    call unpack(eigen_solver%V(:,i),this%V(i))
  enddo

  select case(this%which_modes)
  case(SIGN_KERNEL_LOW_MODES)
    do i=1,NEV
      this%EV(i) = REAL(1.0_DP/eigen_solver%eval(i),kind=DP)
    enddo
  case(SIGN_KERNEL_HIGH_MODES)
    do i=1,NEV
      this%EV(i) = REAL(eigen_solver%eval(i),kind=DP)
    enddo
  end select

  call delete(eigen_solver)
  deallocate(v,w)
  deallocate(eigen_solver)

  if (this%debug) then
    do i=1,NEV
    do j=i,NEV
      ctmp = prod(this%V(i),this%V(j))
      if (abs(ctmp) > 1.d-10) then
        if (0 == nodeid) write(iout,'("%EORTH",2I3,2E24.15)')i,j,ctmp
      endif
    enddo
    enddo
  endif

  if (0 == nodeid) then
    do i=1,NEV
      write(iout,'("%EVAL",I3,E24.15)')i,this%EV(i)
    end do
  endif

  call toc(solver_time)

  etime = get_elapse(solver_time)  

  if (0 == nodeid) then
    if (this%debug) write(   *,'("# ",A,"_ETIME_=",E24.16)')TRIM(ADJUSTL(MODE_STR(this%which_modes))),etime
    write(iout,'("# ",A,"_ETIME_=",E24.16)')TRIM(ADJUSTL(MODE_STR(this%which_modes))),etime
  endif

  return
end subroutine
Subroutine :
this :type(quark_ovlp_kern_eigmodes), intent(inout)
iout :integer, intent(in)
tol :real(DP), intent(in)
iiter :integer, intent(out)
u :type(vfield_gluon_wg), intent(in)

Compute low/high modes of

Hw = gamma_5 Dw/(Dw*a5/2+1)

[Source]

subroutine get_sign_kernel_eigen_modes(this,iout,tol,iiter,u)
!                    
! Compute low/high modes of 
!                    
! Hw  = gamma_5 Dw/(Dw*a5/2+1)
!
  use ks_class
  use f95_lapack, only : LA_HESV
  implicit none
  type(quark_ovlp_kern_eigmodes), intent(inout) :: this
  integer,  intent(in)  :: iout
  real(DP), intent(in)  :: tol
  integer,  intent(out) :: iiter
  type(vfield_gluon_wg), intent(in)    ::  u

  type(field_quark_wg),  allocatable :: w,v
  type(ks_alg),          allocatable :: eigen_solver
  type(timer) :: solver_time
  complex(DP) :: ctmp
  real(DP) :: etime, rtmp
  real(DP) :: a5,m0
  character(CHARLEN) :: str
  integer, parameter :: INV_MAXITER=10000
  integer :: iter_high
  integer :: NS,nev,i,j,jiter,nnev

  call new(solver_time)
  call tic(solver_time)

  if (.not.this%is_initialized) then
    call error_stop("get_sign_kernel_eigen_modes:type(qurark_ovlp_kern_eigmodes) is not initialized.")
  endif

  a5 = this%a5
  m0 = this%m0

  allocate(w,v)
  call new(w)
  call new(v)

  allocate(eigen_solver)

  iiter = 0
  nev = this%NEV
  write(iout,'("# NEV=",I3)')nev
  call new(eigen_solver,n=NSITE,m=nev*3,nev=nev,maxiter=400,tol=tol,debug=this%debug,mode=ESLV_MODE_HERMITE)

  do
    call solve(eigen_solver)
    select case(get_status(eigen_solver))
    case (ESLV_OP_NOP)
      cycle
    case (ESLV_OP_DO_MATVEC)
      call unpack(eigen_solver%src_vec,w)
      jiter = INV_MAXITER

      select case(this%which_modes)
      case(SIGN_KERNEL_LOW_MODES)

        !
        ! for low-modes, apply  1/Hw
        !
        call assign_inv_mult_sign_kernel(this,iout,tol,jiter,v,w,u)

      case(SIGN_KERNEL_HIGH_MODES)

        !
        ! for high-modes, apply  Hw
        !
        call     assign_mult_sign_kernel(this,iout,tol,jiter,v,w,u)

      end select

      iiter = iiter + jiter
      call pack(v,eigen_solver%dst_vec)
    case (ESLV_OP_CONVERGED)
      exit
    case (ESLV_OP_FAIL)
      write(str,'("Error in eigen solver. ITER:",I5,1X,A,I5)') get_current_iteration(eigen_solver),__FILE__,__LINE__
      call error_stop(TRIM(str))
    end select
  enddo

  nnev = get_num_converged(eigen_solver)

  if (0 == nodeid) write(iout,'("# NNEV=",I3)')nnev

  !--------------------
  ! store eigen modes
  !--------------------
  do i=1,NEV
    call unpack(eigen_solver%V(:,i),this%V(i))
  enddo

  select case(this%which_modes)
  case(SIGN_KERNEL_LOW_MODES)
    do i=1,NEV
      this%EV(i) = REAL(1.0_DP/eigen_solver%eval(i),kind=DP)
    enddo
  case(SIGN_KERNEL_HIGH_MODES)
    do i=1,NEV
      this%EV(i) = REAL(eigen_solver%eval(i),kind=DP)
    enddo
  end select

  call delete(eigen_solver)
  deallocate(v,w)
  deallocate(eigen_solver)

  if (this%debug) then
    do i=1,NEV
    do j=i,NEV
      ctmp = prod(this%V(i),this%V(j))
      if (abs(ctmp) > 1.d-10) then
        if (0 == nodeid) write(iout,'("%EORTH",2I3,2E24.15)')i,j,ctmp
      endif
    enddo
    enddo
  endif

  if (0 == nodeid) then
    do i=1,NEV
      write(iout,'("%EVAL",I3,E24.15)')i,this%EV(i)
    end do
  endif

  call toc(solver_time)

  etime = get_elapse(solver_time)  

  if (0 == nodeid) then
    if (this%debug) write(   *,'("# ",A,"_ETIME_=",E24.16)')TRIM(ADJUSTL(MODE_STR(this%which_modes))),etime
    write(iout,'("# ",A,"_ETIME_=",E24.16)')TRIM(ADJUSTL(MODE_STR(this%which_modes))),etime
  endif

  return
end subroutine
Subroutine :
this :type(quark_ovlp_kern_eigmodes), intent(inout)

[Source]

subroutine new_lowmode(this)
  implicit none
  type(quark_ovlp_kern_eigmodes), intent(inout) :: this
  integer :: i

  if (this%NEV == 0) then
    call error_stop("Set subspace size NEV tepe(quark_ovlp_kern_eigmodes).")
  endif

  if (allocated(this%V))  deallocate(this%V)
  if (allocated(this%EV)) deallocate(this%EV)

  allocate(this%V(1:this%NEV))
  do i=1,this%NEV
    call new(this%V(i))
  enddo
  allocate(this%EV(this%NEV))

  this%EV(:) = 0.0_DP
  this%m0 = 1.0_DP
  this%a5 = 1.0_DP
  this%which_modes = SIGN_KERNEL_LOW_MODES
  this%is_initialized = .TRUE.

  return
end subroutine
Subroutine :
this :type(quark_ovlp_kern_eigmodes), intent(inout)

[Source]

subroutine new_lowmode(this)
  implicit none
  type(quark_ovlp_kern_eigmodes), intent(inout) :: this
  integer :: i

  if (this%NEV == 0) then
    call error_stop("Set subspace size NEV tepe(quark_ovlp_kern_eigmodes).")
  endif

  if (allocated(this%V))  deallocate(this%V)
  if (allocated(this%EV)) deallocate(this%EV)

  allocate(this%V(1:this%NEV))
  do i=1,this%NEV
    call new(this%V(i))
  enddo
  allocate(this%EV(this%NEV))

  this%EV(:) = 0.0_DP
  this%m0 = 1.0_DP
  this%a5 = 1.0_DP
  this%which_modes = SIGN_KERNEL_LOW_MODES
  this%is_initialized = .TRUE.

  return
end subroutine
quark_ovlp_kern_eigmodes
Derived Type :
m0 :real(DP)
: kernel negaitive mass paremter, (m0 > 0).
a5 :real(DP)
: kernel mebius parameter
V(:) :type(field_quark_wg), allocatable
: Hw V = V lambda, V : eigen vectors
EV(:) :real(DP), allocatable
: lambda : eigen values
NEV = 0 :integer
: subspace dimension
which_modes = SIGN_KERNEL_LOW_MODES :integer
: high mode or low mode
debug = .FALSE. :logical
is_initialized = .FALSE. :logical

eigen mode container

quark_ovlp_kern_eigmodes
Derived Type :
m0 :real(DP)
: kernel negaitive mass paremter, (m0 > 0).
a5 :real(DP)
: kernel mebius parameter
V(:) :type(field_quark_wg), allocatable
: Hw V = V lambda, V : eigen vectors
EV(:) :real(DP), allocatable
: lambda : eigen values
NEV = 0 :integer
: subspace dimension
which_modes = SIGN_KERNEL_LOW_MODES :integer
: high mode or low mode
debug = .FALSE. :logical
is_initialized = .FALSE. :logical

eigen mode container