3
0
mirror of https://github.com/triqs/dft_tools synced 2024-08-06 20:40:00 +02:00

finished SOC tutorial

This commit is contained in:
aichhorn 2020-05-05 10:29:13 +02:00
parent 4b09b97ce4
commit 0ce04d34a7
4 changed files with 170 additions and 45 deletions

View File

@ -0,0 +1,112 @@
# Import the modules:
from triqs_dft_tools.sumk_dft import *
from pytriqs.gf import *
from pytriqs.archive import HDFArchive
from pytriqs.operators.util import *
from pytriqs.operators.util.U_matrix import *
from triqs_cthyb import *
import pytriqs.utility.mpi as mpi
# Init the SumK class:
filename = 'Sr2MgOsO6_SOC.h5'
SK = SumkDFT(hdf_file=filename,use_dft_blocks=True)
# Find diagonal local basis set:
SK.calculate_diagonalization_matrix(prop_to_be_diagonal='eal',calc_in_solver_blocks=True)
###########################
# Now we pick the orbitals:
# BE CAREFUL: THIS NEEDS TO BE DONE PROPERLY
# AND IS DIFFERENT FORM CASE TO CASE!
SK.block_structure.pick_gf_struct_solver([{'ud_0': [0,1,2],'ud_1': [0,1,2]}])
###########################
# Now we set up the U matrix, first in cubic (Wien2k) convention:
U = 2.0
J = 0.2
U_sph = U_matrix(l=2, U_int=U, J_hund=J)
U_sph = np.kron(np.reshape(np.eye(2),(1,2,1,2)),np.kron(np.reshape(np.eye(2),(2,1,2,1)),U_sph))
U_mat = transform_U_matrix(U_sph, SK.T[0].conjugate())
# Now we set up the interaction Hamiltonian
h_sumk = h_int_slater(['ud'], range(10), U_mat, off_diag=True, complex=True)
# And convert it to the solver structure
h_int = SK.block_structure.convert_operator(h_sumk)
# Solver Init:
beta = 40.0
S = Solver(beta=beta, gf_struct=SK.block_structure.gf_struct_solver_list[0])
# Solver parameters:
p = {}
# solver
p["random_seed"] = 123 * mpi.rank + 567
p["length_cycle"] = 200
p["n_warmup_cycles"] = 100000
p["n_cycles"] = 3000000
# tail fit
p["perform_tail_fit"] = True
p["fit_max_moment"] = 4
p["fit_min_w"] = 4.0
p["fit_max_w"] = 8.0
# double counting correction:
dc_type = 0 # FLL
# DMFT loops:
n_loops = 1
#for first iteration, add the outout group:
if mpi.is_master_node():
with HDFArchive(filename) as ar:
if (not ar.is_group('dmft_output')):
ar.create_group('dmft_output')
for iteration_number in range(1,n_loops+1):
mpi.report("Iteration = %s"%iteration_number)
SK.set_Sigma([ S.Sigma_iw ]) # put Sigma into the SumK class
chemical_potential = SK.calc_mu( precision = 0.01 ) # find the chemical potential for given density
S.G_iw << SK.extract_G_loc()[0]
if (iteration_number==1):
# Put Hartree energy on Re Sigma
dm = S.G_iw.density()
SK.calc_dc(dm, U_interact = U, J_hund = J, orb = 0, use_dc_formula = dc_type)
S.Sigma_iw << SK.block_structure.convert_matrix(SK.dc_imp[0],space_from='sumk',space_to='solver')['ud_0'][0,0]
mpi.report("Orbital densities of local Green function :")
for s,gf in S.G_iw:
mpi.report("Orbital %s: %s"%(s,dm[s].real))
mpi.report("Total charge of Gloc : %.6f"%S.G_iw.total_density().real)
# Calculate new G0_iw to input into the solver:
S.G0_iw << S.Sigma_iw + inverse(S.G_iw)
S.G0_iw << inverse(S.G0_iw)
# Solve the impurity problem:
S.solve(h_int=h_int, **p)
# Solved. Now do post-solution stuff:
dm = S.G_iw.density()
mpi.report("Orbital densities of impurity Green function:")
for s,gf in S.G_iw:
mpi.report("Orbital %s: %s"%(s,dm[s].real))
mpi.report("Total charge of impurity problem : %.6f"%S.G_iw.total_density().real)
# Write the final Sigma and G to the hdf5 archive:
if mpi.is_master_node():
with HDFArchive(filename) as ar:
ar['dmft_output']['iterations'] = iteration_number
ar['dmft_output']['G_0'] = S.G0_iw
ar['dmft_output']['G_tau'] = S.G_tau
ar['dmft_output']['G_iw'] = S.G_iw
ar['dmft_output']['Sigma_iw_%s'%iteration_number] = S.Sigma_iw
# Set the new double counting:
dm = S.G_iw.density() # compute the density matrix of the impurity problem
SK.calc_dc(dm, U_interact = U, J_hund = J, orb = 0, use_dc_formula = dc_type)
# Save stuff into the user_data group of hdf5 archive in case of rerun:
SK.save(['chemical_potential','dc_imp','dc_energ'])

Binary file not shown.

After

Width:  |  Height:  |  Size: 26 KiB

View File

@ -136,7 +136,7 @@ In the last line we use the Wien2k convention to write the U matrix in the cubic
h_int = SK.block_structure.convert_operator(h_sumk) h_int = SK.block_structure.convert_operator(h_sumk)
h_int = h_int.real h_int = h_int.real
Note that we needed to set up the interaction matrix for the full set of five *d* orbitals. The :meth:`convert_operator` method then takes care of rotating and picking the relevant orbitals. In the last line above we made the Hamiltonian real, since we know it this case that there are only real numbers in the interaction Hamiltonian. Note that this is not generally the case! Note that we needed to set up the interaction Hamiltonian for the full set of five *d* orbitals. The :meth:`convert_operator` method then takes care of rotating and picking the relevant orbitals. In the last line above we made the Hamiltonian real, since we know it this case that there are only real numbers in the interaction Hamiltonian. Note that this is not generally the case!
Now we have the interaction Hamiltonian for the solver, which we set up next:: Now we have the interaction Hamiltonian for the solver, which we set up next::

View File

@ -69,8 +69,6 @@ This reads all the data, and stores everything that is necessary for the DMFT ca
The DMFT calculation The DMFT calculation
==================== ====================
[TBC]
Rotating the basis Rotating the basis
------------------ ------------------
@ -81,63 +79,62 @@ We first initialise the SumK class::
from triqs_dft_tools.sumk_dft import * from triqs_dft_tools.sumk_dft import *
SK = SumkDFT(hdf_file='Sr2MgOsO6_SOC.h5',use_dft_blocks=True) SK = SumkDFT(hdf_file='Sr2MgOsO6_SOC.h5',use_dft_blocks=True)
The flag *use_dft_blocks=True* determines, as usual, the smallest possible blocks with non-zero entries, and initialises them as *solver* block structure. In order to disentangle the :math:`d_{x^2-y^2}` and the :math:`d_{xy}` orbitals, and pick out only the partially filled one, we do a transformation into a basis where the local Hamiltonian is diagonal:: The flag *use_dft_blocks=True* determines, as usual, the smallest possible blocks with non-zero entries, and initialises them as *solver* block structure. We want to work in the basis that diagonalises the local Hamiltonian, which is often referred to as numerical j-basis. The transformation into this basis is calculated by::
mat = SK.calculate_diagonalization_matrix(prop_to_be_diagonal='eal',calc_in_solver_blocks=True) mat = SK.calculate_diagonalization_matrix(prop_to_be_diagonal='eal',calc_in_solver_blocks=True)
We can look at the diagonalisation matrix, it is:: This transformation is stored in the SK.block_structure class and can be used to transform Greens function or operators. The next step is actually not needed for a DMFT calculation, but it is good to do this check to see what the transformation does to the local Hamiltonian. Note that in case of SOC, although spin is not a good quantum number any more, there are two blocks of size 5x5, each corresponding to negative/positive :math:`m_j` values::
>>> print mat[0]['down']
[[ 1. +0.j 0. +0.j 0. +0.j 0. +0.j 0. +0.j ]
[ 0. +0.j -0.985+0.j 0. -0.173j 0. +0.j 0. +0.j ]
[ 0. +0.j 0.173+0.j 0. -0.985j 0. +0.j 0. +0.j ]
[ 0. +0.j 0. +0.j 0. +0.j 1. +0.j 0. +0.j ]
[ 0. +0.j 0. +0.j 0. +0.j 0. +0.j 1. +0.j ]]
>>>
This transformation is already stored in the SK.block_structure class. The next step is actually not needed for a DMFT calculation, but it is good to do this check to see what the transformation does to the local Hamiltonian. We can calculate it before rotation, rotate it, and look at the 2x2 block with and without off-diagonals::
eal = SK.eff_atomic_levels() eal = SK.eff_atomic_levels()
eal2 = SK.block_structure.convert_matrix(eal[0],space_from='sumk', space_to='solver') eal2 = SK.block_structure.convert_matrix(eal[0],space_from='sumk', space_to='solver')
print eal[0]['up'][1:3,1:3] # prints the 2x2 block with offiagonals print eal2['ud_0'].diagonal().real. # diagonal of the first 5x5 block
[[ 0.391+0.j -0. -0.815j] [0.071 0.131 0.608 4.572 5.128]
[-0. +0.815j 4.884-0.j ]]
print eal2['up_1'] # prints the 2x2 block after rotation print eal2['ud_1'].diagonal().real. # diagonal of the second 5x5 block
[[0.247-0.j 0. +0.j] [0.071 0.131 0.608 4.572 5.128]
[0. -0.j 5.028+0.j]]
So the local Hamiltonian has been diagonalised. From the other entries we can see that the *up_0* block and the [1,1] entry of the *up_1* block correspond to :math:`e_g`-like orbitals, and the others are the We see that the orbitals are ordered from low to high energy, which makes picking the orbitals around the Fermi level quite easy. We just take indices 0 to 2::
:math:`t_{2g}` orbitals that we want to keep. So we pick the according orbitals in the block structure::
SK.block_structure.pick_gf_struct_solver([{'up_1': [0],'up_2': [0],'up_3': [0],'down_1': [0],'down_2': [0],'down_3': [0]}]) SK.block_structure.pick_gf_struct_solver([{'ud_0': [0,1,2],'ud_1': [0,1,2]}])
We can now look at the final result:: We can now look at the final result::
print SK.block_structure.convert_matrix(eal[0],space_from='sumk',space_to='solver') eal3 = SK.block_structure.convert_matrix(eal[0],space_from='sumk', space_to='solver')
{'up_2': array([[0.156-0.j]]), 'up_3': array([[0.156-0.j]]), 'up_1': array([[0.247-0.j]]), 'down_3': array([[0.156-0.j]]), 'down_2': array([[0.156-0.j]]), 'down_1': array([[0.247-0.j]])}
print eal3['ud_0']
We see that we arrived at a structure with 3 orbitals per spin only, and blocks of size 1x1. [[ 0.071+0.j 0. +0.j -0. -0.j]
[-0. -0.j 0.131+0.j -0. -0.j]
[-0. +0.j 0. +0.j 0.608-0.j]]
print eal3['ud_1']
[[ 0.071+0.j -0. -0.j 0. +0.j]
[-0. +0.j 0.131-0.j 0. +0.j]
[ 0. -0.j 0. -0.j 0.608+0.j]]
We see that we arrived at a diagonal structure with two blocks of size 3x3, and we have picked the orbitals around the Fermi level.
The interaction Hamiltonian The interaction Hamiltonian
--------------------------- ---------------------------
We now set up the interaction Hamiltonian. Since we want to rotate the interaction matrix into the local basis, we are using the Slater convention for it:: We now set up the interaction Hamiltonian. Since we want to rotate the interaction matrix into the local basis, we are using the Slater convention for it. We use *l=2* for *d* orbitals. Also, for SOC calculations, we need to inflate the resulting matrix to size 10x10::
from pytriqs.operators.util import * from pytriqs.operators.util import *
from pytriqs.operators.util.U_matrix import * from pytriqs.operators.util.U_matrix import *
U = 2.0 U = 2.0
J = 0.2 J = 0.2
U_mat = U_matrix(l=2,U_int=U,J_hund=J,basis='other', T=SK.T[0].conjugate()) U_sph = U_matrix(l=2, U_int=U, J_hund=J)
U_sph = np.kron(np.reshape(np.eye(2),(1,2,1,2)),np.kron(np.reshape(np.eye(2),(2,1,2,1)),U_sph)) # inflating the matrix
U_mat = transform_U_matrix(U_sph, SK.T[0].conjugate())
In the last line we use the Wien2k convention to write the U matrix in the cubic harmonics. Next, we want to rotate it, and pick out only the relevant :math:`t_{2g}` orbitals. At this step we need to give the indices of the wanted orbitals in the *sumk* basis:: In the last line we use the Wien2k convention to write the U matrix in the cubic harmonics. Next, we want to set up a Hamiltonian and rotate it into the *solver* basis::
U_trans = transform_U_matrix(U_mat, SK.block_structure.transformation[0]['up'].conjugate()) h_sumk = h_int_slater(['ud'], range(2*(2*l+1)), U_mat, off_diag=True, complex=True)
indices_to_pick_sumk = [1,3,4] # these are the relevant indices in the sumk basis h_int = SK.block_structure.convert_operator(h_sumk)
U_red = subarray(U_trans,len(U_trans.shape)*[indices_to_pick_sumk])
h_int = h_int_slater(['up','down'], indices_to_pick_sumk, U_red, map_operator_structure=SK.block_structure.sumk_to_solver[0], off_diag=False, complex=False) Note that we needed to set up the interaction Hamiltonian first for the full set of *d* orbitals. The :meth:`convert_operator` method then takes care of rotating and picking the relevant orbitals.
Now we have the interaction Hamiltonian for the solver, which we set up next:: Now we have the interaction Hamiltonian for the solver, which we set up next::
@ -157,16 +154,16 @@ Now we have the interaction Hamiltonian for the solver, which we set up next::
# tail fit # tail fit
p["perform_tail_fit"] = True p["perform_tail_fit"] = True
p["fit_max_moment"] = 4 p["fit_max_moment"] = 4
p["fit_min_n"] = 40 p["fit_min_w"] = 4.0
p["fit_max_n"] = 100 p["fit_max_w"] = 8.0
The DMFT loop with automatic basis rotations The DMFT loop with automatic basis rotations
-------------------------------------------- --------------------------------------------
After these initialisation steps, the formal DMFT cycle is very similar to a calculation without these basis rotations, since these rotations are done automatically, once the :class:`BlockStructure` property *transformation* is set, see :ref:`basisrotation`. After these initialisation steps, the formal DMFT cycle is very similar to a calculation without SOC, since the rotations are done automatically, once the :class:`BlockStructure` property *transformation* is set, see :ref:`basisrotation`.
The DMFT loop itself looks very much the same as in :ref:`SrVO3`:: The DMFT loop itself looks very much the same as in :ref:`SrVO3` or :ref:`Sr2MgOsO6_noSOC`::
# double counting correction: # double counting correction:
dc_type = 0 # FLL dc_type = 0 # FLL
@ -177,7 +174,6 @@ The DMFT loop itself looks very much the same as in :ref:`SrVO3`::
mpi.report("Iteration = %s"%iteration_number) mpi.report("Iteration = %s"%iteration_number)
SK.symm_deg_gf(S.Sigma_iw) # symmetrizing Sigma
SK.set_Sigma([ S.Sigma_iw ]) # put Sigma into the SumK class SK.set_Sigma([ S.Sigma_iw ]) # put Sigma into the SumK class
chemical_potential = SK.calc_mu( precision = 0.01 ) # find the chemical potential for given density chemical_potential = SK.calc_mu( precision = 0.01 ) # find the chemical potential for given density
S.G_iw << SK.extract_G_loc()[0] S.G_iw << SK.extract_G_loc()[0]
@ -186,7 +182,7 @@ The DMFT loop itself looks very much the same as in :ref:`SrVO3`::
# Put Hartree energy on Re Sigma # Put Hartree energy on Re Sigma
dm = S.G_iw.density() dm = S.G_iw.density()
SK.calc_dc(dm, U_interact = U, J_hund = J, orb = 0, use_dc_formula = dc_type) SK.calc_dc(dm, U_interact = U, J_hund = J, orb = 0, use_dc_formula = dc_type)
S.Sigma_iw << SK.block_structure.convert_matrix(SK.dc_imp[0],space_from='sumk',space_to='solver')['up_1'][0,0] S.Sigma_iw << SK.block_structure.convert_matrix(SK.dc_imp[0],space_from='sumk',space_to='solver')['ud_0'][0,0]
# Calculate new G0_iw to input into the solver: # Calculate new G0_iw to input into the solver:
S.G0_iw << S.Sigma_iw + inverse(S.G_iw) S.G0_iw << S.Sigma_iw + inverse(S.G_iw)
@ -203,12 +199,29 @@ The DMFT loop itself looks very much the same as in :ref:`SrVO3`::
# Save stuff into the user_data group of hdf5 archive in case of rerun: # Save stuff into the user_data group of hdf5 archive in case of rerun:
SK.save(['chemical_potential','dc_imp','dc_energ']) SK.save(['chemical_potential','dc_imp','dc_energ'])
The only difference to the other example is in the initialisation of the real part of the self energy. We cannot just take an element of the *dc_imp* array, since this array is stored in the *sumk* structure. Therefore, we first need to transform this matrix into *solver* space, and then take the appropriate matrix element. After the first iteration (here done with 24e6 MC sweeps), you should get self energies like this: The only difference to the other example is in the initialisation of the real part of the self energy. We cannot just take an element of the *dc_imp* array, since this array is stored in the *sumk* structure. Therefore, we first need to transform this matrix into *solver* space, and then take the appropriate matrix element. After the first iteration (here done with 18e6 MC sweeps), you should get this output at the end of the run::
.. image:: images_scripts/Sr2MgOsO6_noSOC_Sigmas.png Total number of measures: 18000000
Average sign: (0.884535,-4.11253e-06)
Orbital densities of impurity Green function:
Orbital ud_0:
[[ 5.20045070e-01 -8.24863778e-10 5.95348202e-12]
[-8.24863778e-10 4.30734642e-01 -1.29359496e-03]
[ 5.95348202e-12 -1.29359496e-03 4.80477133e-02]]
Orbital ud_1:
[[ 5.24181422e-01 2.22991244e-09 -8.16290063e-10]
[ 2.22991244e-09 4.30431196e-01 2.19004569e-03]
[-8.16290063e-10 2.19004569e-03 4.77161009e-02]]
Total charge of impurity problem : 2.001156
We see that there is a small sign problem due to the off-diagonal elements in the hybridisation function. However, this sign problem is treatable, since we have rotated into the local basis where these off-diagonal elements are minimised.
The imaginary part of the self energy matrix of the *ud_0* block looks like this:
.. image:: images_scripts/Sr2MgOsO6_SOC_Sigmas.png
:width: 600 :width: 600
:align: center :align: center
The two :math:`d_{xz}` and :math:`d_{yz}` orbitals are degenerate (blocks *up_2* and *up_3*), whereas the :math:`d_{xy}`-like orbital is different. Plotted on the same scale, the off-diagonal elements are very small, only the *(1,2)* and *(2,1)* elements are visibly different from zero.
A complete python script for this tutorial, including some more input/output, is available (:download:`Sr2MgOsO6_noSOC.py <images_scripts/Sr2MgOsO6_noSOC.py>`). When running the script, you will encounter warnings during the transformation from the *sumk* to the *solver* basis. These warnings just reflect that the off-diagonal elements of the full Greens function are not zero at all frequencies, although the local Hamiltonian is. In that sense, we still do an approximation when restricting ourselves to the :math:`t_{2g}`-like orbitals. A complete python script for this tutorial, including some more input/output, is available (:download:`Sr2MgOsO6_SOC.py <images_scripts/Sr2MgOsO6_SOC.py>`). When running the script, you will encounter warnings during the transformation from the *sumk* to the *solver* basis. These warnings just reflect that the off-diagonal elements of the full Greens function are not zero at all frequencies, although the local Hamiltonian is. In that sense, we still do an approximation when restricting ourselves to the low-energy subset.