ALF  dev.
A QMC Code for fermionic models
Functions/Subroutines | Variables
global_mod Module Reference

Handles global updates and parallel tempering. More...

Functions/Subroutines

subroutine exchange_step (Phase, GR, udvr, udvl, Stab_nt, udvst, N_exchange_steps, Tempering_calc_det)
 Handles parrallel tempering
This subroutine is called only if the tempering flag is switched on
Requires MPI. More...
 
subroutine global_updates (Phase, GR, udvr, udvl, Stab_nt, udvst, N_Global)
 Carries out N_global global updates as defeined in the Global_move routine of the Hamiltonian module. More...
 
complex(kind=kind(0.d0)) function compute_ratio_global (Phase_Det_old, Phase_Det_new, Det_vec_old, Det_vec_new, nsigma_old, T0_Proposal_ratio, Ratio)
 This fucntion computes ratio of weights for global moves. More...
 
subroutine compute_fermion_det (Phase_det, Det_Vec, udvl, udvst, Stab_nt, storage)
 Computes the fermion determinant from scratch or from the storage. More...
 
integer function npbc_tempering (n, Isize)
 Periodic boundary conditions required to defined master and slave for the tempering. More...
 
subroutine global_tempering_setup
 Allocates the Tempering_acceptance (Type Obser_vec) variable that monitors the exchange acceptance. More...
 
subroutine global_tempering_init_obs
 Initializes Tempering_acceptance. More...
 
subroutine global_tempering_obser (toggle)
 Measures Tempering_acceptance. More...
 
subroutine global_tempering_pr
 Prints Tempering_acceptance. More...
 

Variables

type(obser_vec), private tempering_acceptance
 

Detailed Description

Handles global updates and parallel tempering.

Author
ALF-project

Function/Subroutine Documentation

◆ compute_fermion_det()

subroutine global_mod::compute_fermion_det ( complex (kind=kind(0.d0)), dimension(:), intent(out)  Phase_det,
real (kind=kind(0.d0)), dimension(:,:), intent(out)  Det_Vec,
class(udv_state), dimension(:), intent(inout), allocatable  udvl,
class(udv_state), dimension(:,:), intent(inout), allocatable  udvst,
integer, dimension(:), intent(in), allocatable  Stab_nt,
character (len=64), intent(in)  storage 
)

Computes the fermion determinant from scratch or from the storage.

Author
 Finite temperature: If the storage is full then computes 
                     det( 1 +  VL*DL*UL^{dag} )= Phase_det*e^{ \sum_{n=1}^{ndim} Det_vec(n)}
                     Note that  Udvl%U = UL, Udvl%D = DL and  Udvl%V = VL^{dag}
                     If storage is empty one first computes  Udvl and in doing so fills up the storage.
 Zero temperature:   If storage is full then computes  det \Psi_L U(\theta_tot,0) \Psi_R  =  Phase_det*e^{ \sum_{n=1}^{N_part} Det_vec(n) }
                     For the projective code the required scales which one can throw away for the Green function are in udvst.
                     If the storage is not full, then things will be computed from scratch and intermediate results will be stored in udvst.
Parameters
[in,out]udvl,udvrClass(UDV_State)
  If storage=Full  : On entry and exit udvl contains the last udv decomposition such that G=(1 + VL*DL*UL^{dag})^{-1}
  If storage=Empty : On exit udv1 will contain the last udv decomposition such that G=(1 + VL*DL*UL^{dag})^{-1}
[out]Phase_detComplex, Dimension(N_FL)
  The phase, per flavor
[out]Det_vecReal, Dimension(:,N_FL)
[in]Stab_ntInteger, Dimension(:)
  List of time slices for stabilization.
[in]StorageCharacter
  storage = "Full".  Compute the fermion det  with the use of the storage.
  storage = "Empty". Compute the fermion det from scratch and  in doing so, fill the storage.

Local variables

Definition at line 780 of file Global_mod.F90.

◆ compute_ratio_global()

complex (kind=kind(0.d0)) function global_mod::compute_ratio_global ( complex (kind=kind(0.d0)), dimension(:), intent(in), allocatable  Phase_Det_old,
complex (kind=kind(0.d0)), dimension(:), intent(in), allocatable  Phase_Det_new,
real (kind=kind(0.d0)), dimension(:,:), intent(in), allocatable  Det_vec_old,
real (kind=kind(0.d0)), dimension(:,:), intent(in), allocatable  Det_vec_new,
type (fields), intent(in)  nsigma_old,
real (kind=kind(0.d0)), intent(in)  T0_Proposal_ratio,
complex (kind=kind(0.d0)), dimension(2), intent(out)  Ratio 
)

This fucntion computes ratio of weights for global moves.

Author
The ALF Project contributors
  Ratio_Global = T0(nsigma--> nsigma_old) W(nsigma)/T0(nsigma_old--> nsigma) W(nsigma_old) =
                 T0_Proposal_ratio   W(nsigma)/  W(nsigma_old)
  On entry the old and new fermion determinants read: Phase_det*e^{ \sum_{n=1}^{ndim} Det_vec(n) }
  as obtained from the Compute_Fermion_det routine.
  On exit Ratio_Global = Ratio(1)*exp(Ratio(2))
Parameters
[in]Phase_det_newComplex, Dimension(N_FL)
[in]Phase_det_oldComplex, Dimension(N_FL)
[in]Det_vec_newReal, Dimension(:,N_FL)
[in]Det_vec_oldReal, Dimension(:,N_FL)
[in]T0_proposal_ratioReal
[in]nsigma_oldType(Fields)
  Old configuration. The new configuration is stored in nsigma. nsigma is a globale variable
  contained in the Hamiltonian module.

Definition at line 671 of file Global_mod.F90.

◆ exchange_step()

subroutine global_mod::exchange_step ( complex (kind=kind(0.d0)), intent(inout)  Phase,
complex (kind=kind(0.d0)), dimension(:,:,:), intent(inout), allocatable  GR,
class(udv_state), dimension(:), intent(inout), allocatable  udvr,
class(udv_state), dimension(:), intent(inout), allocatable  udvl,
integer, dimension(:), intent(in), allocatable  Stab_nt,
class(udv_state), dimension(:, :), intent(inout), allocatable  udvst,
integer, intent(in)  N_exchange_steps,
logical, intent(in)  Tempering_calc_det 
)

Handles parrallel tempering
This subroutine is called only if the tempering flag is switched on
Requires MPI.

Author
ALF-project

On entry and on exit the left storage is full, the Green function is on time slice 0 and the phase is set. The same holds on exit but with updated quantities.

Parameters
[in,out]PhaseComplex
  Is updated upon acceptance
[in,out]udvl,udvrClass(UDV_State)
  On entry udvl contains the last udv decomposition such that G=(1 + VL*DL*UL^{dag})^{-1}
  udvr  is used as working space
[in,out]GRComplex
  Green function. Is updated upon acceptance.
[in,out]udvstClass(UDV_State)
  Storage. Is updated with left propagation upon acceptance
[in]Stab_ntInteger
  List of time slices for stabilization.
[in]N_exchange_stepsInteger
 Number of exchange steps
[in]Tempering_calc_detLogical
  If true then the fermion determinant is computed in the calulcation of the exchange weight.

Local variables.

Store old configuration

Send result of acceptance/rejection decision form master to slave

Definition at line 109 of file Global_mod.F90.

◆ global_tempering_init_obs()

subroutine global_mod::global_tempering_init_obs ( )

Initializes Tempering_acceptance.

Author
The ALF Project contributors

Definition at line 995 of file Global_mod.F90.

◆ global_tempering_obser()

subroutine global_mod::global_tempering_obser ( logical, intent(in)  toggle)

Measures Tempering_acceptance.

Author
The ALF Project contributors

Definition at line 1007 of file Global_mod.F90.

◆ global_tempering_pr()

subroutine global_mod::global_tempering_pr ( )

Prints Tempering_acceptance.

Author
The ALF Project contributors

Definition at line 1023 of file Global_mod.F90.

◆ global_tempering_setup()

subroutine global_mod::global_tempering_setup ( )

Allocates the Tempering_acceptance (Type Obser_vec) variable that monitors the exchange acceptance.

Author
The ALF Project contributors

Definition at line 979 of file Global_mod.F90.

◆ global_updates()

subroutine global_mod::global_updates ( complex (kind=kind(0.d0)), intent(inout)  Phase,
complex (kind=kind(0.d0)), dimension(:,:,:), intent(inout), allocatable  GR,
class (udv_state), dimension(:), intent(inout), allocatable  udvr,
class (udv_state), dimension(:), intent(inout), allocatable  udvl,
integer, dimension(:), intent(in), allocatable  Stab_nt,
class(udv_state), dimension(:,:), intent(inout), allocatable  udvst,
integer, intent(in)  N_Global 
)

Carries out N_global global updates as defeined in the Global_move routine of the Hamiltonian module.

Author
The ALF Project contributors

On entry and on exit the left storage is full, the Green function is on time slice 0 and the phase is set. The same holds on exit but with updated quantities.

Parameters
[in,out]PhaseComplex
  Is updated upon acceptance
[in,out]udvl,udvrClass(UDV_State)
  On entry udvl contains the last udv decomposition such that G=(1 + VL*DL*UL^{dag})^{-1}
  udvr  is used as working space
[in,out]GRComplex
  Green function. Is updated upon acceptance.
[in,out]udvstClass(UDV_State)
  Storage. Is updated with left propagation upon acceptance
[in]Stab_ntInteger
  List of time slices for stabilization.
[in]N_GlobalInteger
  Number of global moves that will be caried out

Definition at line 451 of file Global_mod.F90.

◆ npbc_tempering()

integer function global_mod::npbc_tempering ( integer, intent(in)  n,
integer, intent(in)  Isize 
)

Periodic boundary conditions required to defined master and slave for the tempering.

Author
The ALF Project contributors

Definition at line 961 of file Global_mod.F90.