Class | hmc_truncd_overlap_class |
In: |
HmcTruncdOverlapClass/hmc_truncd_overlap_class.F90
|
Action : $ S_Q $
$ S_Q $ : Pseudo-fermion HMC action, $ S_Q = |D_{OV}^{-1} phi|^{2} $
where
$ D_{OV} = e^{dag} P^{dag} D_{DW}(1)^{-1} D_{DW}(m_q) P e $
and
$ D_{OV}^{-1} = e^{dag} P^{dag} D_{DW}(m_q)^{-1} D_{DW}(1) P e $
$ P $ : premutation and chiral projection in 5-Dim.
$ e $ : prolongation $ e $ or restriction $ e^{dag} $ operaotr between 4-D and 5-D.
$ D_{DW}(m) $ : Domainwall operator with mass $ m $ .
$Id: hmc_truncd_overlap_class.F90,v 1.30 2011/08/11 02:53:03 ishikawa Exp $
Subroutine : | |
this : | type(hmc_truncd_dw_overlap), intent(inout) |
delete HMC quark
subroutine delete_tdwov_hmc(this) ! ! \delete HMC quark ! implicit none type(hmc_truncd_dw_overlap), intent(inout) :: this deallocate(this%y) NULLIFY(this%y) call delete(this%target_phys_param) call delete(this%phys_param) call delete(this%pv_solver_log) call delete(this%phys_solver_log) call delete(this%reweighting_log) call delete(this%lowmode_log) deallocate(this%target_phys_param) deallocate(this%phys_param) return end subroutine
Subroutine : | |||
this : | type(hmc_truncd_dw_overlap),intent(inout)
| ||
BB : | type(vfield_gluon_wg), intent(inout)
| ||
idepth : | integer, intent(in)
| ||
u : | type(vfield_gluon_wg), intent(in)
|
Calc MD force from hopping part
Action : $ S_Q $
$ S_Q $ : Pseudo-fermion HMC action, $ S_Q = |D_{OV}^{-1} phi|^{2} $
where
$ D_{OV} = e^{dag} P^{dag} D_{DW}(1)^{-1} D_{DW}(m_q) P e $
$ D_{OV}^{-1} = e^{dag} P^{dag} D_{DW}(m_q)^{-1} D_{DW}(1) P e $
$ P $ : premutation and chiral projection in 5-Dim.
$ e $ : prolongation $ e $ or restriction $ e^{dag} $ operaotr between 4-D and 5-D.
$ D_{DW}(m) $ : Domainwall operator with mass $ m $ .
Force computation:
$ delta S_Q = - X^{dag} D_{OV}^{-1} delta D_{OV} X + \{ H. C.\}, $
with $ X equiv D_{OV}^{-1}phi. $
This reduces to
\[ delta S_Q = -frac{1}{2} fx^{dag} delta M^{(NS)}_{hop} fy + \{ H. C.\}, \]
with
$ fx equiv left( D_{DW}(1)^{dag} right)^{-1} P e left( D_{OV}^{dag} right)^{-1} X, $
$ fy equiv left( D_{DW}(1)^{-1} D_{DW}(m_q) - 1 right) P e X, $
and
$ M^{(NS)}_{hop} $ : 5-dim hopping matrix (diagonal in 5-dim).
subroutine force_tdwov_hmc(this,BB,idepth,u) ! ! Calc MD \force from hopping part ! ! Action : $ S_Q $ ! ! $ S_Q $ : Pseudo-fermion HMC action, $ S_Q = |D_{OV}^{-1} \phi|^{2} $ ! ! where ! ! $ D_{OV} = e^{\dag} P^{\dag} D_{DW}(1)^{-1} D_{DW}(m_q) P e $ ! ! ! $ D_{OV}^{-1} = e^{\dag} P^{\dag} D_{DW}(m_q)^{-1} D_{DW}(1) P e $ ! ! $ P $ : premutation and chiral projection in 5-Dim. ! ! $ e $ : prolongation $ e $ or restriction $ e^{\dag} $ operaotr between 4-D and 5-D. ! ! $ D_{DW}(m) $ : Domainwall operator with mass $ m $ . ! !--------------------------------------------------------------------------------- ! ! \Force computation: ! ! $ \delta S_Q = - X^{\dag} D_{OV}^{-1} \delta D_{OV} X + \{ H. C.\}, $ ! ! with $ X \equiv D_{OV}^{-1}\phi. $ ! ! This reduces to ! !\[ ! \delta S_Q = -\frac{1}{2} fx^{\dag} \delta M^{(NS)}_{hop} fy + \{ H. C.\}, !\] ! ! with ! ! $ fx \equiv \left( D_{DW}(1)^{\dag} \right)^{-1} P e \left( D_{OV}^{\dag} \right)^{-1} X, $ ! ! $ fy \equiv \left( D_{DW}(1)^{-1} D_{DW}(m_q) - 1 \right) P e X, $ ! ! and ! ! $ M^{(NS)}_{hop} $ : 5-dim hopping matrix (diagonal in 5-dim). ! ! implicit none type(hmc_truncd_dw_overlap),intent(inout) :: this ! quark action/algorithm parameters type(vfield_gluon_wg), intent(inout) :: BB ! pre-force contribution from hopping part integer, intent(in) :: idepth ! MD depth for recursive MD integrator type(vfield_gluon_wg), intent(in) :: u ! gauge field type(field_dw_quark_wg), allocatable :: fx,fy type(field_dw_quark_wg), allocatable :: pw,qw type(field_quark_wg), allocatable :: w real(DP) :: err,rtmp,tol integer :: iout,jout integer :: iiter,NS character(len=CHARLEN) :: str if (idepth == this%depth_hopping) then #ifdef _MD_DEPTH_ if (nodeid==0) write(*,'("depth=",I3," : force quark tdwov hmc hopping")')idepth #endif write(str,'("# FORCE SOLVER")') call print(this%phys_solver_log,TRIM(ADJUSTL(str))) call print( this%pv_solver_log,TRIM(ADJUSTL(str))) iout = get_file_unit(this%phys_solver_log) jout = get_file_unit(this%pv_solver_log) NS = get_NS(this%phys_param) allocate(fx,fy,pw,w) call new(fx,NS) call new(fy,NS) call new(pw,NS) call new(w) !====================================== ! fy <= ( (B+CM5(m)) Ddw(m) \ Ddw(1) - (B+CM5(1)) ) P y !====================================== ! ! pw <= P y ! call assign_proj_4dto5d(pw,this%y) ! ! fy <= Ddw(1) pw ! call set_mass_pauli_villars(this%phys_param) call assign_mult_dwf(this%phys_param,fy,pw,u) call set_mass_physical(this%phys_param) ! ! fx <= Ddw(m) \ fy ! tol = this%force_phys_inv%tol iiter = MAX_ITER call assign_inv_mult_dwf(this%phys_param,iout,tol,iiter,fx,fy,u,cron=this%cron1) call inc(this%force_phys_inv,iiter) !====================================== ! w <= Dov \ y ! = P' Ddw(m) \ Ddw(1) P y = P' fx !====================================== call assign_proj_5dto4d(w,fx) ! fy <= (B+CM(m))fx call assign_mult_impterm_dwf(this%phys_param,fy,fx) ! fx <= (B+CM(1))pw call set_mass_pauli_villars(this%phys_param) call assign_mult_impterm_dwf(this%phys_param,fx,pw) call set_mass_physical(this%phys_param) ! ! fy <= fy - fx ! call accum_sub(fy,fx) !====================================== ! fx <= Ddw(m)' \ P w ! = Gamma5 Ddw(m) \ Gamma5 P w !====================================== call assign_proj_4dto5d(pw,w) deallocate(w) tol = this%force_phys_inv%tol iiter = MAX_ITER call assign_inv_mult_dwf(this%phys_param,iout,tol,iiter,fx,pw,u,dagger=.true.,cron=this%cron2) call inc(this%force_phys_inv,iiter) call delete(pw) deallocate(pw) !====================================== ! force from hopping matrix ! ! delta S = coef fx \dalta M fy + {H.c} ! !====================================== rtmp = +0.5_DP call force_hmc_hopping(BB,rtmp,fx,fy) call delete(fx) call delete(fy) deallocate(fx,fy) endif return end subroutine
Subroutine : | |||
this : | type(hmc_truncd_dw_overlap),intent(inout)
| ||
ifi : | integer, intent(in)
| ||
u : | type(vfield_gluon_wg), intent(in)
| ||
action : | real(DP), intent(out)
|
Calculate Hamiltonian (Action) for HMC quark
Action :
$ S_Q $ : Pseudo-fermion HMC action, $ S_Q = |D_{OV}^{-1} phi|^{2} $
where
$ D_{OV} = e^{dag} P^{dag} D_{DW}(1)^{-1} D_{DW}(m_q) P e $
$ D_{OV}^{-1} = e^{dag} P^{dag} D_{DW}(m_q)^{-1} D_{DW}(1) P e $
$ P $ : premutation and chiral projection in 5-Dim.
$ e $ : prolongation $ e $ or restriction $ e^{dag} $ operaotr between 4-D and 5-D.
$ D_{DW}(m) $ : Domainwall operator with mass $ m $ .
subroutine hamil_tdwov_hmc(this,ifi,u,action) ! ! Calculate Hamiltonian (Action) for HMC quark ! ! Action : ! ! $ S_Q $ : Pseudo-fermion HMC action, $ S_Q = |D_{OV}^{-1} \phi|^{2} $ ! ! where ! ! $ D_{OV} = e^{\dag} P^{\dag} D_{DW}(1)^{-1} D_{DW}(m_q) P e $ ! ! ! $ D_{OV}^{-1} = e^{\dag} P^{\dag} D_{DW}(m_q)^{-1} D_{DW}(1) P e $ ! ! $ P $ : premutation and chiral projection in 5-Dim. ! ! $ e $ : prolongation $ e $ or restriction $ e^{\dag} $ operaotr between 4-D and 5-D. ! ! $ D_{DW}(m) $ : Domainwall operator with mass $ m $ . ! use comlib use hmc_status_class implicit none type(hmc_truncd_dw_overlap),intent(inout) :: this ! quark action/algorithm parameters integer, intent(in) :: ifi ! target switch : 0 for initial generating PF field, 1 for final type(vfield_gluon_wg), intent(in) :: u ! gauge field real(DP), intent(out) :: action ! output action value real(DP) :: Spf type(field_quark_wg), allocatable :: w,v type(field_dw_quark_wg), allocatable :: s,ds character(CHARLEN) :: str real(DP) :: tol,rtmp0,rtmp1 integer :: iiter,iout,i,NEV,NS type(hmc_status) :: status Spf = 0.0_DP if (0 /= this%depth_hopping) then select case(ifi) case(0) ! INITIAL if (this%ncron > 0) then !--------------------------------- ! set up for chronological guess !--------------------------------- allocate(this%cron1) allocate(this%cron2) NS = get_NS(this%phys_param) call new(this%cron1,NSIZE=NS*CRON_NSITE,NMAX=this%ncron) call new(this%cron2,NSIZE=NS*CRON_NSITE,NMAX=this%ncron) endif allocate(w) call new(w) !===================== ! Heatbath !===================== call set_gaussian_noise(w) !===================== ! calc. action !===================== Spf=abs2(w) !===================== ! ! y <= Dov w !===================== write(str,'("# HAMIL SOLVER (INITIAL)")') call print( this%pv_solver_log,TRIM(ADJUSTL(str))) iout = get_file_unit(this%pv_solver_log) iiter = MAX_ITER tol = this%hamil_pv_inv%tol call assign_mult_ov(this%phys_param,iout,tol,iiter,this%y,w,u) call inc(this%hamil_pv_inv,iiter) deallocate(w) case(1:2) ! FINAL or REVCHECK allocate(w) call new(w) !============================ ! ! w <= Dov \ y !============================ write(str,'("# HAMIL SOLVER (FINAL)")') call print(this%phys_solver_log,TRIM(ADJUSTL(str))) iout = get_file_unit(this%phys_solver_log) iiter = MAX_ITER tol = this%hamil_phys_inv%tol call assign_inv_mult_ov(this%phys_param,iout,tol,iiter,w,this%y,u,this%cron1) call inc(this%hamil_phys_inv,iiter) !=================================== ! action = |Dov \ y|^2 = |w|^2 !=================================== Spf=abs2(w) deallocate(w) if (this%ncron > 0) then !--------------------------------- ! delete chronological guess !--------------------------------- call delete(this%cron1) call delete(this%cron2) deallocate(this%cron1) deallocate(this%cron2) endif #ifdef _REVCHECK if (1 == ifi) then if (this%ncron > 0) then !--------------------------------- ! set up for chronological guess !--------------------------------- allocate(this%cron1) allocate(this%cron2) NS = get_NS(this%phys_param) call new(this%cron1,NSIZE=NS*CRON_NSITE,NMAX=this%ncron) call new(this%cron2,NSIZE=NS*CRON_NSITE,NMAX=this%ncron) endif endif #endif end select endif action=Spf write(str,'("@",I8,I2," Quark: ID:",I3," SQPF:",E24.16)') get_trajectory_number(status),ifi,this%myid,Spf call print_log_action(status,str) return end subroutine
Function : | |
has_gmp : | logical |
this : | type(hmc_truncd_dw_overlap), intent(in) |
function has_gmp_tdwov_hmc(this) result(has_gmp) implicit none type(hmc_truncd_dw_overlap), intent(in) :: this logical :: has_gmp has_gmp = this%has_global_metropolis_test return end function
Function : | |
has_rew : | logical |
this : | type(hmc_truncd_dw_overlap), intent(in) |
function has_reweighting_tdwov_hmc(this) result(has_rew) implicit none type(hmc_truncd_dw_overlap), intent(in) :: this logical :: has_rew has_rew = this%has_reweighting return end function
Derived Type : | |||
target_phys_param => NULL() : | type(quark_truncd_dw_overlap), public, pointer
| ||
phys_param => NULL() : | type(quark_truncd_dw_overlap), public, pointer
|
Subroutine : | |||
this : | type(hmc_truncd_dw_overlap),intent(inout)
| ||
u1 : | type(vfield_gluon_wg), intent(in)
| ||
u0 : | type(vfield_gluon_wg), intent(in)
| ||
ihit : | integer, intent(inout)
|
Nothing is done
subroutine metropolis_tdwov_hmc(this,u1,u0,ihit) ! ! Nothing is done ! implicit none type(hmc_truncd_dw_overlap),intent(inout) :: this ! quark action/algorithm parameters type(vfield_gluon_wg), intent(in) :: u1 ! MD final gauge field type(vfield_gluon_wg), intent(in) :: u0 ! MD inital gauge field integer, intent(inout) :: ihit ! 1 for accept, 0 for reject return end subroutine
Subroutine : | |
this : | type(hmc_truncd_dw_overlap), intent(inout) |
id : | integer, intent(in) |
initialize/create
subroutine new_tdwov_hmc(this,id) ! ! initialize/create ! implicit none type(hmc_truncd_dw_overlap), intent(inout) :: this integer, intent(in) :: id character(len=CHARLEN) :: str NULLIFY(this%y) this%myid = id allocate(this%target_phys_param) allocate(this%phys_param) call new(this%target_phys_param,id) call new(this%phys_param,id) call new(this%force_pv_inv) call new(this%force_phys_inv) call new(this%hamil_pv_inv) call new(this%hamil_phys_inv) this%depth_hopping = 1 ! ! append action id to the log file names ! write(str,'(I2.2)')id this%pv_solver_log_fname = TRIM(this%pv_solver_log_fname)//TRIM(str) this%phys_solver_log_fname = TRIM(this%phys_solver_log_fname)//TRIM(str) this%reweighting_log_fname = TRIM(this%reweighting_log_fname)//TRIM(str) this%lowmode_log_fname = TRIM(this%lowmode_log_fname)//TRIM(str) call new( this%pv_solver_log, TRIM(ADJUSTL(this%pv_solver_log_fname))) call new(this%phys_solver_log,TRIM(ADJUSTL(this%phys_solver_log_fname))) call new(this%reweighting_log,TRIM(ADJUSTL(this%reweighting_log_fname))) call new( this%lowmode_log, TRIM(ADJUSTL(this%lowmode_log_fname))) return end subroutine
Subroutine : | |
this : | type(hmc_truncd_dw_overlap), intent(inout) |
print out HMC
subroutine print_tdwov_hmc(this) ! ! \print out HMC ! implicit none type(hmc_truncd_dw_overlap), intent(inout) :: this if (nodeid==0) then write(*,'("**** Quark # ",I3," HMC Truncated-DW-Overlap")')this%myid write(*,'("***** Physical parameters *****")') if (this%has_reweighting) then write(*,'("*** Reweighting Target Action paramter")') call print(this%target_phys_param) write(*,'("** Reweighting paramter")') write(*,'(" Number of Noise :",I4)')this%num_noise_rew write(*,'(" Number of Low-modes :",I4)')this%NEV endif write(*,'("*** HMC/MD Action parameter")') call print(this%phys_param) write(*,'("***** HMC Algorithm parameters *****")') write(*,'(" Force Stopping Condition :",E16.8)')this%force_phys_inv%tol write(*,'(" Hamil Stopping Condition :",E16.8)')this%hamil_phys_inv%tol write(*,'(" MD force depth PSfermion part :",I3)')this%depth_hopping write(*,'(" Chronological History :",I3)')this%ncron write(*,'(80("-"))') endif return end subroutine
Subroutine : | |
this : | type(hmc_truncd_dw_overlap), intent(inout) |
print out HMC quark
subroutine print_stat_tdwov_hmc(this) ! ! \print out HMC quark ! implicit none type(hmc_truncd_dw_overlap), intent(inout) :: this real(DP) :: fph_ave,hph_ave real(DP) :: fpv_ave,hpv_ave fph_ave = get_average(this%force_phys_inv) hph_ave = get_average(this%hamil_phys_inv) fpv_ave = get_average(this%force_pv_inv) hpv_ave = get_average(this%hamil_pv_inv) if (nodeid==0) then write(*,'("**** Quark # ",I3," HMC Truncated-DW-Overlap")')this%myid if (this%has_reweighting) then write(*,'("*** Reweighting Target Action paramter")') call print(this%target_phys_param) write(*,'("** Reweighting paramter")') write(*,'(" Number of Noise :",I4)')this%num_noise_rew write(*,'(" Number of Low-modes :",I4)')this%NEV endif write(*,'("*** HMC/MD Action parameter")') call print(this%phys_param) write(*,'(" HMC algorithm ")') write(*,'(" Force Solver (Phys mode) iter")') write(*,'(" (averaged) :",F12.4)') fph_ave write(*,'(" Force Solver (PV mode) iter")') write(*,'(" (averaged) :",F12.4)') fpv_ave write(*,'(" Hamil Solver (Phys mode) iter")') write(*,'(" (averaged) :",F12.4)') hph_ave write(*,'(" Hamil Solver (PV mode) iter")') write(*,'(" (averaged) :",F12.4)') hpv_ave write(*,'(" Chronological History :",I3)')this%ncron endif return end subroutine
Subroutine : | |
this : | type(hmc_truncd_dw_overlap), intent(inout) |
iout : | integer, intent(in) |
read HMC quark
subroutine read_tdwov_hmc(this,iout) ! ! \read HMC quark ! use comlib implicit none type(hmc_truncd_dw_overlap), intent(inout) :: this integer, intent(in) :: iout integer :: reweighting_flag allocate(this%y) call new(this%y) call new(this%hamil_pv_inv) call new(this%hamil_phys_inv) call new(this%force_pv_inv) call new(this%force_phys_inv) !------------------------------------------ ! read reweighting job parmeters !------------------------------------------ if (nodeid==0) then read(iout,*) ! skip comment line endif if (nodeid==0) then read(iout,*)reweighting_flag read(iout,*)this%num_noise_rew read(iout,*)this%NEV read(iout,*) ! skip comment line endif #ifndef _singlePU call comlib_bcast(reweighting_flag,0) call comlib_bcast(this%num_noise_rew,0) call comlib_bcast(this%NEV,0) #endif if (reweighting_flag == 0) then this%has_reweighting = .false. else this%has_reweighting = .true. endif !------------------------------------------ ! read reweighting target quark parameters !------------------------------------------ call read(this%target_phys_param,iout) !------------------------------------------ ! read HMC quark parameters !------------------------------------------ if (nodeid==0) then read(iout,*) ! skip comment line endif call read(this%phys_param,iout) call set_coef(this%phys_param) ! set/pre-compute DW fermion action parameters for HMC quark !------------------------------------------ ! read HMC algorithm parameters !------------------------------------------ if (nodeid==0) then read(iout,*) ! skip comment line read(iout,*)this%hamil_phys_inv%tol read(iout,*)this%force_phys_inv%tol read(iout,*)this%depth_hopping read(iout,*)this%ncron endif #ifndef _singlePU call comlib_bcast(this%hamil_phys_inv%tol,0) call comlib_bcast(this%force_phys_inv%tol,0) call comlib_bcast(this%depth_hopping,0) call comlib_bcast(this%ncron,0) #endif this%hamil_pv_inv%tol = this%hamil_phys_inv%tol this%force_pv_inv%tol = this%force_phys_inv%tol return end subroutine
Subroutine : | |||
this : | type(hmc_truncd_dw_overlap),intent(inout)
| ||
u : | type(vfield_gluon_wg), intent(in) |
Compute reweighting factor for quark determinant
Extra dimension size NS is enlarged by reweighting
\[ W = exp( - |D_{OV}(target) \ D_{OV}(HMC) eta|^2 + |eta|^2) \] with gaussian noise $ eta $ .
subroutine reweighting_tdwov_hmc(this,u) ! ! Compute reweighting factor for quark determinant ! ! Extra dimension size NS is enlarged by reweighting ! !\[ ! W = \exp( - |D_{OV}(target) \ D_{OV}(HMC) \eta|^2 + |\eta|^2) !\] ! with gaussian noise $ \eta $ . ! ! use timer_class use hmc_status_class implicit none type(hmc_truncd_dw_overlap),intent(inout) :: this ! quark action/algorithm parameters type(vfield_gluon_wg), intent(in) :: u type(field_quark_wg), allocatable :: eta,v,w type(dwf_low_modes), allocatable :: LV integer :: iout,iiter,inoise,itraj real(DP), parameter :: EIGEN_TOL = 1.0e-13_DP real(DP) :: tol,rtmp0,rtmp1,rtmp2,mq, eval_min real(DP) :: ave_dS, ave_exp complex(DP) :: ctmp character(CHARLEN) :: str logical :: is_lowprojd type(timer) :: etime type(hmc_status) :: status if (.not.this%has_reweighting) return call new(etime) call tic(etime) if (this%NEV > 0) then is_lowprojd = .true. else is_lowprojd = .false. endif if (is_lowprojd) then !--------------------------------------------------- ! overlap kernel lowmode comptation ! kernel = Hw for Domain-Wall fermion !--------------------------------------------------- allocate(LV) write(str,'("# REWEIGHTING SOLVER")') call print(this%lowmode_log,TRIM(ADJUSTL(str))) iout = get_file_unit(this%lowmode_log) tol = EIGEN_TOL iiter = MAX_ITER LV%NEV = this%NEV call new(LV) call get_low_modes(LV,this%target_phys_param,iout,tol,iiter,u) else !--------------------------------------------------- ! overlap kernel lowmode comptation ! kernel = Hw for Domain-Wall fermion !--------------------------------------------------- allocate(LV) write(str,'("# REWEIGHTING SOLVER")') call print(this%lowmode_log,TRIM(ADJUSTL(str))) iout = get_file_unit(this%lowmode_log) tol = EIGEN_TOL iiter = MAX_ITER LV%NEV = 6 call new(LV) call get_low_modes(LV,this%target_phys_param,iout,tol,iiter,u) call delete(LV) deallocate(LV) endif call set_coef(this%target_phys_param) allocate(eta) call new(eta) itraj = get_trajectory_number(status) ave_dS = 0.0_DP ave_exp = 0.0_DP do inoise = 1, this%num_noise_rew call set_gaussian_noise(eta) allocate(v,w) call new(v) call new(w) !------------------------------------ ! v = Dov(HMC) eta !------------------------------------ write(str,'("# REWEIGHTING SOLVER")') call print(this%pv_solver_log,TRIM(ADJUSTL(str))) iout = get_file_unit(this%pv_solver_log) iiter = MAX_ITER tol = this%hamil_pv_inv%tol call assign_mult_ov(this%phys_param,iout,tol,iiter,v,eta,u) call inc(this%hamil_pv_inv,iiter) !------------------------------------ ! w = Dov(target) \ v !------------------------------------ write(str,'("# REWEIGHTING SOLVER")') call print(this%phys_solver_log,TRIM(ADJUSTL(str))) iout = get_file_unit(this%phys_solver_log) iiter = MAX_ITER tol = this%hamil_phys_inv%tol if (is_lowprojd) then call assign_inv_mult_ov(this%target_phys_param,LV,iout,tol,iiter,w,v,u) else call assign_inv_mult_ov(this%target_phys_param, iout,tol,iiter,w,v,u) endif call inc(this%hamil_phys_inv,iiter) !------------- ! v = w + eta ! w = w - eta !------------- call assign_add(v,w,eta) call accum_sub(w,eta) !------------------------------------------------------ ! exp_weight = |Dov(target)\Dov(HMC) eta|^2 - |eta|^2 ! = Re(w'*v) ! weight = exp(-exp_weight) !------------------------------------------------------ rtmp1 = real(prod(w,v),kind=DP) rtmp2 = exp(-rtmp1) write(str,'(I7," dS_=",E24.15," exp(-dS)_=",E24.15," NOISE_=",I3)')itraj,rtmp1,rtmp2,inoise call print(this%reweighting_log,str) ave_dS = ave_dS + rtmp1 ave_exp = ave_exp + rtmp2 deallocate(v,w) !---------------------- ! check GW relation !---------------------- if (is_lowprojd) then call check_GW_relation_tdwov(this,LV=LV,eta=eta,u=u) else call check_GW_relation_tdwov(this, eta=eta,u=u) endif enddo ave_dS = ave_dS /this%num_noise_rew ave_exp = ave_exp/this%num_noise_rew write(str,'("@",I6," dS_=",E24.15," exp(-dS)_=",E24.15," NOISE_AVERAGED")')itraj,ave_dS,ave_exp call print(this%reweighting_log,str) deallocate(eta) if (is_lowprojd) then call delete(LV) deallocate(LV) endif call toc(etime) rtmp0 = get_elapse(etime) write(str,'("# ETIME_=",E24.15)')rtmp0 call print(this%reweighting_log,str) return end subroutine
Subroutine : | |||
this : | type(hmc_truncd_dw_overlap), intent(in) | ||
iout : | integer, intent(in)
|
save HMC Quark
subroutine save_config_tdwov_hmc(this,iout) ! ! \save HMC Quark ! implicit none type(hmc_truncd_dw_overlap), intent(in) :: this integer, intent(in) :: iout ! configuration file unit id call save_config(this%target_phys_param,iout) call save_config(this%phys_param,iout) return end subroutine
Subroutine : | |
this : | type(hmc_truncd_dw_overlap), intent(inout) |
id : | integer, intent(in) |
set quark # id
subroutine set_id_tdwov_hmc(this,id) ! ! set quark # id ! implicit none type(hmc_truncd_dw_overlap), intent(inout) :: this integer, intent(in) :: id this%myid = id call set_id(this%target_phys_param,id) call set_id(this%phys_param,id) return end subroutine