Examples
In these section, we will use some typical examples to show you how to apply the Fortran binding and explicit Fortran types to control the atomic eigenvalue problem solver, CT-HYB and HF-QMC quantum impurity solvers.
Case 1: JASMINE component (atomic eigenvalue problem solver)
In the following, we will use an example to show you how to use API to control the jasmine code. When you want to compile your code, you have to ensure that japi.mod, libMM.a and libatomic.a are in correct PATH. Or else the compiler will complain that it can not find them.
NOTE:
You can find japi.mod, libMM.a, and libatomic.a in iqist/src/capi, iqist/src/base, and iqist/src/tools/jasmine directories, respectively.
(1) Import API support
use japi
(2) Create T_jasmine
type (T_jasmine) :: I_solver ! define I_solver
...
I_solver%ibasis = 1 ! setup I_solver
I_solver%ictqmc = 1
...
I_solver%Uc = 4.0
I_solver%Uv = 4.0
I_solver%Jz = 0.0
I_solver%Js = 0.0
I_solver%Jp = 0.0
...
NOTE:
Every parameter for the atomic eigenvalue problem solver must be initialized here.
(3) Initialize the atomic eigenvalue problem solver
call cat_init_atomic(I_solver)
(4) Start the atomic eigenvalue problem solver
call cat_exec_atomic()
(5) Close the atomic eigenvalue problem solver
call cat_stop_atomic()
(6) Access the computational results
You have to write your Fortran codes manually to access the results.
Case 2: CT-HYB quantum impurity solvers
In the following, we will use the AZALEA code as an example to show how to use API to control it. When you want to compile your code, you have to ensure that capi.mod, libMM.a and libctqmc.a are in correct PATH. Or else the compiler will complain that it can not find them.
NOTE:
You can find capi.mod, libMM.a, and libctqmc.a in iqist/src/capi, iqist/src/base, and iqist/src/ctqmc/azalea directories, respectively.
(1) Import API support
use capi
(2) create T_mpi
type (T_mpi) :: I_mpi ! define I_mpi
...
call mp_init() ! init mpi environment
call mp_comm_rank(I_mpi%myid) ! init I_mpi structure
call mp_comm_size(I_mpi%nprocs) ! init I_mpi structure
NOTE:
The above codes need MPI support. Namely, you have to import the MPI support explicitly, such as,
use mmpi ! import mpi support
(3) Create T_segment_azalea
type (T_segment_azalea) :: I_solver ! define I_solver
...
I_solver%isscf = 1 ! setup I_solver
I_solver%issun = 2
I_solver%isspn = 1
I_solver%isbin = 1
I_solver%nband = 1
I_solver%nspin = 2
I_solver%norbs = 2
I_solver%ncfgs = 4
I_solver%niter = 20
I_solver%mkink = 1024
I_solver%mfreq = 8193
I_solver%nfreq = 128
I_solver%ntime = 1024
I_solver%nflip = 10000
I_solver%ntherm = 20000
I_solver%nsweep = 20000000
I_solver%nwrite = 2000000
I_solver%nclean = 20000
I_solver%nmonte = 100
I_solver%ncarlo = 100
I_solver%U = 4.0
I_solver%Uc = 4.0
I_solver%Uv = 4.0
I_solver%Jz = 0.0
I_solver%Js = 0.0
I_solver%Jp = 0.0
I_solver%mune = 2.0
I_solver%beta = 10.0
I_solver%part = 0.50
I_solver%alpha = 0.50
NOTE:
If you want to use the other solvers, instead of the AZALEA code, please choose suitable solver type.
Every parameter for quantum impurity solver must be initialized here, or else the solver will not work properly.
(4) Initialize the CT-HYB quantum impurity solver
call cat_init_ctqmc(I_mpi, I_solver)
(5) Setup hybf, symm, eimp, and ktau
For examples:
integer :: size_t
complex(dp) :: hybf(size_t)
...
call cat_set_hybf(size_t, hybf) ! setup hybridization function: hybf
NOTE:
This step is optional, because the CT-HYB quantum impurity solvers will provide default values for hybf, symm, eimp, and ktau or read them from external disk files.
(6) Start the CT-HYB quantum impurity solver
call cat_exec_ctqmc(i)
Here is the current iteration number.
(7) Retrieve the calculated results
Through this API, the user can only access the sigf (i.e., self-energy function ), grnf (i.e., impurity Green's function ), nmat (i.e., impurity occupation number ), and nnmat (i.e., double occupation number ). As for the other physical observables, the user should parse the other output files generated by iQIST.
integer :: size_t
complex(dp) :: grnf(size_t)
call cat_get_grnf(size_t, grnf)
(8) Close the CT-HYB quantum impurity solver
call cat_stop_ctqmc()
(9) Finalize the MPI environment
call mp_barrier()
call mp_finalize()
NOTE:
This step is also optional.
Case 3: HF-QMC quantum impurity solver
In the following, we will use the DAISY code as an example to show how to use API to control it. When you want to compile your code, you have to ensure that dapi.mod, libMM.a and libhfqmc.a are in correct PATH. Or else the compiler will complain that it can not find them.
NOTE:
You can find dapi.mod, libMM.a, and libhfqmc.a in iqist/src/capi, iqist/src/base, and iqist/src/hfqmc/daisy directories, respectively.
(1) Import API support
use dapi
(2) Create T_mpi
type (T_mpi) :: I_mpi ! define I_mpi
...
call mp_init() ! init mpi environment
call mp_comm_rank(I_mpi%myid) ! init I_mpi structure
call mp_comm_size(I_mpi%nprocs) ! init I_mpi structure
NOTE:
The above codes need MPI support. Namely, you have to import the MPI support explicitly, such as,
use mmpi ! import mpi support
(3) Create T_daisy
type (T_daisy) :: I_solver ! define I_solver
...
I_solver%isscf = 1 ! setup I_solver
I_solver%issun = 2
I_solver%isspn = 1
I_solver%isbin = 1
I_solver%nband = 1
I_solver%nspin = 2
I_solver%norbs = 2
I_solver%niter = 20
I_solver%mstep = 16
I_solver%mfreq = 8193
I_solver%nsing = 1
I_solver%ntime = 128
I_solver%ntherm = 100
I_solver%nsweep = 240000
I_solver%nclean = 100
I_solver%ncarlo = 10
I_solver%Uc = 4.0
I_solver%Jz = 0.0
I_solver%mune = 2.0
I_solver%beta = 10.0
I_solver%part = 0.50
I_solver%alpha = 0.50
NOTE:
Every parameter for quantum impurity solver must be initialized here, or else the solver will not work properly.
(4) Init the HF-QMC quantum impurity solver
call cat_init_hfqmc(I_mpi, I_solver)
(5) Setup wssf, symm, eimp, and ktau
For examples:
integer :: size_t
complex(dp) :: wssf(size_t)
...
call cat_set_wssf(size_t, wssf) ! setup bath weiss's function: wssf
NOTE:
This step is optional, because the HF-QMC quantum impurity solver will provide default values for wssf, symm, eimp, and ktau or read them from external disk files.
(6) Start the HF-QMC quantum impurity solver
call cat_exec_hfqmc(i)
Here is the current iteration number.
(7) Retrieve the calculated results
Through this API, the user can only access the sigf (i.e., self-energy function ), grnf (i.e., impurity Green's function ), nmat (i.e., impurity occupation number ), and nnmat (i.e., double occupation number ). As for the other physical observables, the user should parse the other output files generated by iQIST.
integer :: size_t
complex(dp) :: grnf(size_t)
call cat_get_grnf(size_t, grnf)
(8) Close the HF-QMC quantum impurity solver
call cat_stop_hfqmc()
(9) Finalize the MPI environment
call mp_barrier()
call mp_finalize()
NOTE:
This step is also optional.