3
0
mirror of https://github.com/triqs/dft_tools synced 2024-06-24 14:12:21 +02:00

[doc] updates

This commit is contained in:
Priyanka Seth 2015-03-11 14:35:28 +01:00
parent 1a48a8afe1
commit 9f08ec1d80
4 changed files with 227 additions and 175 deletions

View File

@ -1,5 +1,7 @@
.. index:: DFT+DMFT calculation
.. highlight:: python
.. _DFTDMFTmain:
The DFT+DMFT calculation
@ -20,64 +22,93 @@ to get the local quantities used in DMFT. It is initialized by::
from pytriqs.applications.dft.sumk_dft import *
SK = SumkDFT(hdf_file = filename)
The only necessary parameter is the filename of the hdf5 archive. In addition, there are some optional parameters:
* `mu`: The chemical potential at initialization. This value is only used if no other value is found in the hdf5 archive. The default value is 0.0.
* `h_field`: External magnetic field. The default value is 0.0.
* `use_dft_blocks`: If true, the structure of the density matrix is analysed at initialisation, and non-zero matrix elements
are identified. The DMFT calculation is then restricted to these matrix elements, yielding a more efficient solution of the
local interaction problem. Degeneracies in orbital and spin space are also identified and stored for later use. The default value is `False`.
* `dft_data`, `symmcorr_data`, `parproj_data`, `symmpar_data`, `bands_data`: These string variables define the subgroups in the hdf5 archive,
where the corresponding information is stored. The default values are consistent with those in :ref:`interfacetowien`.
This is the reference for the function:
At initialisation, the necessary data is read from the hdf5 file. If a calculation is restarted based on a previous hdf5 file, information on
degenerate shells, the block structure of the density matrix, the chemical potential, and double counting correction is also read in.
.. function:: SumkDFT(hdf_file, h_field = 0.0, use_dft_blocks = False, dft_data = 'dft_input', symmcorr_data = 'dft_symmcorr_input', parproj_data = 'dft_parproj_input', symmpar_data = 'dft_symmpar_input', bands_data = 'dft_bands_input', transp_data = 'dft_transp_input', misc_data = 'dft_misc_input')
The parameters needed to initialise SumkDFT as as follows (with any default variables as shown above):
========================= ==================== ===========================================================================
Name Type Meaning
========================= ==================== ===========================================================================
hdf_file String Name of the main hdf5 file containing converted dft information.
h_field Scalar External magnetic field.
use_dft_blocks Boolean Analyse the structure of the density matrix at initialisation,
and identify non-zero matrix elements.
The DMFT calculation is then restricted to these matrix elements,
yielding a more efficient solution of the local interaction problem.
Degeneracies in orbital and spin space are also identified and stored for later use.
dft_data String hdf5 subgroup containing required DFT data.
symmcorr_data String hdf5 subgroup containing correlated shell symmetry data.
parproj_data String hdf5 subgroup containing partial projector data.
symmpar_data String hdf5 subgroup containing non-correlated shell symmetry data.
bands_data String hdf5 subgroup containing DFT band data.
transp_data String hdf5 subgroup containing optics/transport data.
misc_data String hdf5 subgroup containing miscellaneous DFT data.
========================= ==================== ===========================================================================
.. index:: Multiband solver
Setting up the Multi-Band Solver
--------------------------------
There is a module that helps setting up the multiband CTQMC solver. It is loaded and initialized by::
The next step is to setup the impurity solver.
For more details here, see the `CTHYB <http://ipht.cea.fr/triqs/1.2/applications/cthyb/>`_ documentation.
from pytriqs.applications.dft.solver_multiband import *
S = SolverMultiBand(beta, n_orb, gf_struct = SK.gf_struct_solver[0], map=SK.map[0])
The necessary parameters are the inverse temperature `beta`, the Coulomb interaction `U_interact`, the Hund's rule coupling `J_hund`,
and the number of orbitals `n_orb`. There are again several optional parameters that allow the tailoring of the local Hamiltonian to
specific needs. They are:
* `gf_struct`: The block structure of the local density matrix given in the format calculated by :class:`SumkDFT`.
* `map`: If `gf_struct` is given as parameter, `map` also must be given. This is the mapping from the block structure to a general
up/down structure.
The solver method is called later by this statement::
S.solve(U_interact,J_hund,use_spinflip=False,use_matrix=True,
l=2,T=None, dim_reps=None, irep=None, n_cycles =10000,
length_cycle=200,n_warmup_cycles=1000)
The parameters for the Coulomb interaction `U_interact` and the Hund's coupling `J_hund` are necessary input parameters. The rest are optional
parameters for which default values are set. Generally, they should be reset for the problem at hand. Here is a description of the parameters:
* `use_matrix`: If `True`, the interaction matrix is calculated from Slater integrals, which are computed from `U_interact` and
`J_hund`. Otherwise, a Kanamori representation is used. Attention: We define the intraorbital interaction as
`U_interact`, the interorbital interaction for opposite spins as `U_interact-2*J_hund`, and interorbital for equal spins as
`U_interact-3*J_hund`.
* `T`: The matrix that transforms the interaction matrix from spherical harmonics to a symmetry-adapted basis. Only effective for Slater
parametrisation, i.e. `use_matrix=True`.
* `l`: The orbital quantum number. Again, only effective for Slater parametrisation, i.e. `use_matrix=True`.
* `use_spinflip`: If `True`, the full rotationally-invariant interaction is used. Otherwise, only density-density terms are
kept in the local Hamiltonian.
* `dim_reps`: If only a subset of the full d-shell is used as correlated orbtials, one can specify here the dimensions of all the subspaces
of the d-shell, i.e. t2g and eg. Only effective for Slater parametrisation.
* `irep`: The index in the list `dim_reps` of the subset that is used. Only effective for Slater parametrisation.
* `n_cycles`: Number of CTQMC cycles (a sequence of moves followed by a measurement) per core. The default value of 10000 is the minimum, and generally should be increased.
* `length_cycle`: Number of CTQMC moves per one cycle.
* `n_warmup_cycles`: Number of initial CTQMC cycles before measurements start. Usually of order of 10000, sometimes needs to be increased significantly.
Most of above parameters can be taken directly from the :class:`SumkDFT` class, without defining them by hand. We will see a specific example
at the end of this tutorial.
.. This is initialised as follows::
..
.. from pytriqs.applications.impurity_solvers.cthyb import *
.. beta = 40.0
.. gf_struct = SK.gf_struct_solver[0]
.. S = Solver(beta=beta, gf_struct=gf_struct)
..
.. The solver method is called later by this statement::
..
.. S.solve(h_loc=h_loc, **p)
..
.. where `p` represents the solve parameters.
..
.. There is a module that helps setting up the multiband CTQMC solver. It is loaded and initialized by::
..
.. from pytriqs.applications.dft.solver_multiband import *
.. S = SolverMultiBand(beta, n_orb, gf_struct = SK.gf_struct_solver[0], map=SK.map[0])
..
.. The necessary parameters are the inverse temperature `beta`, the Coulomb interaction `U_interact`, the Hund's rule coupling `J_hund`,
.. and the number of orbitals `n_orb`. There are again several optional parameters that allow the tailoring of the local Hamiltonian to
.. specific needs. They are:
..
.. * `gf_struct`: The block structure of the local density matrix given in the format calculated by :class:`SumkDFT`.
.. * `map`: If `gf_struct` is given as parameter, `map` also must be given. This is the mapping from the block structure to a general
.. up/down structure.
..
.. The solver method is called later by this statement::
..
.. S.solve(U_interact,J_hund,use_spinflip=False,use_matrix=True,
.. l=2,T=None, dim_reps=None, irep=None, n_cycles =10000,
.. length_cycle=200,n_warmup_cycles=1000)
..
.. The parameters for the Coulomb interaction `U_interact` and the Hund's coupling `J_hund` are necessary input parameters. The rest are optional
.. parameters for which default values are set. Generally, they should be reset for the problem at hand. Here is a description of the parameters:
..
.. * `use_matrix`: If `True`, the interaction matrix is calculated from Slater integrals, which are computed from `U_interact` and
.. `J_hund`. Otherwise, a Kanamori representation is used. Attention: We define the intraorbital interaction as
.. `U_interact`, the interorbital interaction for opposite spins as `U_interact-2*J_hund`, and interorbital for equal spins as
.. `U_interact-3*J_hund`.
.. * `T`: The matrix that transforms the interaction matrix from spherical harmonics to a symmetry-adapted basis. Only effective for Slater
.. parametrisation, i.e. `use_matrix=True`.
.. * `l`: The orbital quantum number. Again, only effective for Slater parametrisation, i.e. `use_matrix=True`.
.. * `use_spinflip`: If `True`, the full rotationally-invariant interaction is used. Otherwise, only density-density terms are
.. kept in the local Hamiltonian.
.. * `dim_reps`: If only a subset of the full d-shell is used as correlated orbtials, one can specify here the dimensions of all the subspaces
.. of the d-shell, i.e. t2g and eg. Only effective for Slater parametrisation.
.. * `irep`: The index in the list `dim_reps` of the subset that is used. Only effective for Slater parametrisation.
.. * `n_cycles`: Number of CTQMC cycles (a sequence of moves followed by a measurement) per core. The default value of 10000 is the minimum, and generally should be increased.
.. * `length_cycle`: Number of CTQMC moves per one cycle.
.. * `n_warmup_cycles`: Number of initial CTQMC cycles before measurements start. Usually of order of 10000, sometimes needs to be increased significantly.
..
.. Most of above parameters can be taken directly from the :class:`SumkDFT` class, without defining them by hand. We will see a specific example
.. at the end of this tutorial.
.. index:: DFT+DMFT loop, one-shot calculation
@ -91,18 +122,16 @@ set up the loop over DMFT iterations and the self-consistency condition::
n_loops = 5
for iteration_number in range(n_loops) : # start the DMFT loop
SK.put_Sigma(Sigma_imp = [ S.Sigma ]) # Put self energy to the SumK class
chemical_potential = SK.calc_mu() # calculate the chemical potential for the given density
S.G << SK.extract_G_loc()[0] # extract the local Green function
S.G0 << inverse(S.Sigma + inverse(S.G)) # finally get G0, the input for the Solver
SK.put_Sigma(Sigma_imp = [ S.Sigma ]) # Put self energy to the SumK class
chemical_potential = SK.calc_mu() # calculate the chemical potential for the given density
S.G_iw << SK.extract_G_loc()[0] # extract the local Green function
S.G0_iw << inverse(S.Sigma_iw + inverse(S.G_iw)) # finally get G0, the input for the Solver
S.solve(U_interact,J_hund,use_spinflip=False,use_matrix=True, # now solve the impurity problem
l=2,T=None, dim_reps=None, irep=None, n_cycles =10000,
length_cycle=200,n_warmup_cycles=1000)
S.solve(h_loc=h_loc, **p) # now solve the impurity problem
dm = S.G.density() # Density matrix of the impurity problem
SK.calc_dc( dm, U_interact = U, J_hund = J, use_dc_formula = 0) # Set the double counting term
SK.save(['chemical_potential','dc_imp','dc_energ']) # Save data in the hdf5 archive
dm = S.G_iw.density() # Density matrix of the impurity problem
SK.calc_dc(dm, U_interact=U, J_hund=J, orb=0, use_dc_formula=dc_type) # Set the double counting term
SK.save(['chemical_potential','dc_imp','dc_energ']) # Save data in the hdf5 archive
These basic steps are enough to set up the basic DMFT Loop. For a detailed description of the :class:`SumkDFT` routines,
see the reference manual. After the self-consistency steps, the solution of the Anderson impurity problem is calculation by CTQMC.
@ -120,8 +149,8 @@ At the end of the calculation, we can save the Greens function and self energy i
import pytriqs.utility.mpi as mpi
if mpi.is_master_node():
ar = HDFArchive("YourDFTDMFTcalculation.h5",'w')
ar["G"] = S.G
ar["Sigma"] = S.Sigma
ar["G"] = S.G_iw
ar["Sigma"] = S.Sigma_iw
This is it!

View File

@ -9,9 +9,11 @@ First, we load the necessary modules::
from pytriqs.applications.dft.sumk_dft import *
from pytriqs.applications.dft.converters.wien2k_converter import *
from pytriqs.applications.dft.solver_multiband import *
from pytriqs.gf.local import *
from pytriqs.archive import *
from pytriqs.archive import HDFArchive
from pytriqs.operators.hamiltonians import *
from pytriqs.applications.impurity_solvers.cthyb import *
Then we define some parameters::
@ -19,20 +21,20 @@ Then we define some parameters::
U = 2.7
J = 0.65
beta = 40
loops = 10 # Number of DMFT sc-loops
mix = 0.8 # Mixing factor of Sigma after solution of the AIM
Delta_mix = 1.0 # Mixing factor of Delta as input for the AIM
loops = 10 # Number of DMFT sc-loops
sigma_mix = 0.8 # Mixing factor of Sigma after solution of the AIM
delta_mix = 1.0 # Mixing factor of Delta as input for the AIM
dc_type = 1 # DC type: 0 FLL, 1 Held, 2 AMF
use_blocks = True # use bloc structure from DFT input
use_matrix = False # True: Slater parameters, False: Kanamori parameters U+2J, U, U-J
use_spinflip = False # use the full rotational invariant interaction?
prec_mu = 0.0001
qmc_cycles = 20000
length_cycle = 200
warming_iterations = 2000
# Solver parameters
p = {}
p["length_cycle"] = 200
p["n_warmup_cycles"] = 2000
p["n_cycles"] = 20000
Most of these parameters are self-explaining. The first, `dft_filename`, gives the filename of the input files.
Most of these parameters are self-explanatory. The first, `dft_filename`, gives the filename of the input files.
The next step, as described in the previous section, is to convert the input files::
Converter = Wien2kConverter(filename=dft_filename, repacking=True)
@ -46,106 +48,119 @@ from scratch::
previous_runs = 0
previous_present = False
if mpi.is_master_node():
ar = HDFArchive(dft_filename+'.h5','a')
if 'iterations' in ar:
previous_present = True
previous_runs = ar['iterations']
del ar
f = HDFArchive(dft_filename+'.h5','a')
if 'dmft_output' in f:
ar = f['dmft_output']
if 'iterations' in ar:
previous_present = True
previous_runs = ar['iterations']
else:
f.create_group('dmft_output')
del f
previous_runs = mpi.bcast(previous_runs)
previous_present = mpi.bcast(previous_present)
# if previous runs are present, no need for recalculating the bloc structure:
calc_blocs = use_blocks and (not previous_present)
Now we can use all this information to initialise the :class:`SumkDFT` class::
SK=SumkDFT(hdf_file=dft_filename+'.h5',use_dft_blocks=calc_blocs)
SK = SumkDFT(hdf_file=dft_filename+'.h5',use_dft_blocks=use_blocks)
If there was a previous run, we know already about the block structure, and therefore `UseDFTBlocs` is set to `False`.
The next step is to initialise the Solver::
Norb = SK.corr_shells[0]['dim']
n_orb = SK.corr_shells[0]['dim']
l = SK.corr_shells[0]['l']
S = SolverMultiBand(beta=beta,n_orb=Norb,gf_struct=SK.gf_struct_solver[0],map=SK.map[0])
As we can see, many options of the solver are set by properties of the :class:`SumkDFT` class, so we don't have
to set them manually.
spin_names = ["up","down"]
orb_names = [i for i in range(n_orb)]
# Use GF structure determined by DFT blocks
gf_struct = SK.gf_struct_solver[0]
# Construct U matrix for density-density calculations
Umat, Upmat = U_matrix_kanamori(n_orb=n_orb, U_int=U, J_hund=J)
# Construct Hamiltonian and solver
h_loc = h_loc_density(spin_names, orb_names, map_operator_structure=SK.sumk_to_solver[0], U=Umat, Uprime=Upmat, H_dump="H.txt")
S = Solver(beta=beta, gf_struct=gf_struct)
If there are previous runs stored in the hdf5 archive, we can now load the self energy
of the last iteration::
if (previous_present):
if (mpi.is_master_node()):
ar = HDFArchive(dft_filename+'.h5','a')
S.Sigma << ar['SigmaImFreq']
del ar
S.Sigma = mpi.bcast(S.Sigma)
if previous_present:
if mpi.is_master_node():
S.Sigma_iw << HDFArchive(dft_filename+'.h5','a')['dmft_output']['Sigma_iw']
chemical_potential,dc_imp,dc_energ = SK.load(['chemical_potential','dc_imp','dc_energ'])
S.Sigma_iw << mpi.bcast(S.Sigma_iw)
SK.set_mu(chemical_potential)
SK.set_dc(dc_imp,dc_energ)
The last command is the broadcasting of the self energy from the master node to the slave nodes.
Now we can go to the definition of the self-consistency step. It consists again of the basic steps discussed in the
previous section, with some additional refinement::
The self-energy is broadcast from the master node to the slave nodes. Also, the
last saved chemical potential and double counting values are read in and set.
for iteration_number in range(1,loops+1) :
SK.symm_deg_gf(S.Sigma,orb=0) # symmetrise Sigma
SK.put_Sigma(Sigma_imp = [ S.Sigma ]) # put Sigma into the SumK class:
chemical_potential = SK.calc_mu( precision = prec_mu ) # find the chemical potential
S.G << SK.extract_G_loc()[0] # calculation of the local Green function
mpi.report("Total charge of Gloc : %.6f"%S.G.total_density())
if ((iteration_number==1)and(previous_present==False)):
# Init the DC term and the real part of Sigma, if no previous run was found:
dm = S.G.density()
SK.calc_dc( dm, U_interact = U, J_hund = J, orb = 0, use_dc_formula = dc_type)
S.Sigma << SK.dc_imp[0]['up'][0,0]
S.G0 << inverse(S.Sigma + inverse(S.G))
# Solve the impurity problem:
S.solve(U_interact=U,J_hund=J,use_spinflip=use_spinflip,use_matrix=use_matrix,
l=l,T=SK.T[0], dim_reps=SK.dim_reps[0], irep=2, n_cycles=qmc_cycles,
length_cycle=length_cycle,n_warmup_cycles=warming_iterations)
# solution done, do the post-processing:
mpi.report("Total charge of impurity problem : %.6f"%S.G.total_density())
S.Sigma <<(inverse(S.G0)-inverse(S.G))
# Solve the impurity problem:
S.solve(U_interact=U,J_hund=J,use_spinflip=use_spinflip,use_matrix=use_matrix,
l=l,T=SK.T[0], dim_reps=SK.dim_reps[0], irep=2, n_cycles=qmc_cycles,
length_cycle=length_cycle,n_warmup_cycles=warming_iterations)
# solution done, do the post-processing:
mpi.report("Total charge of impurity problem : %.6f"%S.G.total_density())
S.Sigma <<(inverse(S.G0)-inverse(S.G))
# Now mix Sigma and G with factor Mix, if wanted:
if ((iteration_number>1) or (previous_present)):
if (mpi.is_master_node()):
ar = HDFArchive(dft_filename+'.h5','a')
mpi.report("Mixing Sigma and G with factor %s"%mix)
S.Sigma << mix * S.Sigma + (1.0-mix) * ar['Sigma']
S.G << mix * S.G + (1.0-mix) * ar['GF']
del ar
S.G = mpi.bcast(S.G)
S.Sigma = mpi.bcast(S.Sigma)
# Write the final Sigma and G to the hdf5 archive:
if (mpi.is_master_node()):
ar = HDFArchive(dft_filename+'.h5','a')
ar['iterations'] = previous_runs + iteration_number
ar['Sigma'] = S.Sigma
ar['GF'] = S.G
del ar
Now we can go to the definition of the self-consistency step. It consists again
of the basic steps discussed in the previous section, with some additional
refinement::
# Now set new double counting:
dm = S.G.density()
SK.calc_dc( dm, U_interact = U, J_hund = J, orb = 0, use_dc_formula = dc_type)
for iteration_number in range(1,loops+1):
if mpi.is_master_node(): print "Iteration = ", iteration_number
#Save stuff:
SK.save(['chemical_potential','dc_imp','dc_energ'])
SK.symm_deg_gf(S.Sigma_iw,orb=0) # symmetrise Sigma
SK.put_Sigma(Sigma_imp = [ S.Sigma_iw ]) # put Sigma into the SumK class
chemical_potential = SK.calc_mu( precision = prec_mu ) # find the chemical potential for given density
S.G_iw << SK.extract_G_loc()[0] # calc the local Green function
mpi.report("Total charge of Gloc : %.6f"%S.G_iw.total_density())
# Init the DC term and the real part of Sigma, if no previous runs found:
if (iteration_number==1 and previous_present==False):
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.dc_imp[0]['up'][0,0]
# Calculate new G0_iw to input into the solver:
if mpi.is_master_node():
# We can do a mixing of Delta in order to stabilize the DMFT iterations:
S.G0_iw << S.Sigma_iw + inverse(S.G_iw)
ar = HDFArchive(dft_filename+'.h5','a')['dmft_output']
if (iteration_number>1 or previous_present):
mpi.report("Mixing input Delta with factor %s"%delta_mix)
Delta = (delta_mix * delta(S.G0_iw)) + (1.0-delta_mix) * ar['Delta_iw']
S.G0_iw << S.G0_iw + delta(S.G0_iw) - Delta
ar['Delta_iw'] = delta(S.G0_iw)
S.G0_iw << inverse(S.G0_iw)
del ar
S.G0_iw << mpi.bcast(S.G0_iw)
# Solve the impurity problem:
S.solve(h_loc=h_loc, **p)
# Solved. Now do post-processing:
mpi.report("Total charge of impurity problem : %.6f"%S.G_iw.total_density())
# Now mix Sigma and G with factor sigma_mix, if wanted:
if (iteration_number>1 or previous_present):
if mpi.is_master_node():
ar = HDFArchive(dft_filename+'.h5','a')['dmft_output']
mpi.report("Mixing Sigma and G with factor %s"%sigma_mix)
S.Sigma_iw << sigma_mix * S.Sigma_iw + (1.0-sigma_mix) * ar['Sigma_iw']
S.G_iw << sigma_mix * S.G_iw + (1.0-sigma_mix) * ar['G_iw']
del ar
S.G_iw << mpi.bcast(S.G_iw)
S.Sigma_iw << mpi.bcast(S.Sigma_iw)
# Write the final Sigma and G to the hdf5 archive:
if mpi.is_master_node():
ar = HDFArchive(dft_filename+'.h5','a')['dmft_output']
if previous_runs: iteration_number += previous_runs
ar['iterations'] = iteration_number
ar['G_0'] = S.G0_iw
ar['G_tau'] = S.G_tau
ar['G_iw'] = S.G_iw
ar['Sigma_iw'] = S.Sigma_iw
del ar
# 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 dft_output group of hdf5 archive in case of rerun:
SK.save(['chemical_potential','dc_imp','dc_energ'])
This is all we need for the DFT+DMFT calculation. At the end, all results are stored in the hdf5 output file.

View File

@ -2,13 +2,14 @@
The interface
=============
The basic function of the interface to the Wien2k program package is
to take the output of the program that constructs the projected local
orbitals (:program:`dmftproj`, for documentation see :download:`TutorialDmftproj.pdf <TutorialDmftproj.pdf>`), and to store all the necessary information into
an hdf5 file. This latter file is then used to do the DMFT calculation. The
reason for this structure is that this enables the user to have everything
that is necessary to reproduce the calculation in one single hdf5 archive.
The basic function of the interface to the Wien2k program package is to take
the output of the program that constructs the projected local orbitals
(:program:`dmftproj`, for documentation see
:download:`TutorialDmftproj.pdf <TutorialDmftproj.pdf>`),
and to store all the necessary information into an hdf5 file. This latter file
is then used to do the DMFT calculation. The reason for this structure is that
this enables the user to have everything that is necessary to reproduce the
calculation in one single hdf5 archive.
.. index:: Interface to Wien2k
@ -31,17 +32,20 @@ calculation for TiO, the :program:`Wien2k` naming convention is that all files a
:file:`TiO.*`, so you would give `filename = "TiO"`. The constructor opens
an hdf5 archive, named :file:`material_of_interest.h5`, where all the data is stored.
There are three optional parameters to the Constructor:
These are the parameters to the Constructor:
* `dft_subgrp`: We store all data in subgroups of the hdf5 archive. For the main data
that is needed for the DMFT loop, we use the subgroup specified by this optional parameter.
The default value `dft_input` is used as the subgroup name.
* `symmcorr_subgrp`: In this subgroup we store all the data for applying the symmetry
operations in the DMFT loop. The default value is `dft_symmcorr_input`.
* `repacking`: If true, and the hdf5 file already exists, the system command :program:`h5repack`
is invoked. This command ensures a minimal file size of the hdf5
file. The default value is `False`. If you wish to use this, ensure
that :program:`h5repack` is in your path variable!
========================= ============================ ===========================================================================
Name Type, Default Meaning
========================= ============================ ===========================================================================
filename String Material being studied, corresponding to the :program:`Wien2k` file names.
The constructor stores the data in the hdf5 archive :file:`material_of_interest.h5`.
dft_subgrp String, dft_input hdf5 subgroup containing required DFT data
symmcorr_subgrp String, dft_symmcorr_input hdf5 subgroup containing all necessary data to apply
the symmetry operations in the DMFT loop
repacking Boolean, False Does the hdf5 file already exist and should the :program:`h5repack` be
invoked to ensures a minimal archive file size?
Note that the :program:`h5repack` must be in your path variable!
========================= ============================ ===========================================================================
After initialising the interface module, we can now convert the input text files into the
hdf5 archive by::
@ -70,8 +74,12 @@ of :program:`Wien2k`, you have to use::
This reads the files :file:`material_of_interest.parproj` and :file:`material_of_interest.sympar`.
Again, there are two optional parameters
* `parproj_subgrp`: The subgroup for partial projectors data. The default value is `dft_parproj_input`.
* `symmpar_subgrp`: The subgroup for symmetry operations data. The default value is `dft_symmpar_input`.
========================= ============================ ===========================================================================
Name Type, Default Meaning
========================= ============================ ===========================================================================
parproj_subgrp String, dft_parproj_input hdf5 subgroup containing partial projectors data.
symmpar_subgrp String, dft_symmpar_input hdf5 subgroup containing symmetry operations data.
========================= ============================ ===========================================================================
Another routine of the class allows to read the input for plotting the momentum-resolved
spectral function. It is done by::

View File

@ -9,7 +9,7 @@ All the data is stored using the hdf5 standard, as described also in the documen
hdf5 data format
----------------
In order to be used with the DMFT routines, the following data needs to be provided in the hdf5 file. It contains a lot of information in order to perform DMFT calculations for all kinds of situations, e.g. d-p Hamiltonians, more than one correlated atomic shell, or using symmetry operations for the k-summation. We store all data in subgroups of the hdf5 arxive:
In order to be used with the DMFT routines, the following data needs to be provided in the hdf5 file. It contains a lot of information in order to perform DMFT calculations for all kinds of situations, e.g. d-p Hamiltonians, more than one correlated atomic shell, or using symmetry operations for the k-summation. We store all data in subgroups of the hdf5 archive:
:program:`Main data`: There needs to be one subgroup for the main data of the calculation. The default name of this group is `dft_input`. Its contents are