Class hmc_qcd_class
In: HmcQcdClass/hmc_qcd_class.F90
lattice_class error_class timer_class counter_class hmc_logfile_class hmc_gluon_class hmc_quark_class hmc_status_class comlib print_status_class file_tools_class random_class hmc_identity sf_measure_class field_gauge_class hmc_qcd_class dot/f_0.png

HMC algorithm with QCD action

Version

 $Id: hmc_qcd_class.F90,v 1.32 2011/12/01 12:48:46 ishikawa Exp $

Methods

Included Modules

lattice_class error_class timer_class counter_class hmc_logfile_class hmc_gluon_class hmc_quark_class hmc_status_class comlib print_status_class file_tools_class random_class hmc_identity sf_measure_class field_gauge_class

Public Instance methods

Subroutine :
this :type(qcd_hmc), intent(inout)

Delelete HMC QCD action

[Source]

subroutine delete_qcd_hmc(this)
!
! Delelete HMC QCD action
!
  implicit none
  type(qcd_hmc), intent(inout) :: this

  call delete(this%gluon)
  if (SW_ON == this%switch_dynamical_quark) then
    call delete(this%quark)
    deallocate(this%quark)
  endif
  call delete(this%plaq_log)
  call delete(this%force_log)

  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
: HMC QCD action
ifi :integer, intent(in)
: select target
0
for initial Hamiltonian with generating P.
1
for final Hamiltonian without generating P.
2
for initial Hamiltonian without generating P.

Compute Total Hamiltonian for HMC QCD action

[Source]

subroutine hamil_qcd_hmc(this,ifi)
!
! Compute Total Hamiltonian for HMC QCD action
! 
  implicit none
  type(qcd_hmc), intent(inout) :: this  ! HMC QCD action
  integer,       intent(in)    :: ifi ! select target
                                      ! [0] for initial Hamiltonian with generating P.
                                      ! [1] for final Hamiltonian without generating P.
                                      ! [2] for initial Hamiltonian without generating P.

  type(tfield_gluon_wg), allocatable :: wl34,wl98
  real(DP) :: H
  real(DP) :: SU,SP
  real(DP) :: SQ_all
  real(DP) :: SQ(SIZE(this%quark))
  character(CHARLEN) :: str

  H      = 0.0_DP
  SU     = 0.0_DP
  SP     = 0.0_DP
  SQ_all = 0.0_DP
  SQ(:)  = 0.0_DP

  if (0 == ifi) then
    this%mult0 = get_count(mult_iter)
    if (.not. associated(this%p) ) then
      allocate(this%p)
      call new(this%p)
    endif
  endif

  !========================================
  ! calc. Hamiltonian for gluon action
  ! if (ifi==0)  p is updated
  !========================================
  allocate(wl34,wl98)
  call new(wl34)
  call new(wl98)
  call make_two_links(this%u,wl34,wl98)
  call hamil(this%gluon,ifi,this%u,this%p,SU,SP,wl34,wl98)

#ifndef _CLOVER
  deallocate(wl34)
  deallocate(wl98)
#endif

  !========================================
  ! calc. Hamiltonian for quark action
  !========================================
  if (SW_ON == this%switch_dynamical_quark) then
#ifdef _CLOVER
    call set_two_links(this%quark,wl34,wl98)
#endif
    call hamil(this%quark,ifi,this%u,SQ)
  endif

#ifdef _CLOVER
  call remove_two_links(this%quark)
  deallocate(wl34)
  deallocate(wl98)
#endif

  SQ_all=sum(SQ)
  H=SU+SP+SQ_all
!  Hamil%H=H
!  Hamil%SU=SU
!  Hamil%SP=SP
!  Hamil%SQ=SQ_all

  select case(ifi)
  case (0) 
    this%H0 = H
    this%WLP0 = this%gluon%WLP
  case (1)
    this%H1 = H
    this%WLP1 = this%gluon%WLP
  case (2)
    this%H2 = H
    this%WLP2 = this%gluon%WLP
  end select

  if (nodeid == 0) then
    write(str,'("@",I8,I2,"   H:",E24.16,"  SU:",E24.16,"  SP:",E24.16,"  SQ:",E24.16)') get_trajectory_number(this%hmc_status), ifi,H,SU,SP,SQ_all
    call print_log_action(this%hmc_status,str)
    write(str,*)
    call print_log_action(this%hmc_status,str)
  endif

  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
iflg :integer, intent(in)
fname :character(CHARLEN), intent(in)

Initialize/load field configuration according to HMC status

[Source]

subroutine initialize_config_qcd_hmc(this,iflg,fname)
!
! Initialize/load field configuration according to HMC status
!
  use print_status_class
  use random_class
  implicit none
  type(qcd_hmc),   intent(inout) :: this
  integer,            intent(in) :: iflg
  character(CHARLEN), intent(in) :: fname


  select case (iflg)
  case (0)   ! cold start with new random number
    call set_identity(this%u)
    this%u0 = this%u
  case default
    call read_cont_config(this,fname)
    this%u0 = this%u
  end select

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

Do measurment on the fly during HMC algorithm.

This is called within the MD integrator.

ifdef _SF

[Source]

subroutine measurement_onfly_qcd_hmc(this)
!
! Do measurment on the fly during HMC algorithm.
!
! This is called within the MD integrator.
!
#ifdef _SF
  use sf_measure_class
#endif
  implicit none
  type(qcd_hmc), intent(inout) :: this
  integer(8) :: mult1
  integer    :: itraj
  character(CHARLEN) :: str
#ifdef _SF
  type(sf_measurement) :: sf_obs
#endif

  itraj = get_trajectory_number(this%hmc_status)

  !------------------------
  ! save plaquette history
  !------------------------
  if (nodeid == 0) then
    mult1 = get_count(mult_iter) - this%mult0
    write(str,'(I6,10E24.16)') #ifdef _RGGAUGE #else #endif
    call print(this%plaq_log,TRIM(str))
  endif

  !-------------------
  ! print status
  !-------------------
  if (nodeid == 0) then
    write(*,'("# ",I8,"  PLQ: ",F19.16,2I12)') itraj,this%WLP1%PLQ,get_count(mult_iter),mult1
  endif

  !-------------------------------
  ! compute reweighting factor
  !-------------------------------
  if (SW_ON == this%switch_dynamical_quark) then
    call reweighting(this%quark,this%u)
  endif

  !-------------------------------
  ! measurements
  !-------------------------------
  if (SW_ON == this%switch_measure_meson) then
#ifdef _SF
    call new(sf_obs,"Params_SF_PCAC",this%hmc_status)
    call measure(sf_obs,this%gluon%gluon,this%u)
    call delete(sf_obs)
#endif
  endif

  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
: HMC QCD action
ihit :integer, intent(out)
: 1 for accepted, 0 for rejected.

Do all Metropolis tests for HMC QCD action

This routine contains HMC Metropolis test and other Globarl Metropolis tests. This routine update configuration according to the result of the Metropolis test.

If accepted the trial configuration U is copied to U0, If rejceed the original configuration U0 is copied to U.

 After calling this Metropolis test, U is the new configuration.

[Source]

subroutine metropolis_all_qcd_hmc(this,ihit)
!
! Do all Metropolis tests for HMC QCD action
!
! This routine contains HMC Metropolis test and other Globarl Metropolis tests.
! This routine update configuration according to the result of the Metropolis test.
!
! If accepted the trial configuration U  is copied to U0,
! If rejceed  the original configuration U0 is copied to U.
!
!  After calling this Metropolis test, U is the new configuration.
! 
  implicit none
  type(qcd_hmc), intent(inout) :: this   ! HMC QCD action
  integer,       intent(out)   :: ihit   ! 1 for accepted, 0 for rejected.
  real(DP) :: dH

  if (associated(this%p)) then
    deallocate(this%p)
  endif

  dH = this%H1 - this%H0
  call metropolis_qcd_hmc(this,dH,ihit)

  if (SW_ON == this%switch_dynamical_quark) then
    call metropolis(this%quark,this%u0,this%u,ihit)
  endif

  this%mult1 = get_count(mult_iter)

  if (ihit==1) then
    this%u0   = this%u
    this%WLP0 = this%WLP1
  else
    this%u    = this%u0
    this%WLP1 = this%WLP0
  endif

  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
status :type(hmc_status), intent(in), target

Initialize and Create HMC QCD action

[Source]

subroutine new_qcd_hmc(this,status)
!
! Initialize and Create HMC QCD action
!
  use comlib
  implicit none
  type(qcd_hmc), intent(inout) :: this
  type(hmc_status), intent(in), target :: status

  call new(m_lattice)      ! create global lattice

  allocate(this%u,this%u0)
  call new(this%u)
  call new(this%u0)

  this%hmc_status => status

  call new(this%force_log,TRIM(ADJUSTL(this%force_log_fname)))
  call new( this%plaq_log, TRIM(ADJUSTL(this%plaq_log_fname)))

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

Print out HMC QCD parameters on display

[Source]

subroutine print_qcd_hmc(this)
!
! Print out HMC QCD parameters on display
!
  implicit none
  type(qcd_hmc), intent(inout) :: this

  if (nodeid == 0) then
    write(*,'("**** Lattice ****")')
  endif
  call print(m_lattice) ! print global lattice

  if (nodeid == 0) then
    write(*,*)
    write(*,'("**** Gluon Action ****")')
  endif
  call print(this%gluon)

  select case (this%switch_dynamical_quark)
  case (SW_OFF)
    if (nodeid == 0) then
      write(*,*)
      write(*,'("**** Quenched QCD ****")')
    endif
  case (SW_ON)
    if (nodeid == 0) then
      write(*,*)
      write(*,'("**** Full QCD ****")')
    endif
    call print(this%quark)
  end select

  select case (this%switch_dynamical_quark)
  case (SW_OFF)
  case (SW_ON)
    if (nodeid == 0) then
      write(*,*)
      write(*,'("**** DO measurement at config save ****")')
    endif
  end select
  
  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
elapse_time :real(DP), intent(in)

Print out HMC QCD job statistics on display

[Source]

subroutine print_stat_qcd_hmc(this,elapse_time)
!
! Print out HMC QCD job statistics on display
!
  use comlib
  implicit none
  type(qcd_hmc), intent(inout) :: this
  real(DP), intent(in) :: elapse_time
  real(DP) :: mult_traj
  real(DP) :: copy_gluon_time
  real(DP) :: copy_quark_time

  mult_traj = REAL(get_count(mult_iter),kind=DP)/get_sweep_count(this%hmc_status)
#ifndef _singlePU
  call comlib_sumcast(mult_traj)
  mult_traj=mult_traj/NPU
#endif

  copy_gluon_time = get_elapse(copy_sfg_time)
  copy_quark_time = get_elapse(copy_fq_time)

  call print_copy_sfg_statistics()

  if (nodeid == 0) then
    write(*,'("                copy_u time (sec) :",F12.4," (",F8.4," %)")') copy_gluon_time, 100*copy_gluon_time/elapse_time
  endif

  call print_copy_fq_statistics()

  if (nodeid == 0) then
    write(*,'("             copy_y(c) time (sec) :",F12.4," (",F8.4," %)")') copy_quark_time, 100*copy_quark_time/elapse_time
    write(*,'("           # of Mult (par 1sweep) :",F12.4)') mult_traj
  endif

  if (SW_ON == this%switch_dynamical_quark) then
    if (nodeid == 0) then
      write(*,'(80("-"))')
      write(*,'("  Quark Statistics")')
    endif
    call print_statistics(this%quark)
  endif

  if (nodeid == 0) then
    write(*,'(80("-"))')
  endif

  return
end subroutine
qcd_hmc
Derived Type :

QCD HMC program parameter

Subroutine :
this :type(qcd_hmc), intent(inout)
iout :integer, intent(in)

Read HMC QCD job parameters form file indicated by unit id

[Source]

subroutine read_qcd_hmc(this,iout)
!
! Read HMC QCD job parameters form file indicated by unit id 
!
  use comlib
  use print_status_class
  implicit none
  type(qcd_hmc), intent(inout) :: this
  integer,       intent(in)    :: iout
  integer :: iqk,nqk

  allocate(this%gluon)
  call new(this%gluon,id=0)

  if (nodeid==0) then
    read(iout,'(A)')this%gluon%file_name
    read(iout,*)this%number_of_quark
  endif

#ifndef _singlePU
  call comlib_bcast(this%gluon%file_name,0)
  call comlib_bcast(this%number_of_quark,0)
#endif

  nqk=this%number_of_quark
  allocate(this%quark(1:nqk))
  call new(this%quark)

  if (nodeid==0) then
    do iqk=1,nqk
      read(iout,'(A)')this%quark(iqk)%file_name
    enddo
    read(iout,*)this%switch_dynamical_quark
    read(iout,*)this%switch_measure_meson
  endif

#ifndef _singlePU
  do iqk=1,nqk
    call comlib_bcast(this%quark(iqk)%file_name,0)
  enddo
  call comlib_bcast(this%switch_dynamical_quark,0)
  call comlib_bcast(this%switch_measure_meson,0)
#endif

  !==================================
  ! read gauge parameter
  !==================================
  call read(this%gluon)

  if (SW_ON == this%switch_dynamical_quark) then
    !==================================
    ! read quark parameter
    !==================================
    call read(this%quark)
  endif

  if (nodeid == 0) write(*,'(80("-"))')

  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
: HMC QCD action
ifi :integer, intent(in)
: switch, 0 for initial, 1 for final, 2 for rev check.
 MD Reversivility check routine

 Save, restore, and check the differneces of configurations

[Source]

subroutine rev_check_qcd_hmc(this,ifi)
!
!  MD Reversivility check routine
!
!  Save, restore, and check the differneces of configurations
!
  implicit none
  type(qcd_hmc), intent(inout) :: this   ! HMC QCD action
  integer,       intent(in)    :: ifi    ! switch, 0 for initial, 1 for final, 2 for rev check.
  type(vfield_gluon_wg),  allocatable :: duu
  type(vfield_gluon_wog), allocatable :: dpp
  real(DP) :: diff_u,diff_p

  select case(ifi)
  case(0)
  !================================
  ! Keep initial configuration
  !          (p(t=0),u(t=0))
  !================================
    if (.not. associated(m_u0)) allocate(m_u0)
    if (.not. associated(m_p0)) allocate(m_p0)
    call new(m_u0)
    call new(m_p0)
    m_u0 = this%u0
    m_p0 = this%p

  case(1)
  !================================
  ! Keep final configuration
  !         (p(t=tau),u(t=tau))
  !================================
    if (.not. associated(m_u1)) allocate(m_u1)
    if (.not. associated(m_p1)) allocate(m_p1)
    call new(m_u1)
    call new(m_p1)
    m_u1 = this%u
    m_p1 = this%p

  case(2)
  !==========================================
  ! compare inital and reversed
  !                 configuration
  ! (p(0),u(0)) vs (p(tau-tau),u(tau-tau))
  !==========================================

    !==================================================
    ! diff_u = |u(t=0) - u(t=tau-tau)|/COL/Sqrt(VOL*T)
    ! diff_p = |p(t=0) - p(t=tau-tau)|/COL/Sqrt(VOL*T)
    !==================================================
    diff_u = delta_u(m_u0,this%u)
    diff_p = delta_p(m_p0,this%p)

    if (nodeid==0) then
      write(*,'("%   DU,DP: ",2E24.16)') diff_u,diff_p
      write(*,'("% DH,DH/H: ",2E24.16)')  this%H0-this%H2 , (this%H0-this%H2)/this%H0
    endif

    !==========================
    ! restore configuration
    !         (p(tau),u(tau))
    !==========================
    this%u = m_u1
    this%p = m_p1

    if (associated(m_u0)) deallocate(m_u0)
    if (associated(m_p0)) deallocate(m_p0)
    if (associated(m_u1)) deallocate(m_u1)
    if (associated(m_p1)) deallocate(m_p1)
  end select

  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
iout :integer, intent(in)

Save HMC QCD job parameter in a file with unit id iout.

[Source]

subroutine save_qcd_hmc(this,iout)
!
! Save HMC QCD job parameter in a file with unit id iout.
!
  implicit none
  type(qcd_hmc), intent(inout) :: this
  integer,       intent(in)    :: iout

  integer :: iqk,nqk

  if (nodeid == 0) then
    write(iout,'(a)') TRIM(this%gluon%file_name)
    write(iout,'(I3)')this%number_of_quark
    nqk=this%number_of_quark
    do iqk=1,nqk
      write(iout,'(a)')TRIM(this%quark(iqk)%file_name)
    enddo
    write(iout,'(I3)')this%switch_dynamical_quark
    write(iout,'(I3)')this%switch_measure_meson
  endif

  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
iout :integer, intent(in)

Save configuration for measurment to a file indicated by an unit id.

[Source]

subroutine save_config_qcd_hmc(this,iout)
!
! Save configuration for measurment to a file indicated by an unit id.
!
  use hmc_identity
  use print_status_class
  use comlib
  implicit none
  type(qcd_hmc), intent(inout) :: this
  integer,       intent(in)    :: iout
  real(DP) :: plq

  call plaquette(this%u,plq)

!**************************************
!* write configuration
!**************************************
  write(iout) this%u
!**************************************
!* configuration information
!**************************************
  write(iout) NTX,NTY,NTZ,NTT
  write(iout) NDIMX,NDIMY,NDIMZ
  write(iout) plq
  call save_config(this%gluon,iout)
  if (SW_ON == this%switch_dynamical_quark) then
    write(iout) this%number_of_quark
    call save_config(this%quark,iout)
  else
    write(iout) 0
  endif

  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
iout :integer, intent(in)

Save configuration for continuation to a file indicated by an unit id.

[Source]

subroutine save_contg_qcd_hmc(this,iout)
!
! Save configuration for continuation to a file indicated by an unit id.
!
  use comlib
  use print_status_class
  use random_class
  implicit none
  type(qcd_hmc), intent(inout) :: this
  integer,       intent(in)    :: iout

  real(DP) :: plq


  call plaquette(this%u,plq)

  write(iout)this%u
  call save(g_rand,iout)
  write(iout)plq

  if (nodeid == 0) then
    write(*,'(80("-"))')
    write(*,'(" PLQ :",F24.16)')plq
  endif

  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
idepth :integer, intent(in)

[Source]

subroutine set_current_md_depth_qh(this,idepth)
  implicit none
  type(qcd_hmc), intent(inout) :: this
  integer, intent(in) :: idepth
  this%idepth = idepth
  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
: HMC QCD action
dt :real(DP), intent(in)
: time step for MD

Update Momentum (P)

p = p + dt F

[Source]

subroutine update_p_qcd_hmc(this,dt)
!
! Update Momentum (P)
!
! p = p + dt F
!
  implicit none
  type(qcd_hmc), intent(inout) :: this  ! HMC QCD action
  real(DP),      intent(in)    :: dt    ! time step for MD

  call force(this)

  call update_p(this%p,dt,this%F)

  if (associated(this%F)) then
    deallocate(this%F)
  endif

  return
end subroutine
Subroutine :
this :type(qcd_hmc), intent(inout)
: HMC QCD action
dt :real(DP), intent(in)
: time step for MD

Update Corrdinate (Q)

q = q + dt p

[Source]

subroutine update_q_qcd_hmc(this,dt)
!
! Update Corrdinate (Q)
!
! q = q + dt p
!
  implicit none
  type(qcd_hmc), intent(inout) :: this ! HMC QCD action
  real(DP),      intent(in)    :: dt   ! time step for MD

  call update_u(this%u,dt,this%p)

  return
end subroutine
wloops_obj
Derived Type :
PLQ =0.0_DP :real(DP)
: PLQ = Sum[Re[1x1],vol,direction]/direction/Vol
RCT =0.0_DP :real(DP)
: RCT = Sum[Re[2x1+1x2],vol,direction]/direction/Vol
    = 2 Sum[Re[1x2],vol,direction]/direction/Vol

contains Wilson loop/plaquette values

Original external subprogram is hmc_gluon_class#wloops_obj