abio Package¶
abio
Package¶
This sub-package provides objects and methods to generate Input files.
abivars
Module¶
This module contains lookup table with the name of the ABINIT variables.
- abipy.abio.abivars.is_abiunit(s: str) bool [source]¶
True if string is one of the units supported by the ABINIT parser
- class abipy.abio.abivars.AbinitInputFile(filepath: str)[source]¶
Bases:
TextFile
,Has_Structure
,NotebookWriter
This object parses the Abinit input file, stores the variables in dict-like objects (Datasets) and build Structure objects from the input variables. Mainly used for inspecting the structure declared in the Abinit input file.
- classmethod from_string(string: str) AbinitInputFile [source]¶
Build the object from string.
- structure()[source]¶
The structure defined in the input file.
If the input file contains multiple datasets AND the datasets have different structures, this property returns None. In this case, one has to access the structure of the individual datasets. For example:
input.datasets[0].structure
gives the structure of the first dataset.
- yield_figs(**kwargs)[source]¶
This function generates a predefined list of matplotlib figures with minimal input from the user.
- class abipy.abio.abivars.AbinitInputParser[source]¶
Bases:
object
- verbose = 0¶
- abipy.abio.abivars.structure_from_abistruct_fmt(string)[source]¶
Parse geometrical information given in the structure:abivars format return Structure object
A typical input file in “structure:abivars” format looks like:
# MgB2 lattice structure. natom 3 acell 2*3.086 3.523 Angstrom rprim 0.866025403784439 0.5 0.0 -0.866025403784439 0.5 0.0 0.0 0.0 1.0 # Atomic positions xred_symbols 0.0 0.0 0.0 Mg 1/3 2/3 0.5 B 2/3 1/3 0.5 B
i.e. the zcucl, typat and ntypat vars are implicitly defined via the xred_symbols table.
abivars_db
Module¶
Database with the names of the input variables used in Abinit and in other main programs.
- abipy.abio.abivars_db.get_abinit_variables()[source]¶
Returns the database with the description of the ABINIT variables.
- abipy.abio.abivars_db.get_anaddb_variables()[source]¶
Returns the database with the description of the ANADDB variables.
- abipy.abio.abivars_db.docvar(varname, executable='abinit')[source]¶
Return the Variable object associated to this name.
decorators
Module¶
Decorators for AbinitInput or MultiDataset objects.
- exception abipy.abio.decorators.InputDecoratorError[source]¶
Bases:
Exception
Error class raised by
AbinitInputDecorator
.
- class abipy.abio.decorators.AbinitInputDecorator[source]¶
Bases:
MSONable
An AbinitInputDecorator adds new options to an existing
AbinitInput
or an existingMultiDataset
without altering its structure. This is an abstract Base class.Example
decorator = MyDecorator(arguments)
new_abinit_input = decorator(abinit_input) new_multidataset = decorator(multidataset)
Note that a decorator does not modify the object on which it acts.
Warning
Please avoid introducing decorators acting on the structure (in particular the lattice) since the initial input may use the initial structure to compute important variables. For instance, the list of k-points for band structure calculation depend on the bravais lattice and a decorator that changes it should recompute the path. This should not represent a serious limitation because it’s always possible to change the structure with its methods and then call the factory function without having to decorate an already existing object.
- Error¶
alias of
InputDecoratorError
- class abipy.abio.decorators.SpinDecorator(spinmode, kptopt_ifspinor=4)[source]¶
Bases:
AbinitInputDecorator
This decorator changes the spin polarization.
- classmethod from_dict(d: dict) SpinDecorator [source]¶
- Parameters:
d – Dict representation.
- Returns:
MSONable class.
- class abipy.abio.decorators.SmearingDecorator(smearing)[source]¶
Bases:
AbinitInputDecorator
This decorator changes the electronic smearing.
- classmethod from_dict(d: dict) SmearingDecorator [source]¶
- Parameters:
d – Dict representation.
- Returns:
MSONable class.
- class abipy.abio.decorators.XcDecorator(ixc: int)[source]¶
Bases:
AbinitInputDecorator
Change the exchange-correlation functional.
- classmethod from_dict(d: dict) XcDecorator [source]¶
- Parameters:
d – Dict representation.
- Returns:
MSONable class.
- class abipy.abio.decorators.LdaUDecorator(symbols_luj, usepawu=1, unit='eV')[source]¶
Bases:
AbinitInputDecorator
This decorator adds LDA+U parameters to an
AbinitInput
object.- classmethod from_dict(d: dict) LdaUDecorator [source]¶
- Parameters:
d – Dict representation.
- Returns:
MSONable class.
- class abipy.abio.decorators.LexxDecorator(symbols_lexx, exchmix=None)[source]¶
Bases:
AbinitInputDecorator
This decorator add local exact exchange to a
AbinitInput
object.- classmethod from_dict(d: dict) LexxDecorator [source]¶
- Parameters:
d – Dict representation.
- Returns:
MSONable class.
enums
Module¶
This module defines enumerators associated to important Abinit input variables
- class abipy.abio.enums.EnumMixin[source]¶
Bases:
object
Mixin for enums provides extra capabilities.
- class abipy.abio.enums.StrEnum(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]¶
- class abipy.abio.enums.RUNL(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]¶
-
Values of optdriver corresponding to the different run-levels. See defs_basis.F90
- GSTATE = 0¶
- RESPFN = 1¶
- SCREENING = 3¶
- SIGMA = 4¶
- NONLINEAR = 5¶
- GWR = 6¶
- EPH = 7¶
- WFK = 8¶
- RTTDDFT = 9¶
- GWLS = 66¶
- BSE = 99¶
- LONGWAVE = 10¶
- class abipy.abio.enums.WFK_TASK(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]¶
-
Integer flags defining the task to be performed in wfk_analyze. See defs_basis.F90
- NONE = 0¶
- FULLBZ = 1¶
- CLASSIFY = 2¶
- PAW_AEPSI = 3¶
- EINTERP = 4¶
- DDK = 5¶
- DDK_DIAGO = 6¶
- OPTICS_FULLBZ = 7¶
- KPTS_ERANGE = 8¶
- CHECK_SYMTAB = 9¶
- class abipy.abio.enums.GWR_TASK(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]¶
-
String flags defining the task to be performed in the GWR code.
- HDIAGO = 'HDIAGO'¶
- HDIAGO_FULL = 'HDIAGO_FULL'¶
- CC4S = 'CC4S'¶
- CC4S_FULL = 'CC4S_FULL'¶
- G0W0 = 'G0W0'¶
- G0V = 'G0V'¶
- EGEW = 'EGEW'¶
- EGW0 = 'EGW0'¶
- G0EW = 'G0EW'¶
- RPA_ENERGY = 'RPA_ENERGY'¶
- GAMMA_GW = 'GAMMA_GW'¶
- CHI0 = 'CHI0'¶
factories
Module¶
Factory functions for Abinit input files
- abipy.abio.factories.gs_input(structure: Structure, pseudos, kppa=None, ecut=None, pawecutdg=None, scf_nband=None, accuracy='normal', spin_mode='polarized', smearing='fermi_dirac:0.1 eV', charge=0.0, scf_algorithm=None) AbinitInput [source]¶
Returns an
abipy.abio.inputs.AbinitInput
for ground-state calculation.- Parameters:
structure –
abipy.core.structure.Structure
object.pseudos – List of filenames or list of
pymatgen.io.abinit.pseudos.Pseudo
objects orpymatgen.io.abinit.pseudos.PseudoTable
object.kppa – Defines the sampling used for the SCF run. Defaults to 1000 if not given.
ecut – cutoff energy in Ha (if None, ecut is initialized from the pseudos according to accuracy)
pawecutdg – cutoff energy in Ha for PAW double-grid (if None, pawecutdg is initialized from the pseudos according to accuracy)
scf_nband – Number of bands for SCF run. If scf_nband is None, nband is automatically initialized from the list of pseudos, the structure and the smearing option.
accuracy – Accuracy of the calculation.
spin_mode – Spin polarization.
smearing – Smearing technique.
charge – Electronic charge added to the unit cell.
scf_algorithm – Algorithm used for solving of the SCF cycle.
- abipy.abio.factories.ebands_input(structure: Structure, pseudos, kppa=None, nscf_nband=None, ndivsm=15, ecut=None, pawecutdg=None, scf_nband=None, accuracy='normal', spin_mode='polarized', smearing='fermi_dirac:0.1 eV', charge=0.0, scf_algorithm=None, dos_kppa=None) MultiDataset [source]¶
Returns a
abipy.abio.inputs.MultiDataset
object for band structure calculations.- Parameters:
structure –
abipy.core.structure.Structure
object.pseudos – List of filenames or list of
pymatgen.io.abinit.pseudos.Pseudo
objects orpymatgen.io.abinit.pseudos.PseudoTable
object.kppa – Defines the sampling used for the SCF run. Defaults to 1000 if not given.
nscf_nband – Number of bands included in the NSCF run. Set to scf_nband + 10 if None.
ndivsm – Number of divisions used to sample the smallest segment of the k-path. if 0, only the GS input is returned in multi[0].
ecut – cutoff energy in Ha (if None, ecut is initialized from the pseudos according to accuracy)
pawecutdg – cutoff energy in Ha for PAW double-grid (if None, pawecutdg is initialized from the pseudos according to accuracy)
scf_nband – Number of bands for SCF run. If scf_nband is None, nband is automatically initialized from the list of pseudos, the structure and the smearing option.
accuracy – Accuracy of the calculation.
spin_mode – Spin polarization.
smearing – Smearing technique.
charge – Electronic charge added to the unit cell.
scf_algorithm – Algorithm used for solving the SCF cycle.
dos_kppa – Scalar or List of integers with the number of k-points per atom to be used for the computation of the DOS (None if DOS is not wanted).
- abipy.abio.factories.phonons_from_gsinput(gs_inp, ph_ngqpt=None, qpoints=None, with_ddk=True, with_dde=True, with_bec=False, ph_tol=None, ddk_tol=None, dde_tol=None, wfq_tol=None, qpoints_to_skip=None, qptopt=1, manager=None)[source]¶
Returns a list of inputs in the form of a MultiDataset to perform phonon calculations, based on a ground state
abipy.abio.inputs.AbinitInput
. It will determine if WFQ files should be calculated for some q points and add the NSCF AbinitInputs to the set. The inputs have the following tags, according to their function: “ddk”, “dde”, “nscf”, “ph_q_pert”. All of them have the tag “phonon”.- Parameters:
gs_inp – an
abipy.abio.inputs.AbinitInput
representing a ground state calculation, likely the SCF performed to get the WFK.ph_ngqpt – a list of three integers representing the gamma centered q-point grid used for the calculation. If None and qpoint is None None, the ngkpt value present in the gs_input will be used. Incompatible with qpoints.
qpoints – a list of coordinates of q points in reduced coordinates for which the phonon perturbations will be calculated. Incompatible with ph_ngqpt.
with_ddk – If True, if Gamma is included in the list of qpoints it will add inputs for the calculations of the DDK.
with_dde – If True, if Gamma is included in the list of qpoints it will add inputs for the calculations of the DDE. Automatically sets with_ddk=True.
with_bec – If Truem if Gamma is included in the list of qpoints the DDE will be calculated in the same input as the phonons. This will allow to determine the BECs. Automatically sets with_ddk=True and with_dde=False.
ph_tol – a dictionary with a single key defining the type of tolerance used for the phonon calculations and its value. Default: {“tolvrs”: 1.0e-10}.
ddk_tol – a dictionary with a single key defining the type of tolerance used for the DDK calculations and its value. Default: {“tolwfr”: 1.0e-22}.
dde_tol – a dictionary with a single key defining the type of tolerance used for the DDE calculations and its value. Default: {“tolvrs”: 1.0e-10}.
wfq_tol – a dictionary with a single key defining the type of tolerance used for the NSCF calculations of the WFQ and its value. Default {“tolwfr”: 1.0e-22}.
qpoints_to_skip – a list of coordinates of q points in reduced coordinates that will be skipped. Useful when calculating multiple grids for the same system to avoid duplicate calculations. If a DDB needs to be extended with more q points use e.g. ddb.qpoints.to_array().
qptopt – Option for the generation of q-points. Default: 1
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- abipy.abio.factories.g0w0_with_ppmodel_inputs(structure, pseudos, kppa, nscf_nband, ecuteps, ecutsigx, ecut=None, pawecutdg=None, shifts=(0.0, 0.0, 0.0), accuracy='normal', spin_mode='polarized', smearing='fermi_dirac:0.1 eV', ppmodel='godby', charge=0.0, scf_algorithm=None, inclvkb=2, scr_nband=None, sigma_nband=None, gw_qprange=1) MultiDataset [source]¶
Returns a
abipy.abio.inputs.MultiDataset
object that performs G0W0 calculations with the plasmon pole approximation.- Parameters:
structure –
abipy.core.structure.Structure
object.pseudos – List of filenames or list of
pymatgen.io.abinit.pseudos.Pseudo
objects orpymatgen.io.abinit.pseudos.PseudoTable
object.kppa – Defines the sampling used for the SCF run.
nscf_nband – Number of bands included in the NSCF run.
ecuteps – Cutoff energy [Ha] for the screening matrix.
ecutsigx – Cutoff energy [Ha] for the exchange part of the self-energy.
ecut – cutoff energy in Ha (if None, ecut is initialized from the pseudos according to accuracy)
pawecutdg – cutoff energy in Ha for PAW double-grid (if None, pawecutdg is initialized from the pseudos according to accuracy)
shifts – Shifts for k-mesh.
accuracy – Accuracy of the calculation.
spin_mode – Spin polarization.
smearing – Smearing technique.
ppmodel – Plasmonpole technique.
charge – Electronic charge added to the unit cell.
scf_algorithm – Algorithm used for solving of the SCF cycle.
inclvkb – Treatment of the dipole matrix elements (see abinit variable).
scr_nband – Number of bands used to compute the screening (default is nscf_nband)
sigma_nband – Number of bands used to compute the self-energy (default is nscf_nband)
gw_qprange – Option for the automatic selection of k-points and bands for GW corrections. See Abinit docs for more detail. The default value makes the code compute the QP energies for all the point in the IBZ and one band above and one band below the Fermi level.
- abipy.abio.factories.g0w0_convergence_inputs(structure, pseudos, kppa, nscf_nband, ecuteps, ecutsigx, scf_nband, ecut, accuracy='normal', spin_mode='polarized', smearing='fermi_dirac:0.1 eV', response_models=None, charge=0.0, scf_algorithm=None, inclvkb=2, gw_qprange=1, gamma=True, nksmall=None, extra_abivars=None) MultiDataset [source]¶
Returns a
abipy.abio.inputs.MultiDataset
object to generate a G0W0 work for the given the material. See also [Setten2017].- Parameters:
structure –
abipy.core.structure.Structure
objectpseudos – List of
pymatgen.io.abinit.pseudos.Pseudo
objects.kppa – k points per reciprocal atom.
scf_nband – number of scf bands
ecut – ecut for all calcs that that are not ecut convergence cals at scf level
run. (scf Defines the sampling used for the SCF)
nscf_nband – a list of number of bands included in the screening and sigmaruns. The NSCF run will be done with the maximum.
ecuteps – list of Cutoff energy [Ha] for the screening matrix.
ecutsigx – Cutoff energy [Ha] for the exchange part of the self-energy.
accuracy – Accuracy of the calculation.
spin_mode – Spin polarization.
smearing – Smearing technique.
charge – Electronic charge added to the unit cell.
scf_algorithm – Algorithm used for solving of the SCF cycle.
inclvkb – Treatment of the dipole matrix elements (see abinit variable).
response_models – List of response models
gw_qprange – selectpr for the qpoint mesh
gamma – is true a gamma centered mesh is enforced
nksmall – Kpoint division for additional band and dos calculations
extra_abivars – Dictionary with extra variables passed to ABINIT for all tasks.
extra abivars that are provided with _s appended will be take as a list of values to be tested a scf level
- abipy.abio.factories.bse_with_mdf_inputs(structure: Structure, pseudos, scf_kppa, nscf_nband, nscf_ngkpt, nscf_shiftk, ecuteps, bs_loband, bs_nband, mbpt_sciss, mdf_epsinf, ecut=None, pawecutdg=None, exc_type='TDA', bs_algo='haydock', accuracy='normal', spin_mode='polarized', smearing='fermi_dirac:0.1 eV', charge=0.0, scf_algorithm=None) MultiDataset [source]¶
Returns a
abipy.abio.inputs.MultiDataset
object that performs a GS + NSCF + Bethe-Salpeter calculation. The self-energy corrections are approximated with the scissors operator. The screening is modeled with the model dielectric function.- Parameters:
structure –
abipy.core.structure.Structure
object.pseudos – List of filenames or list of
pymatgen.io.abinit.pseudos.Pseudo
objects orpymatgen.io.abinit.pseudos.PseudoTable
object.scf_kppa – Defines the sampling used for the SCF run.
nscf_nband – Number of bands included in the NSCF run.
nscf_ngkpt – Divisions of the k-mesh used for the NSCF and the BSE run.
nscf_shiftk – Shifts used for the NSCF and the BSE run.
ecuteps – Cutoff energy [Ha] for the screening matrix.
bs_loband – Index of the first occupied band included the e-h basis set (ABINIT convention i.e. first band starts at 1). Can be scalar or array of shape (nsppol,)
bs_nband – Highest band idex used for the construction of the e-h basis set.
mbpt_sciss – Scissor energy in Hartree.
mdf_epsinf – Value of the macroscopic dielectric function used in expression for the model dielectric function.
ecut – cutoff energy in Ha (if None, ecut is initialized from the pseudos according to accuracy)
pawecutdg – cutoff energy in Ha for PAW double-grid (if None, pawecutdg is initialized from the pseudos according to accuracy)
exc_type – Approximation used for the BSE Hamiltonian (Tamm-Dancoff or coupling).
bs_algo – Algorith for the computatio of the macroscopic dielectric function.
accuracy – Accuracy of the calculation.
spin_mode – Spin polarization.
smearing – Smearing technique.
charge – Electronic charge added to the unit cell.
scf_algorithm – Algorithm used for solving the SCF cycle.
- abipy.abio.factories.ion_ioncell_relax_input(structure, pseudos, kppa=None, nband=None, ecut=None, pawecutdg=None, accuracy='normal', spin_mode='polarized', smearing='fermi_dirac:0.1 eV', charge=0.0, scf_algorithm=None, shift_mode='Monkhorst-pack') MultiDataset [source]¶
Returns a
abipy.abio.inputs.MultiDataset
for a structural relaxation. The first dataset optmizes the atomic positions at fixed unit cell. The second datasets optimizes both ions and unit cell parameters.- Parameters:
structure –
abipy.core.structure.Structure
object.pseudos – List of filenames or list of
pymatgen.io.abinit.pseudos.Pseudo
objects orpymatgen.io.abinit.pseudos.PseudoTable
object.kppa – Defines the sampling used for the Brillouin zone.
nband – Number of bands included in the SCF run.
accuracy – Accuracy of the calculation.
spin_mode – Spin polarization.
smearing – Smearing technique.
charge – Electronic charge added to the unit cell.
scf_algorithm – Algorithm used for the solution of the SCF cycle.
- abipy.abio.factories.ion_ioncell_relax_and_ebands_input(structure, pseudos, kppa=None, nband=None, ecut=None, pawecutdg=None, accuracy='normal', spin_mode='polarized', smearing='fermi_dirac:0.1 eV', charge=0.0, scf_algorithm=None) MultiDataset [source]¶
Returns a
abipy.abio.inputs.MultiDataset
for a structural relaxation followed by a band structure run. The first dataset optimizes the atomic positions at fixed unit cell. The second datasets optimizes both ions and unit cell parameters. The other datasets perform a band structure calculation.Warning
Client code is responsible for propagating the relaxed structure obtained with the second dataset to the inputs used for the band structure calculation.
- Parameters:
structure –
abipy.core.structure.Structure
object.pseudos – List of filenames or list of
pymatgen.io.abinit.pseudos.Pseudo
objects orpymatgen.io.abinit.pseudos.PseudoTable
object.kppa – Defines the sampling used for the Brillouin zone.
nband – Number of bands included in the SCF run.
accuracy – Accuracy of the calculation.
spin_mode – Spin polarization.
smearing – Smearing technique.
charge – Electronic charge added to the unit cell.
scf_algorithm – Algorithm used for solving of the SCF cycle.
Returns:
abipy.abio.inputs.MultiDataset
object
- abipy.abio.factories.scf_phonons_inputs(structure, pseudos, kppa, ecut=None, pawecutdg=None, scf_nband=None, accuracy='normal', spin_mode='polarized', smearing='fermi_dirac:0.1 eV', charge=0.0, scf_algorithm=None, qptopt=1) list[AbinitInput] [source]¶
Returns a list of input files for performing phonon calculations. GS input + the input files for the phonon calculation.
- Parameters:
structure –
abipy.core.structure.Structure
object.pseudos – List of filenames or list of
pymatgen.io.abinit.pseudos.Pseudo
objects orpymatgen.io.abinit.pseudos.PseudoTable
object.kppa – Defines the sampling used for the SCF run.
ecut – cutoff energy in Ha (if None, ecut is initialized from the pseudos according to accuracy)
pawecutdg – cutoff energy in Ha for PAW double-grid (if None, pawecutdg is initialized from the pseudos according to accuracy)
scf_nband – Number of bands for SCF run. If scf_nband is None, nband is automatically initialized from the list of pseudos, the structure and the smearing option.
accuracy – Accuracy of the calculation.
spin_mode – Spin polarization.
smearing – Smearing technique.
charge – Electronic charge added to the unit cell.
scf_algorithm – Algorithm used for solving of the SCF cycle.
qptopt – Option for the generation of q-points. Default: 1
- abipy.abio.factories.piezo_elastic_inputs_from_gsinput(gs_inp, ddk_tol=None, rf_tol=None, ddk_split=False, rf_split=False, manager=None) MultiDataset [source]¶
Returns a
abipy.abio.inputs.MultiDataset
for performing elastic and piezoelectric constants calculations. GS input + the input files for the elastic and piezoelectric constants calculation.- Parameters:
gs_inp – Ground State input to build piezo elastic inputs from.
ddk_tol – Tolerance for the DDK calculation (i.e. {“tolwfr”: 1.0e-20}).
rf_tol – Tolerance for the Strain RF calculations (i.e. {“tolvrs”: 1.0e-12}).
ddk_split – Whether to split the DDK calculations.
rf_split – whether to split the RF calculations.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- abipy.abio.factories.scf_piezo_elastic_inputs(structure, pseudos, kppa, ecut=None, pawecutdg=None, scf_nband=None, accuracy='normal', spin_mode='polarized', smearing='fermi_dirac:0.1 eV', charge=0.0, scf_algorithm=None, ddk_tol=None, rf_tol=None, ddk_split=False, rf_split=False) MultiDataset [source]¶
Returns a
abipy.abio.inputs.MultiDataset
for performing elastic and piezoelectric constants calculations. GS input + the input files for the elastic and piezoelectric constants calculation.- Parameters:
structure –
abipy.core.structure.Structure
object.pseudos – List of filenames or list of
pymatgen.io.abinit.pseudos.Pseudo
objects orpymatgen.io.abinit.pseudos.PseudoTable
object.kppa – Defines the sampling used for the SCF run.
ecut – cutoff energy in Ha (if None, ecut is initialized from the pseudos according to accuracy)
pawecutdg – cutoff energy in Ha for PAW double-grid (if None, pawecutdg is initialized from the pseudos according to accuracy)
scf_nband – Number of bands for SCF run. If scf_nband is None, nband is automatically initialized from the list of pseudos, the structure and the smearing option.
accuracy – Accuracy of the calculation.
spin_mode – Spin polarization.
smearing – Smearing technique.
charge – Electronic charge added to the unit cell.
scf_algorithm – Algorithm used for solving of the SCF cycle.
ddk_tol – Tolerance for the Ddk calculation (i.e. {“tolwfr”: 1.0e-20}).
rf_tol – Tolerance for the Strain RF calculations (i.e. {“tolvrs”: 1.0e-12}).
ddk_split – Whether to split the ddk calculations.
rf_split – whether to split the RF calculations.
- abipy.abio.factories.scf_for_phonons(structure, pseudos, kppa=None, ecut=None, pawecutdg=None, nband=None, accuracy='normal', spin_mode='polarized', smearing='fermi_dirac:0.1 eV', charge=0.0, scf_algorithm=None, shift_mode='Symmetric') AbinitInput [source]¶
- abipy.abio.factories.ddkpert_from_gsinput(gs_input, ddk_pert, nband=None, use_symmetries=False, ddk_tol=None, manager=None) AbinitInput [source]¶
Returns an
abipy.abio.inputs.AbinitInput
to perform a DDK calculations for a specific perturbation based on a ground stateabipy.abio.inputs.AbinitInput
.- Parameters:
gs_input – an
abipy.abio.inputs.AbinitInput
representing a ground state calculation, likely the SCF performed to get the WFK.ddk_pert – dict with the Abinit variables defining the perturbation Example: {‘idir’: 1, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]},
use_symmetries – boolean that determines if the irreducible components of the perturbation are used. Default to False. (TODO: Should be implemented)
ddk_tol – a dictionary with a single key defining the type of tolerance used for the DDK calculations and its value. Default: {“tolvrs”: 1.0e-22}.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- abipy.abio.factories.ddepert_from_gsinput(gs_input, dde_pert, use_symmetries=True, dde_tol=None, manager=None) AbinitInput [source]¶
Returns an
abipy.abio.inputs.AbinitInput
to perform a DDE calculations for a specific perturbation based on a ground stateabipy.abio.inputs.AbinitInput
.- Parameters:
gs_input – an
abipy.abio.inputs.AbinitInput
representing a ground state calculation, likely the SCF performed to get the WFK.dde_pert – dict with the Abinit variables defining the perturbation Example: {‘idir’: 1, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]},
use_symmetries – boolean that determines if the irreducible components of the perturbation are used. Default to True. Should be set to False for nonlinear coefficients calculation.
dde_tol – a dictionary with a single key defining the type of tolerance used for the DDE calculations and its value. Default: {“tolvrs”: 1.0e-22}.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- abipy.abio.factories.dtepert_from_gsinput(gs_input, dte_pert, manager=None) AbinitInput [source]¶
Returns an
abipy.abio.inputs.AbinitInput
to perform a DTE calculations for a specific perturbation based on a ground stateabipy.abio.inputs.AbinitInput
.- Parameters:
gs_input – an
abipy.abio.inputs.AbinitInput
representing a ground state calculation, likely the SCF performed to get the WFK.dte_pert – dict with the Abinit variables defining the perturbation Example: {‘idir’: 1, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]},
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- abipy.abio.factories.dte_from_gsinput(gs_input, use_phonons=True, ph_tol=None, ddk_tol=None, dde_tol=None, skip_dte_permutations=False, manager=None) MultiDataset [source]¶
Returns a list of inputs in the form of a
abipy.abio.inputs.MultiDataset
to perform calculations of non-linear properties, based on a ground state AbinitInput.The inputs have the following tags, according to their function: “ddk”, “dde”, “ph_q_pert” and “dte”. All of them have the tag “dfpt”.
- Parameters:
gs_input – an
abipy.abio.inputs.AbinitInput
representing a ground state calculation, likely the SCF performed to get the WFK.use_phonons – determine wether the phonon perturbations at gamma should be included or not
ph_tol – a dictionary with a single key defining the type of tolerance used for the phonon calculations and its value. Default: {“tolvrs”: 1.0e-22}.
ddk_tol – a dictionary with a single key defining the type of tolerance used for the DDK calculations and its value. Default: {“tolwfr”: 1.0e-22}.
dde_tol – a dictionary with a single key defining the type of tolerance used for the DDE calculations and its value. Default: {“tolvrs”: 1.0e-22}.
skip_dte_permutations – Since the current version of abinit always performs all the permutations of the perturbations, even if only one is asked, if True avoids the creation of inputs that will produce duplicated outputs.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- abipy.abio.factories.dfpt_from_gsinput(gs_inp, ph_ngqpt=None, qpoints=None, do_ddk=True, do_dde=True, do_strain=True, do_dte=False, ph_tol=None, ddk_tol=None, dde_tol=None, wfq_tol=None, strain_tol=None, skip_dte_permutations=False, manager=None) MultiDataset [source]¶
Returns a list of inputs in the form of a MultiDataset to perform a set of calculations based on DFPT including phonons, elastic and non-linear properties. Requires a ground state
abipy.abio.inputs.AbinitInput
as a starting point.It will determine if WFQ files should be calculated for some q points and add the NSCF AbinitInputs to the set. The original input is included and the inputs have the following tags, according to their function: “scf”, “ddk”, “dde”, “nscf”, “ph_q_pert”, “strain”, “dte”, “dfpt”.
N.B. Currently (version 8.8.3) anaddb does not support a DDB containing both 2nd order derivatives with qpoints different from gamma AND 3rd order derivatives. The calculations could be run, but the global DDB will not be directly usable as is.
- Parameters:
gs_inp – an
abipy.abio.inputs.AbinitInput
representing a ground state calculation, likely the SCF performed to get the WFK.ph_ngqpt – a list of three integers representing the gamma centered q-point grid used for the calculation. If None and qpoint==None the ngkpt value present in the gs_input will be used. Incompatible with qpoints.
qpoints – a list of coordinates of q points in reduced coordinates for which the phonon perturbations will be calculated. Incompatible with ph_ngqpt.
do_ddk – If True, if Gamma is included in the list of qpoints it will add inputs for the calculations of the DDK.
do_dde – If True, if Gamma is included in the list of qpoints it will add inputs for the calculations of the DDE. Automatically sets with_ddk=True.
do_strain – If True inputs for the strain perturbations will be included.
do_dte – If True inputs for the non-linear perturbations will be included. The phonon non-linear perturbations will be included only if a phonon calculation at gamma is present. The caller is responsible for adding it. Automatically sets with_dde=True.
ph_tol – a dictionary with a single key defining the type of tolerance used for the phonon calculations and its value. Default: {“tolvrs”: 1.0e-10}.
ddk_tol – a dictionary with a single key defining the type of tolerance used for the DDK calculations and its value. Default: {“tolwfr”: 1.0e-22}.
dde_tol – a dictionary with a single key defining the type of tolerance used for the DDE calculations and its value. Default: {“tolvrs”: 1.0e-10}.
wfq_tol – a dictionary with a single key defining the type of tolerance used for the NSCF calculations of the WFQ and its value. Default {“tolwfr”: 1.0e-22}.
strain_tol – dictionary with a single key defining the type of tolerance used for the strain calculations of and its value. Default {“tolvrs”: 1.0e-12}.
skip_dte_permutations – Since the current version of abinit always performs all the permutations of the perturbations, even if only one is asked, if True avoids the creation of inputs that will produce duplicated outputs.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- abipy.abio.factories.minimal_scf_input(structure: Structure, pseudos) AbinitInput [source]¶
Provides an input for a calculation with the minimum possible requirements. Can be used to execute abinit with minimal requirements when needing files that are produced only after a full calculation completes. In general this will contain 1 kpt, 1 band, very low cutoff, no polarization, no smearing. Disables checks on primitive cell and symmetries. Even for large system it will require small memory allocations and few seconds to execute.
- Parameters:
structure –
abipy.core.structure.Structure
object.pseudos – List of filenames or list of
pymatgen.io.abinit.pseudos.Pseudo
objects orpymatgen.io.abinit.pseudos.PseudoTable
object.
inputs
Module¶
This module defines objects to facilitate the creation of ABINIT input files.
- class abipy.abio.inputs.AbstractInput[source]¶
Bases:
MutableMapping
Abstract class defining the methods that must be implemented by Input classes.
- write(filepath: str = 'run.abi', files_file: bool = False) None [source]¶
Write the input file to file
filepath
.
- make_targz(tarname: str = 'input.tar.gz', **kwargs) str [source]¶
Build targz file with the Abinit input and the associated pseudos. Return path to the targz file.
- set_vars(*args, **kwargs) dict [source]¶
Set the value of the variables. Accept also comment=”string” Return: dict with the variables added to the input.
- set_vars_ifnotin(*args, **kwargs) dict [source]¶
Set the value of the variables only if the variable is not already present. Return dict with the variables added to the input.
Example
input.set_vars(ecut=10, ionmov=3)
- pop_vars(keys: list[str]) dict [source]¶
Remove the variables listed in keys. Return dictionary with the variables that have been removed. Unlike
remove_vars
, no exception is raised if the variables are not in the input.- Parameters:
keys – string or list of strings with variable names.
Example
inp.pop_vars([“ionmov”, “optcell”, “ntime”, “dilatmx”])
- remove_vars(keys: Iterable[str] | str, strict: bool = True) dict [source]¶
Remove the variables listed in keys. Return dictionary with the variables that have been removed.
- Parameters:
keys – string or list of strings with variable names.
strict – If True, KeyError is raised if at least one variable is not present.
- abstract property vars: dict¶
Dictionary with the input variables. Used to implement the dict-like interface.
- generate(**kwargs)[source]¶
This function generates new inputs by replacing the variables specified in kwargs.
- Parameters:
kwargs – keyword arguments with the values used for each variable.
gs_inp = call_function_to_generate_initial_template() # To generate two input files with different values of ecut: for inp_ecut in gs_inp.generate(ecut=[10, 20]): print("do something with inp_ecut %s" % inp_ecut) # To generate four input files with all the possible combinations of ecut and nsppol: for inp_ecut in gs_inpt.generate(ecut=[10, 20], nsppol=[1, 2]): print("do something with inp_ecut %s" % inp_ecut)
- class abipy.abio.inputs.AbiAbstractInput[source]¶
Bases:
AbstractInput
Abstract class defining the methods that must be implemented by Input objects. associated to Abinit executables.
- add_abiobjects(*abi_objects) dict [source]¶
This function receive a list of
AbiVarable
objects and add the corresponding variables to the input.
- abstract abivalidate(workdir=None, manager=None)[source]¶
This method should invoke the executable associated to the input object. to test whether the input variables are correct and consistent. The executable is supposed to implement some sort of –dry-run option that invokes the parser to validate the input and exits immediately.
- Parameters:
workdir – Working directory of the fake task used to compute the IBZ. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- Returns:
retcode: Return code. 0 if OK. output_file: output file of the run. log_file: log file of the Abinit run, use log_file.read() to access its content. stderr_file: stderr file of the Abinit run. use stderr_file.read() to access its content. task: Task object
- Return type:
namedtuple
with the following attributes
- exception abipy.abio.inputs.AbinitInputError[source]¶
Bases:
Exception
Base error class for exceptions raised by
AbinitInput
.
- class abipy.abio.inputs.AbinitInput(structure, pseudos, pseudo_dir=None, comment=None, decorators=None, abi_args=None, abi_kwargs=None, tags=None, enforce_znucl=None, enforce_typat=None)[source]¶
Bases:
AbiAbstractInput
,MSONable
,Has_Structure
This object stores the ABINIT variables for a single dataset.
Inheritance Diagram
The name of methods that invoke the Abinit executable to extact dimensions and additional quantities such as list of q-points, perturbations, etc start with the abiget_ prefix.
- Error¶
alias of
AbinitInputError
- enforce_znucl_and_typat(znucl, typat)[source]¶
These arrays are used to enforce a particular value of znucl and typat when writing the Abinit input file This trick may be useful if the ordering of the types based on the sites in the structure does not match the one used to produce other input file such as DDB/DVDB. The ordering indeed changed …
znucl[ntypat] = typat[natom] = Fortran convention. Start to count from 1.
- variable_checksum() str [source]¶
Return string with sha1 value in hexadecimal format. This method is mainly used in unit tests to check the invariance of the input objects
- property is_multidataset: bool¶
Used to understand if we have an AbinitInput or a MultiDataset in polymorphic APIs.
- property is_input: bool¶
Used to understand if we have an AbinitInput or a MultiDataset in polymorphic APIs.
- classmethod from_dict(d: dict) AbinitInput [source]¶
JSON interface used in pymatgen for easier serialization.
- property runlevel¶
A set of strings defining the calculation type.
- property uses_ktimereversal: bool¶
True if time-reversal symmetry is used to generate k-points in the IBZ.
- to_string(sortmode='section', post=None, with_mnemonics=False, mode='text', with_structure=True, with_pseudos=True, exclude=None, verbose=0, files_file=False) str [source]¶
String representation.
- Parameters:
sortmode – “section” if variables should be grouped by sections. “a” for alphabetical order, None if no sorting is wanted.
with_mnemonics – True if mnemonics should be added.
mode – Either text or html if HTML output with links is wanted.
post – String that will be appended to the name of the variables Note that post is usually autodetected when we have multiple datatasets It is mainly used when we have an input file with a single dataset so that we can prevent the code from adding “1” to the name of the variables (In this case, indeed, Abinit complains if ndtset=1 is not specified and we don’t want ndtset=1 simply because the code will start to add _DS1_ to all the input and output files.
with_structure – False if section with structure variables should not be printed.
with_pseudos – False if JSON section with pseudo data should not be added.
exclude – List of variable names that should be ignored.
files_file – if True a string compatible with the presence of a files file will be generated. Otherwise all the required variables will be added to the string.
- property pseudos_abivars: dict¶
A dictionary with the abivars related to the pseudopotentials extracted from the path of the internal pseudopotentials. Empty if pseudopotential related variables are already defined.
- property structure: Structure¶
The
abipy.core.structure.Structure
object associated to this input.
- replicate(ndtset: int) MultiDataset [source]¶
Helper function to construct a
abipy.abio.inputs.MultiDataset
with ndtset datasets from anabipy.abio.inputs.AbinitInput
.
- set_cutoffs_for_accuracy(accuracy: str) dict [source]¶
Set the value of ecut and pawecutdg (if PAW) using the hints reported in the pseudos. Raises
AbinitInputError
if pseudos do not provide hints. In the case of PAW, pawecutdg is either taken from the hints or computed from the recommended value of ecut using a scaling factor that depends onaccuracy
.- Parameters:
accuracy
Return: Dictionary with variables.
- set_scf_nband_semicond(nsppol=1, nspinor=1, nspden=1, charge=0.0, spinat=None) dict [source]¶
Set electronic pararameters, smearing options and compute
nband
for a GS-SCF run. assuming a semiconductor. See set_scf_nband. Return dict with variables.
- set_scf_nband(nsppol: int, nspinor: int, nspden: int, occopt: int, tsmear: float, charge: float, spinat) dict [source]¶
Set electronic pararameters, smearing options and compute
nband
for a GS-SCF run. Return dict with variables.- Parameters:
nsppol – Number of spins.
nspinor – Number of spinor components.
nspden – Number of spin density components.
occopt – Occoputation option.
tsmear – Electronic smearing.
charge – Extra charge.
spinat – If None and nsppol 2, spinat is automatically computed.
- set_kmesh(ngkpt, shiftk, kptopt: int = 1) dict [source]¶
Set the variables for the sampling of the BZ.
- Parameters:
ngkpt – Monkhorst-Pack divisions
shiftk – List of shifts.
kptopt – Option for the generation of the mesh.
- set_autokmesh(nksmall: int, kptopt: int = 1) dict [source]¶
Set the variables (ngkpt, shift, kptopt) for the sampling of the BZ.
- Parameters:
nksmall – Number of k-points used to sample the smallest lattice vector.
kptopt – Option for the generation of the mesh.
- get_ngkpt_shiftk() tuple [source]¶
Return info on the k-point sampling from the input file, more specifically a tuple with nkgpt and shift. ngkpt is set to None if the BZ sampling cannot be described in terms of three divisions + one shift.
- set_phdos_qmesh(nqsmall, method='tetra', ph_qshift=(0, 0, 0)) dict [source]¶
Set the variables (ngkpt, shift, kptopt) for the computation of the Phonon DOS in Abinit. Remember that the Phdos is computed via Fourier interpolation so there’s no constraint of the q-mesh.
- Parameters:
nqsmall – Number of k-points used to sample the smallest lattice vector.
method – “gaussian” or “tetra”.
ph_qshift – Shift for the mesh.
- set_kpath(ndivsm: int, kptbounds=None, iscf=-2) dict [source]¶
Set the variables for the NSCF computation of the electronic band structure.
- Parameters:
ndivsm – if > 0, it’s the number of divisions for the smallest segment of the path (Abinit variable). if < 0, it’s interpreted as the pymatgen line_density parameter in which the number of points in the segment is proportional to its length. Typical value: -20. This option is the recommended one if the k-path contains two consecutive high symmetry k-points that are very close as ndivsm > 0 may produce a very large number of wavevectors.
kptbounds – k-points defining the path in k-space. If None, we use the default high-symmetry k-path defined in the pymatgen database.
- set_qpath(ndivsm: int, qptbounds=None) dict [source]¶
Set the variables for the computation of the phonon band structure and phonon linewidths in the EPH part.
- Parameters:
ndivsm – Number of divisions for the smallest segment in the q-path.
qptbounds – q-points defining the path in q-space. If None, we use the default high-symmetry q-path defined in the pymatgen database.
- set_kptgw(kptgw, bdgw) dict [source]¶
Set the variables (k-points, bands) and nkptgw for the computation of self-energy matrix elements.
- Parameters:
kptgw – List of k-points in reduced coordinates.
bdgw – Specifies the range of bands for the GW corrections. Accepts iterable that can be reshaped as (nkptgw, 2) or a tuple of two integers if the extrema are the same for each k-point.
- set_spin_mode(spin_mode) dict [source]¶
Set the variables used to the treat the spin degree of freedom. Return dictionary with the variables that have been removed.
- Parameters:
spin_mode –
SpinMode
object or string. Possible values for string are:polarized (-)
unpolarized (-)
afm (-)
spinor (-)
spinor_nomag (-)
- set_autospinat(default=0.6) dict [source]¶
Set the variable spinat for collinear calculations in the format (0, 0, m) with the value of m determined with the following order of preference:
If the site of the structure has a magmom setting, that is used.
If the species on the site has a spin setting, that is used.
If the species itself has a particular setting in the config file, that is used, e.g., Mn3+ may have a different magmom than Mn4+.
The element symbol itself is checked in the config file.
If there are no settings, the default value is used.
- set_spinat_from_symbols(symb2spinat: dict, default=(0, 0, 0)) dict [source]¶
Set spinat parameters from a dictionary mapping chemical simbol to spinat value. If an element in the structure is not present in symb2luj, default is used.
Example
symb2spinat = {“Eu”: [0, 0, 7]} inp.set_spinat_from_symbols(symb2spinat)
- set_usepawu(usepawu, symb2luj, units='eV') dict [source]¶
Set DFT+U parameters for PAW calculations.
- Parameters:
usepawu – Option specifying the DFT+U flavor (Abinit input variable).
symb2luj – Dictionary mapping chemical symbol to the values of lpawu, upawu and jpawu. If an element in the structure is not present in symb2luj, the U+J term is automatically disabled for this element. In other words, only the element on which U+J should be applied must be speficied.
units – Energy units for U and J. Note that defaultis eV although ABINIT uses Hartree by default!
Example
symb2luj = {“Eu”: {“lpawu”: 3, “upawu”: 7, “jpawu”: 0.7} inp.set_luj(symb2luj)
- set_kmesh_nband_and_occ(ngkpt, shiftk, nsppol, occ1k_spin, nspinor=1, kptopt=1, occopt=2) dict [source]¶
Helper function to set occupancies when occopt == 2
- Parameters:
ngkpt – Divisions of the k-mesh.
shiftk – Shifts in reduced coordinates.
nsppol – Spin polarization (Abinit input variable).
occ1k_spin – List of nppols strings with the occupations for a single kpoint and the spin_up, spin_down channels if nsppol is 2.
nspinor
kptopt
occopt
- property pseudos: list[Pseudo]¶
List of
pymatgen.io.abinit.pseudos.Pseudo
objects.
- property num_valence_electrons: int¶
Number of valence electrons computed from the pseudos and the structure.
- property valence_electrons_per_atom: int¶
Number of valence electrons for each atom in the structure.
- linspace(varname, start, stop, num=50, endpoint=True) list[AbinitInput] [source]¶
Generate num input files by changing the value of variable varname over the interval [start, stop]. The endpoint of the interval can optionally be excluded.
Example
input_list = gs_template.linspace(“ecut”, start=10, stop=60)
- Parameters:
start – The starting value of the sequence.
stop – The end value of the sequence, unless endpoint is set to False. In this case, the sequence consists of all but the last of
ndtset + 1
evenly spaced samples, so that stop is excluded. Note that the step size changes when endpoint is False.num (int) – Number of samples to generate. Default is 50.
endpoint (bool) – optional. If True, stop is the last sample. Otherwise, it is not included. Default is True.
- arange(varname: str, start, stop, step) list[AbinitInput] [source]¶
Generate list of input files by changing the value of variable varname over the interval [start, stop]. Values are generated within the half-open interval
[start, stop)
(in other words, the interval including start but excluding stop).When using a non-integer step, such as 0.1, the results will often not be consistent. It is better to use
linspace
for these cases.Example
input_list = gs_template.linspace(“ecut”, start=10, stop=60, step=10)
- Parameters:
start – Start of interval. The interval includes this value. The default start value is 0.
stop – End of interval. The interval does not include this value, except in some cases where step is not an integer and floating point
step – Spacing between values. For any output out, this is the distance between two adjacent values,
out[i+1] - out[i]
. The default step size is 1. If step is specified, start must also be given.
- product(*items) list[AbinitInput] [source]¶
Cartesian product of input iterables. Equivalent to nested for-loops.
inp.product("ngkpt", "tsmear", [[2, 2, 2], [4, 4, 4]], [0.1, 0.2, 0.3])
- new_with_vars(*args, **kwargs) AbinitInput [source]¶
Return a new input with the given variables.
Example
new = input.new_with_vars(ecut=20, tsmear=0.04)
- new_with_structure(new_structure, scdims=None, verbose=1, **abi_vars) AbinitInput [source]¶
Return a new
abipy.abio.inputs.AbinitInput
with a different structure. See notes below for the constraints that must be fulfilled by the new structure- Parameters:
new_structure – Parameters defining the crystalline structure. Accepts
abipy.core.structure.Structure
object file with structure (CIF, netcdf file, …) or dictionary with ABINIT geo variables.scdims – 3 integer giving with the number of cells in the supercell along the three reduced directions. Must be used when new_structure represents a supercell of the initial structure defined in the input file.
verbose – Verbosity level.
abi_vars – Abinit variables added to the Abinit input.
Warning
If
scdims
is None (i.e. no supercell), the two structures must have the same value of natom and typat, they can only differ at the level of the lattice and of the atomic positions. When structure represents a supercell, scdims must be coherent with the new_structure passed as argument.
- new_with_decorators(decorators) AbinitInput [source]¶
This function receives a list of
AbinitInputDecorator
objects or just a single object, applies the decorators to the input and returns a newabipy.abio.inputs.AbinitInput
object. self is not changed.
- news_varname_values(varname_values) list[AbinitInput] [source]¶
Return list of new inputs by adding the variable and values specified in varname_values. Useful for performing convergence studies. Two cases are possible:
Generate e.g. 3 inputs with different number of nband:
varname_values = (“nband”, [8, 12, 14])
- for new_inp in input.news_with_varname_values(varname_values):
print(new_inp)
Take Cartesian product of two or multiple variables:
- var_vals = [
(“nband”, [8, 12]), (“ecut”, [4, 8]),
]
- for new_inp in input.news_with_varname_values(varname_values):
print(new_inp)
- pop_tolerances() dict [source]¶
Remove all the tolerance variables present in the input. Return dictionary with the variables that have been removed.
- pop_irdvars() dict [source]¶
Remove all the ird* variables present in the input. Return dictionary with the variables that have been removed.
- property scf_tolvar: tuple[str, float]¶
Returns the tolerance variable and value relative to the SCF convergence. If more than one is present, an error is raised.
- make_ebands_input(ndivsm=15, tolwfr=1e-20, nscf_nband=None, nb_extra=10, nbdbuf=None, nstep=100) AbinitInput [source]¶
Generate an input file for a NSCF band structure calculation along k-path from a GS SCF input.
- Parameters:
ndivsm – if > 0, it’s the number of divisions for the smallest segment of the path (Abinit variable). if < 0, it’s interpreted as the pymatgen line_density parameter in which the number of points in the segment is proportional to its length. Typical value: -20. This option is the recommended one if the k-path contains two high symmetry k-points that are very close as ndivsm > 0 may produce a very large number of wavevectors.
tolwfr – Tolerance on residuals for NSCF calculation.
nscf_nband – Number of bands for NSCF calculation. If None, use nband + nb_extra
nb_extra – Extra bands to to be added to input nband if nscf_nband is None.
nbdbuf – Number of states in buffer
nstep – Max number of NSCF iterations.
- make_edos_input(ngkpt, shiftk=(0, 0, 0), tolwfr=1e-20, nscf_nband=None, nb_extra=10, nstep=100) AbinitInput [source]¶
Generate an input file for electron DOS calculation from a GS-SCF input.
- Parameters:
ngkpt – Number of divisions for the k-mesh.
shiftk – List of k-mesh shifts.
tolwfr – Tolerance on residuals for NSCF calculation
nscf_nband – Number of bands for NSCF calculation. If None, use nband + nb_extra
nb_extra – Extra bands to to be added to input nband if nscf_nband is None.
nstep – Max number of NSCF iterations.
- make_nscf_kptopt0_input(kpts, tolwfr=1e-20, iscf=-2) AbinitInput [source]¶
Build an input for NSCF calculation from a GS SCF one. Uses explicit list of k-points and kptopt 0.
- Parameters:
kpts – List of k-points in reduced coordinates.
tolwfr – Tolerance on residuals.
- make_dfpt_effmass_inputs(kpts, effmass_bands_f90, ngfft=None, tolwfr=1e-20, iscf=-2) MultiDataset [source]¶
Return a
abipy.abio.inputs.MultiDataset
object with 2 inputs for the calculation of effective masses with DFPT The first input in a standard NSCF run, the second input computes the effective masses.- Parameters:
kpts – List of k-points in reduced coordinates where effective masses are wanted.
effmas_bands_f90 – (nkpt, 2) array with band range for effmas computation. WARNING: Assumes Fortran convention with indices starting from 1.
ngfft – FFT divisions (3 integers). Used to enforce the same FFT mesh in the NSCF run as the one used for GS.
tolwfr – Tolerance on residuals.
- make_ph_inputs_qpoint(qpt, tolerance=None, prtwf=-1, prepgkk=0, manager=None) MultiDataset [source]¶
Builds and returns a
abipy.abio.inputs.MultiDataset
list of input files for the calculation of phonons at the given q-point qpt. This method should be called with an input the represents a GS run.- Parameters:
qpt – q-point in reduced coordinates.
tolerance – dict {varname: value} with the tolerance to be used in the DFPT run. Defaults to {“tolvrs”: 1.0e-10}.
prtwf – 1WFs are only needed for restarting 2nd order DFPT or non-linear response. Since these files are huge, we use prtwf -1 so that the 1WF file is produced only if the calculation is not converged so that AbiPy can restart it.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
Warning
The routine assumes the q-point is such that k + q belongs to the initial GS mesh so that the DFPT run can be started from the WFK file directly without having to generate intermediate WFQ files.
- make_ddkpert_input(perturbation, kptopt=2, only_vk=False, use_symmetries=False, tolerance=None, manager=None) AbinitInput [source]¶
Returns
abipy.abio.inputs.AbinitInput
for the calculation of an electric field perturbation. This function should be called with an input that represents a GS run and an electric field perturbation.- Parameters:
perturbation – dict with the Abinit variables defining the irreducible perturbation Example: {‘idir’: 1, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]},
kptopt – 2 to take into account time-reversal symmetry. Note that kptopt 1 is not available.
only_vk – If only matrix elements of the velocity operator are needed. First-order wavefunctions won’t be converged –> not usable for other DFPT calculations.
use_symmetries – boolean that determines if the irreducible components of the perturbation are used. Default to False. (TODO: Should be implemented)
tolerance – dict {varname: value} with the tolerance to be used in the DFPT run. Defaults to {“tolwfr”: 1.0e-22}.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- make_ddk_inputs(tolerance=None, kptopt=2, only_vk=False, manager=None) MultiDataset [source]¶
Return inputs for performing DDK calculations. This functions should be called with an input the represents a GS run.
- Parameters:
tolerance – dict {varname: value} with the tolerance to be used in the DFPT run. Defaults to {“tolwfr”: 1.0e-22}.
kptopt – 2 to take into account time-reversal symmetry. Note that kptopt 1 is not available.
only_vk – If only matrix elements of the velocity operator are needed. First-order wavefunctions won’t be converged –> not usable for other DFPT calculations.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
Return:
abipy.abio.inputs.MultiDataset
object for DFPT runs.
- make_dkdk_input(tolerance=None, kptopt=2, manager=None) AbinitInput [source]¶
Return an input for performing d2/dkdk calculations. This functions should be called with an input the represents a GS run.
- Parameters:
tolerance – dict {varname: value} with the tolerance to be used in the DFPT run. Defaults to {“tolwfr”: 1.0e-22}.
kptopt – 2 to take into account time-reversal symmetry. Note that kptopt 1 is not available.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
Return:
abipy.abio.inputs.AbinitInput
object.
- make_ddepert_input(perturbation, use_symmetries=True, tolerance=None, manager=None) AbinitInput [source]¶
Returns
abipy.abio.inputs.AbinitInput
for the calculation of an electric field perturbation. This function should be called with an input that represents a GS run and an electric field perturbation.- Parameters:
perturbation – dict with the Abinit variables defining the irreducible perturbation Example: {‘idir’: 1, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]},
use_symmetries – boolean that determines if the irreducible components of the perturbation are used. Default to True. Should be set to False for nonlinear coefficients calculation.
tolerance – dict {varname: value} with the tolerance to be used in the DFPT run. Defaults to {“tolvrs”: 1.0e-22}.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- make_dde_inputs(tolerance=None, use_symmetries=True, manager=None) MultiDataset [source]¶
Return
abipy.abio.inputs.MultiDataset
inputs for the calculation of electric field perturbations. This functions should be called with an input the represents a GS run.- Parameters:
tolerance – dict {varname: value} with the tolerance to be used in the DFPT run. Defaults to {“tolvrs”: 1.0e-22}.
use_symmetries – boolean that computes the irreducible components of the perturbation. Default to True. Should be set to False for nonlinear coefficients calculation.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- make_dtepert_input(perturbation, ixc=None, manager=None) AbinitInput [source]¶
Return
abipy.abio.inputs.AbinitInput
for DTE calculation for a given perturbation. This functions should be called with an input that represents a GS run.- Parameters:
perturbation – dict with the Abinit variables defining the irreducible perturbation Example: {‘idir’: 1, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]},
ixc – Value of ixc variable. Used to overwrite the default value read from pseudos.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- make_dte_inputs(phonon_pert=False, skip_permutations=False, ixc=None, manager=None) MultiDataset [source]¶
Return
abipy.abio.inputs.MultiDataset
inputs for DTE calculation. This functions should be called with an input that represents a GS run.- Parameters:
phonon_pert – is True also the phonon perturbations will be considered. Default False.
skip_permutations – Since the current version of abinit always performs all the permutations of the perturbations, even if only one is asked, if True avoids the creation of inputs that will produce duplicated outputs.
ixc – Value of ixc variable. Used to overwrite the default value read from pseudos.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- make_bec_inputs(tolerance=None, prepalw=0, manager=None) MultiDataset [source]¶
Return
abipy.abio.inputs.MultiDataset
inputs for the calculation of Born effective charges. This functions should be called with an input that represents a GS run.- Parameters:
tolerance – dict {varname: value} with the tolerance to be used in the DFPT run. Defaults to {“tolvrs”: 1.0e-10}.
prepalw – 1 to activate computation of all 3*natom perts. Used to prepare longwave-limit calculation.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- make_quad_input(tolerance=None, kptopt=2, nstep=100, manager=None) AbinitInput [source]¶
Return an
abipy.abio.inputs.AbinitInput
for the calculation of dynamical quadrupoles.. This function should be called with an input that represents a GS run.Note that only selected features are compatible with dynamical quadrupoles. Please consult <https://docs.abinit.org/topics/longwave/>.
- Parameters:
tolerance – dict {varname: value} with the tolerance to be used in the DFPT run. Defaults to {“tolvrs”: 1.0e-10}.
kptopt – 2 to take into account time-reversal symmetry.
nstep – Max number of SCF iterations. Since AbiPy is still not able to restart
calculation (a Quadrupole)
user. (we increase the value if not already provided by the)
- make_strain_perts_inputs(tolerance=None, phonon_pert=True, efield_pert=False, kptopt=2, prepalw=0, manager=None) MultiDataset [source]¶
Return
abipy.abio.inputs.MultiDataset
inputs for strain perturbation calculation. This function should be called with an input that represents a GS run.- Parameters:
tolerance – dict {varname: value} with the tolerance to be used in the DFPT run. Defaults to {“tolvrs”: 1.0e-12}.
phonon_pert – isf True also the phonon perturbations are considered. Default True.
efield_pert – if True, the electric field perturbation at gamma is included. Requires DDK files
kptopt – 2 to take into account time-reversal symmetry.
prepalw – 1 to activate computation of all 3*natom perts. Used to prepare calculation of flexoelectric tensor.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- make_wfk_kerange_inputs(sigma_kerange, sigma_ngkpt, kptopt=1, einterp=(1, 5, 0, 0)) MultiDataset [source]¶
Return a
abipy.abio.inputs.MultiDataset
with two inputs for performing a WFK calculation with the kerange trick. This method should be called with the input associated to the NSCF run that produces the WFK file used for the interpolation.- Parameters:
sigma_kerange – The energy window for the WFK generation.
sigma_ngkpt – The fine grid of kpt inside the sigma_kerange interval.
kptopt – Option for the generation of the mesh. Default: use all symmetries.
einterp – The interpolation used. By default it is a star-function interpolation.
- make_eph_transport_input(ddb_ngqpt, sigma_erange, tmesh, eph_ngqpt_fine=None, kptopt=1, mixprec=1, boxcutmin=1.1, ibte_prep=0, ibte_niter=200, ibte_abs_tol=0, eph_restart=1) AbinitInput [source]¶
Return an
abipy.abio.inputs.AbinitInput
to perform phonon-limited transport calculations. This method is usually called with the input associated to the NSCF run that produces the WFK file used to start the EPH run so that we can directly inherit the k-mesh- Parameters:
ddb_ngqpt – the coarse q-grid used to compute the DDB and the DVDB file.
sigma_erange – Energy window for k-states (see Abinit variable)
tmesh – The mesh of temperatures (in Kelvin)
eph_ngqpt_fine – the fine q-grid used for the Fourier interpolation.
kptopt – Option for the generation of the mesh. Default: use all symmetries.
mixprec – For the last task only, 1 is often used to make the EPH calculation faster. Note that the Abinit default is 0.
boxcutmin – For the last task only, 1.1 is often used to decrease memory and is faster over the Abinit default of 2.
ibte_prep – Set it to 1 to activate the iterative Boltzmann equation. Default is SERTA/MRTA.
ibte_niter – Max number of iterations to solve the Boltzmann equation.
ibte_abs_tol – Stopping criterion for the IBTE. 0 activatess the automatic determination of the threshold.
eph_restart – 1 to allow the restart of a run based on the SIGEPH.nc file.
- make_gwr_qprange_input(gwr_ntau, nband, ecuteps, gw_qprange=0, gwr_task=GWR_TASK.G0W0, **kwargs) AbinitInput [source]¶
Build and return an input file to compute QP corrections with the GWR code.
- Parameters:
gwr_ntau – Number of minimax points.
nband – Number of bands in Green’s function
ecuteps – Cutoff energy for chi0
gap. (gw_qprange = 0 to compute the QP corrections only for the fundamental and the direct)
gwr_task – String defining the GWR task
- abivalidate(workdir=None, manager=None)[source]¶
Run ABINIT in dry-run mode to validate the input file. No change done in the input file.
- Parameters:
workdir – Working directory of the temporary task. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- Returns:
retcode: Return code. 0 if OK. output_file: output file of the run. log_file: log file of the Abinit run, use log_file.read() to access the content. stderr_file: stderr file of the Abinit run. use stderr_file.read() to access the content. task: Task object
- Return type:
namedtuple with the following attributes
- abiget_dryrun_task(workdir, manager, **extra_vars)[source]¶
Run ABINIT in dry-run mode. Change some control variable to prevent Abinit from stopping. Return Task object.
- abiget_spacegroup(tolsym=None, retdict=False, workdir=None, manager=None, verbose=0)[source]¶
This function invokes Abinit to get the space group (as detected by Abinit, not by spglib) It should be called with an input file that contains all the mandatory variables required by ABINIT.
- Parameters:
tolsym – Abinit tolsym input variable. None corresponds to the default value.
retdict – True to return dictionary with space group information instead of Structure.
workdir – Working directory of the temporary task used to compute the IBZ. None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.verbose – Verbosity level.
- Returns:
abipy.core.structure.Structure
object with AbinitSpaceGroup obtained from the main output file if retdict is False else dict with e.g. {‘bravais’: ‘Bravais cF (face-center cubic)’, ‘spg_number’: 227, ‘spg_symbol’: ‘Fd-3m’}.
- abiget_dims_spginfo(workdir=None, manager=None, verbose=0) tuple[dict, dict] [source]¶
Return tuple of dictionaries with dimensions and spacegroup info.
- abiget_ibz(ngkpt=None, shiftk=None, kptopt=None, workdir=None, manager=None, verbose=0)[source]¶
This function computes the list of points in the IBZ with the corresponding weights. It should be called with an input file that contains all the mandatory variables required by ABINIT.
- Parameters:
ngkpt – Number of divisions for the k-mesh (default None i.e. use ngkpt from self)
shiftk – List of shifts (default None i.e. use shiftk from self)
kptopt – Option for k-point generation. If None, the value in self is used.
workdir – Working directory of the fake task used to compute the ibz. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.verbose – verbosity level.
- Returns:
points:
numpy.ndarray
with points in the IBZ in reduced coordinates. weights:numpy.ndarray
with the weights of the points.- Return type:
namedtuple with attributes
- abiget_scr_ibz(ngkpt=None, shiftk=None, kptopt=None, workdir=None, manager=None, verbose=0)[source]¶
This function computes the list of points in the IBZ expected by the SCREENING code It should be called with an input file that contains all the mandatory variables required by ABINIT. This is usually the NSCF input file used to generate the WFK file used by SCR/SIGMA
- Parameters:
ngkpt – Number of divisions for the k-mesh (default None i.e. use ngkpt from self)
shiftk – List of shifts (default None i.e. use shiftk from self)
kptopt – Option for k-point generation. If None, the value in self is used.
workdir – Working directory of the fake task used to compute the ibz. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.verbose – verbosity level.
- Returns:
points:
numpy.ndarray
with q-points in the IBZ in reduced coordinates. weights:numpy.ndarray
with the weights of the points.- Return type:
namedtuple with attributes
- abiget_irred_phperts(qpt=None, ngkpt=None, shiftk=None, kptopt=None, prepgkk=0, prepalw=0, workdir=None, manager=None)[source]¶
This function, computes the list of irreducible perturbations for DFPT. It should be called with an input file that contains all the mandatory variables required by ABINIT.
- Parameters:
qpt – qpoint of the phonon in reduced coordinates. Used to shift the k-mesh if qpt is not passed, self must already contain “qpt” otherwise an exception is raised.
ngkpt – Number of divisions for the k-mesh (default None i.e. use ngkpt from self)
shiftk – Shiftks (default None i.e. use shiftk from self)
kptopt – Option for k-point generation. If None, the value in self is used.
prepalw – 1 to activate computation of all 3*natom perts. Used to prepare longwave-limit calculation.
prepgkk – 1 to activate computation of all 3*natom perts (debugging option).
workdir – Working directory of the fake task used to compute the ibz. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- Returns:
List of dictionaries with the Abinit variables defining the irreducible perturbation
Example
- [{‘idir’: 1, ‘ipert’: 1, ‘qpt’: [0.25, 0.0, 0.0]},
{‘idir’: 2, ‘ipert’: 1, ‘qpt’: [0.25, 0.0, 0.0]}]
- abiget_irred_ddeperts(ngkpt=None, shiftk=None, kptopt=None, workdir=None, manager=None)[source]¶
This function, computes the list of irreducible perturbations for DFPT computations. It should be called with an input file that contains all the mandatory variables required by ABINIT.
- Parameters:
ngkpt – Number of divisions for the k-mesh (default None i.e. use ngkpt from self)
shiftk – Shiftks (default None i.e. use shiftk from self)
kptopt – Option for k-point generation. If None, the value in self is used.
workdir – Working directory of the fake task used to compute the ibz. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- Returns:
List of dictionaries with the Abinit variables defining the irreducible perturbation
Example
- [{‘idir’: 1, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]},
{‘idir’: 2, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]}]
- abiget_irred_dteperts(ngkpt=None, shiftk=None, kptopt=None, ixc=None, workdir=None, manager=None, phonon_pert=False)[source]¶
This function, computes the list of irreducible perturbations for DFPT. It should be called with an input file that contains all the mandatory variables required by ABINIT.
- Parameters:
ngkpt – Number of divisions for the k-mesh (default None i.e. use ngkpt from self)
shiftk – Shiftks (default None i.e. use shiftk from self)
kptopt – Option for k-point generation. If None, the value in self is used.
workdir – Working directory of the fake task used to compute the ibz. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.phonon_pert – if True also the phonon perturbations will be considered. Default False.
- Returns:
List of dictionaries with the Abinit variables defining the irreducible perturbation
Example
- [{‘idir’: 1, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]},
{‘idir’: 2, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]}]
- abiget_irred_strainperts(ngkpt=None, shiftk=None, kptopt=None, phonon_pert=True, efield_pert=False, prepalw=0, workdir=None, manager=None)[source]¶
This function, computes the list of irreducible perturbations for strain perturbations in DFPT. It should be called with an input file that contains all the mandatory variables required by ABINIT.
- Parameters:
ngkpt – Number of divisions for the k-mesh (default None i.e. use ngkpt from self)
shiftk – Shiftks (default None i.e. use shiftk from self)
kptopt – Option for k-point generation. If None, the value in self is used.
phonon_pert – if True the phonon perturbation at gamma is included.
efield_pert – if True, the electric field perturbation at gamma is included. Requires DDK files
prepalw – 1 to activate computation of all 3*natom perts. Used to prepare calculation of flexoelectric tensor.
workdir – Working directory of the fake task used to compute the ibz. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- Returns:
List of dictionaries with the Abinit variables defining the irreducible perturbation
Example
- [{‘idir’: 1, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]},
{‘idir’: 2, ‘ipert’: 4, ‘qpt’: [0.0, 0.0, 0.0]}]
- pop_par_vars(all=False) dict [source]¶
Remove all the variables associated to parallelism from the input file. Useful in case of a restart when we need to remove the parallel variables before rerunning autoparal
- abiget_autoparal_pconfs(max_ncpus, autoparal=1, workdir=None, manager=None, verbose=0)[source]¶
Get all the possible configurations up to
max_ncpus
. Return list of parallel configurations.
- add_tags(tags)[source]¶
Add tags to the input.
- Parameters:
tags – A single tag or list/tuple/set of tags
- remove_tags(tags)[source]¶
Remove tags from the input.
- Parameters:
tags – A single tag or list/tuple/set of tags
- run_in_shell(workdir=None, manager=None, verbose=0) AbinitTask [source]¶
Helper method to quickly execute abinit in a shell.
- Parameters:
workdir – Working directory of the fake task used to compute the ibz. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.verbose – verbosity level.
- Returns:
an
abipy.flowtk.tasks.AbinitTask
of the executed calculation.
- class abipy.abio.inputs.MultiDataset(structure, pseudos, pseudo_dir='', ndtset=1)[source]¶
Bases:
object
This object is essentially a list of
abipy.abio.inputs.AbinitInput
objects. It provides an easy-to-use interface to apply global changes to the inputs.Let’s assume for example that multi contains two
AbinitInput
objects and that we want to set ecut to 1 in both dictionaries. A direct (lengthy) approach would be:- for inp in multi:
inp.set_vars(ecut=1)
or alternatively:
- for i in range(multi.ndtset):
multi[i].set_vars(ecut=1)
MultiDataset provides its own implementaion of __getattr__ so that one can simply use:
multi.set_vars(ecut=1)
multi.get(“ecut”) returns a list of values. It is equivalent to:
[inp[“ecut”] for inp in multi]
Note that if “ecut” is not present in one of the input of multi, the corresponding entry is set to None. A default value can be specified with:
multi.get(“paral_kgb”, 0)
Warning
MultiDataset does not support calculations done with different sets of pseudopotentials. The inputs can have different crystalline structures (as long as the atom types are equal) but each input in MultiDataset must have the same set of pseudopotentials.
- Error¶
alias of
AbinitInputError
- classmethod from_inputs(inputs: list[AbinitInput] | MultiDataset) MultiDataset [source]¶
Build a
abipy.abio.inputs.MultiDataset
from a list ofabipy.abio.inputs.AbinitInput
objects.
- classmethod replicate_input(input: AbinitInput, ndtset: int) MultiDataset [source]¶
Construct a multidataset with ndtset from the
abipy.abio.inputs.AbinitInput
input.
- property is_multidataset: bool¶
Used to understand if we have an AbinitInput or a MultiDataset in polymorphic APIs.
- property is_input: bool¶
Used to understand if we have an AbinitInput or a MultiDataset in polymorphic APIs.
- property pseudos¶
Pseudopotential objects.
- append(abinit_input: AbinitInput) None [source]¶
Add a
abipy.abio.inputs.AbinitInput
to the list.
- extend(abinit_inputs: AbinitInput | MultiDataset) None [source]¶
Extends self with a list of
abipy.abio.inputs.AbinitInput
objects.
- addnew_from(dtindex: int) None [source]¶
Add a new entry in the multidataset by copying the input with index
dtindex
.
- split_datasets() list[AbinitInput] [source]¶
Return list of
abipy.abio.inputs.AbinitInput
objects..
- deepcopy() MultiDataset [source]¶
Deep copy of the MultiDataset.
- to_string(mode='text', verbose=0, with_pseudos=True, files_file=False) str [source]¶
String representation i.e. the ABINIT input file.
- Parameters:
mode – Either
text
orhtml
if HTML output with links is wanted.with_pseudos – False if JSON section with pseudo data should not be added.
files_file – if True a string compatible with the presence of a files file will be generated. Otherwise all the required variables will be added to the string.
- get_vars_dataframe(*varnames) DataFrame [source]¶
Return pandas DataFrame with the value of the variables specified in varnames.
- filter_by_tags(tags=None, exclude_tags=None) MultiDataset [source]¶
Filters the inputs according to the tags
- Parameters:
tags – A single tag or list/tuple/set of tags
exclude_tags – A single tag or list/tuple/set of tags that should be excluded
- Returns:
A
abipy.abio.inputs.MultiDataset
containing the inputs containing all the requested tags.
- add_tags(tags, dtindeces=None) None [source]¶
Add tags to the selected inputs
- Parameters:
tags – A single tag or list/tuple/set of tags
dtindeces – a list of indices to which the tags will be added. None=all the inputs.
- remove_tags(tags, dtindeces=None) None [source]¶
Remove tags from the selected inputs
- Parameters:
tags – A single tag or list/tuple/set of tags
dtindeces – a list of indices from which the tags will be removed. None=all the inputs.
- filter_by_runlevel(runlevel) MultiDataset [source]¶
Return new
abipy.abio.inputs.MultiDataset
object in which only the inputs with the given runlevel are selected.
- exception abipy.abio.inputs.AnaddbInputError[source]¶
Bases:
Exception
Base error class for exceptions raised by AnaddbInput
- class abipy.abio.inputs.AnaddbInput(structure: Structure, comment: str = '', anaddb_args=None, anaddb_kwargs=None, spell_check: bool = True)[source]¶
Bases:
AbiAbstractInput
,MSONable
,Has_Structure
This object stores the anaddb variables.
Inheritance Diagram
- Error¶
alias of
AnaddbInputError
- classmethod from_dict(d: dict) AnaddbInput [source]¶
JSON interface used in pymatgen for easier serialization.
- property vars: dict¶
Dictionary with the input variables. Used to implement the dict-like interface.
- classmethod modes_at_qpoint(structure, qpoint, asr=2, chneut=1, dipdip=1, dipquad=1, quadquad=1, ifcflag=0, lo_to_splitting=False, directions=None, anaddb_args=None, anaddb_kwargs=None, spell_check=False) AnaddbInput [source]¶
Build an
abipy.abio.inputs.AnaddbInput
for the calculation of the phonon frequencies at a given q-point.- Parameters:
structure –
abipy.core.structure.Structure
objectqpoint – Reduced coordinates of the q-point where phonon frequencies and modes are wanted
asr – Anaddb input variable. See official documentation.
chneut – Anaddb input variable. See official documentation.
dipdp – Anaddb input variable. See official documentation.
ifcflag – Anaddb input variable. See official documentation.
dipquad – 1 to include DQ, QQ terms (provided DDB contains dynamical quadrupoles).
quadquad – 1 to include DQ, QQ terms (provided DDB contains dynamical quadrupoles).
lo_to_splitting – if True calculation of the LO-TO splitting will be included if qpoint==Gamma
directions – list of 3D directions along which the LO-TO splitting will be calculated. If None the three cartesian direction will be used
anaddb_args – List of tuples (key, value) with Anaddb input variables (default: empty)
anaddb_kwargs – Dictionary with Anaddb input variables (default: empty)
spell_check – False to disable spell checking for input variables.
- classmethod modes_at_qpoints(structure, qpoints, asr=2, chneut=1, dipdip=1, dipquad=1, quadquad=1, ifcflag=0, lo_to_splitting=False, directions=None, ngqpt=None, q1shft=(0, 0, 0), anaddb_args=None, anaddb_kwargs=None, spell_check=False) AnaddbInput [source]¶
Build an
abipy.abio.inputs.AnaddbInput
for the calculation of the phonon frequencies at a given q-point.- Parameters:
structure –
abipy.core.structure.Structure
objectqpoints – List of reduced coordinates of the q-point where phonon frequencies and modes are wanted
asr – Anaddb input variable. See official documentation.
chneut – Anaddb input variable. See official documentation.
dipdp – Anaddb input variable. See official documentation.
ifcflag – Anaddb input variable. See official documentation.
dipquad – 1 to include DQ, QQ terms (provided DDB contains dynamical quadrupoles).
quadquad – 1 to include DQ, QQ terms (provided DDB contains dynamical quadrupoles).
lo_to_splitting – if True calculation of the LO-TO splitting will be included if qpoint==Gamma
directions – list of 3D directions along which the LO-TO splitting will be calculated. If None the three cartesian direction will be used
ngqpt – Monkhorst-Pack divisions for the phonon Q-mesh. Required if ifcflag=1.
q1shft – Shifts used for the coarse Q-mesh
anaddb_args – List of tuples (key, value) with Anaddb input variables (default: empty)
anaddb_kwargs – Dictionary with Anaddb input variables (default: empty)
spell_check – False to disable spell checking for input variables.
- classmethod piezo_elastic(structure, relaxed_ion=True, stress_correction=False, asr=2, chneut=1, dipdip=1, anaddb_args=None, anaddb_kwargs=None) AnaddbInput [source]¶
Build an
abipy.abio.inputs.AnaddbInput
for the calculation of piezoelectric and elastic tensor calculations.- Parameters:
asr – Anaddb input variable. See official documentation.
chneut – Anaddb input variable. See official documentation.
dipdp – Anaddb input variable. See official documentation.
- classmethod phbands_and_dos(structure, ngqpt, nqsmall, qppa=None, ndivsm=20, line_density=None, q1shft=(0, 0, 0), qptbounds=None, asr=2, chneut=0, dipdip=1, dipquad=1, quadquad=1, dos_method='tetra', lo_to_splitting=False, with_ifc=False, anaddb_args=None, anaddb_kwargs=None, spell_check=False, comment=None) AnaddbInput [source]¶
Build an
abipy.abio.inputs.AnaddbInput
for the computation of phonon bands and phonon DOS.- Parameters:
structure –
abipy.core.structure.Structure
objectngqpt – Monkhorst-Pack divisions for the phonon Q-mesh (coarse one)
nqsmall – Used to generate the (dense) mesh for the DOS. It defines the number of q-points used to sample the smallest lattice vector.
qppa – Defines the homogeneous q-mesh used for the DOS in units of q-points per atom. Overrides nqsmall.
line_density – Defines the a density of k-points per reciprocal atom to plot the phonon dispersion. Overrides ndivsm.
ndivsm – Used to generate a normalized path for the phonon bands. If gives the number of divisions for the smallest segment of the path.
q1shft – Shifts used for the coarse Q-mesh
None (qptbounds Boundaries of the path. If) – depending on the input structure.
database (the path is generated from an internal) – depending on the input structure.
asr – Anaddb input variable. See official documentation.
chneut – Anaddb input variable. See official documentation.
dipdip – 1 to activate the treatment of the dipole-dipole interaction (requires BECS and dielectric tensor).
dipquad – 1 to include DQ, QQ terms (provided DDB contains dynamical quadrupoles).
quadquad – 1 to include DQ, QQ terms (provided DDB contains dynamical quadrupoles).
dos_method – Possible choices: “tetra”, “gaussian” or “gaussian:0.001 eV”. In the later case, the value 0.001 eV is used as gaussian broadening
lo_to_splitting – if True calculation of the LO-TO splitting will be included
with_ifc – if True calculation of the interatomic force constants will be included.
anaddb_args – List of tuples (key, value) with Anaddb input variables (default: empty)
anaddb_kwargs – Dictionary with Anaddb input variables (default: empty)
spell_check – False to disable spell checking for input variables.
comment – Optional string with a comment that will be placed at the beginning of the file.
- classmethod modes(structure, enunit=2, asr=2, chneut=1, anaddb_args=None, anaddb_kwargs=None) AnaddbInput [source]¶
Build an
abipy.abio.inputs.AnaddbInput
for the computation of phonon modes.- Parameters:
Structure –
abipy.core.structure.Structure
objectngqpt – Monkhorst-Pack divisions for the phonon Q-mesh (coarse one)
nqsmall – Used to generate the (dense) mesh for the DOS. It defines the number of q-points used to sample the smallest lattice vector.
q1shft – Shifts used for the coarse Q-mesh
None (qptbounds Boundaries of the path. If) – depending on the input structure.
database (the path is generated from an internal) – depending on the input structure.
asr – Anaddb input variable. See official documentation.
chneut – Anaddb input variable. See official documentation.
dipdp – Anaddb input variable. See official documentation.
anaddb_args – List of tuples (key, value) with Anaddb input variables (default: empty)
anaddb_kwargs – Dictionary with Anaddb input variables (default: empty)
- classmethod ifc(structure, ngqpt, ifcout=None, q1shft=(0, 0, 0), asr=2, chneut=1, dipdip=1, anaddb_args=None, anaddb_kwargs=None) AnaddbInput [source]¶
Build an
abipy.abio.inputs.AnaddbInput
for the computation of interatomic force constants.- Parameters:
structure –
abipy.core.structure.Structure
objectngqpt – Monkhorst-Pack divisions for the phonon Q-mesh (coarse one)
ifcout – Number of neighbouring atoms for which the ifc’s will be output. If None all the atoms in the big box.
q1shft – Shifts used for the coarse Q-mesh
asr – Anaddb input variable. See official documentation.
chneut – Anaddb input variable. See official documentation.
dipdip – Anaddb input variable. See official documentation.
anaddb_args – List of tuples (key, value) with Anaddb input variables (default: empty)
anaddb_kwargs – Dictionary with Anaddb input variables (default: empty)
- classmethod dfpt(structure, ngqpt=None, relaxed_ion=False, piezo=False, dde=False, strain=False, dte=False, raman=False, stress_correction=False, nqsmall=None, qppa=None, ndivsm=20, line_density=None, q1shft=(0, 0, 0), qptbounds=None, asr=2, chneut=1, dipdip=1, ramansr=1, alphon=1, dos_method='tetra', directions=None, anaddb_args=None, anaddb_kwargs=None, comment=None) AnaddbInput [source]¶
Builds an
abipy.abio.inputs.AnaddbInput
to post-process a generic DFPT calculation.- Parameters:
structure –
abipy.core.structure.Structure
object.ngqpt – Monkhorst-Pack divisions for the phonon Q-mesh (coarse one)
stress_correction – True to activate computation of stress correction in elastic tensor. Requires DDB with stress entries.
relaxed_ion – True to activate computation of relaxed-ion elastic and piezoelectric tensors. (assume the DDB has atomic perturbations at Gamma)
piezo – if True the piezoelectric tensor are calculated (requires piezoelectric perturbations)
dde – if True dielectric tensors will be calculated. If phonon band structure is calculated will also enable the calculation of the lo_to splitting (requires the DDE perturbations)
strain – if True the elastic tensors will be calculated (requires the strain perturbations)
dte – if True properties related to the nonlinear tensors will be calculated (requires third orders perturbations)
raman – if True the Raman tensor will be calculated (sets dte to True).
nqsmall – Used to generate the (dense) mesh for the DOS. It defines the number of q-points used to sample the smallest lattice vector.
qppa – Defines the homogeneous q-mesh used for the DOS in units of q-points per reciproval atom. Overrides nqsmall.
line_density – Defines the a density of k-points per reciprocal atom to plot the phonon dispersion. Overrides ndivsm.
ndivsm – Used to generate a normalized path for the phonon bands. If gives the number of divisions for the smallest segment of the path.
q1shft – Shifts used for the coarse Q-mesh
qptbounds – Boundaries of the path. If None, the path is generated from an internal database depending on the input structure.
asr – Anaddb input variable. See official documentation.
chneut – Anaddb input variable. See official documentation.
dipdp – Anaddb input variable. See official documentation.
ramansr – Anaddb input variable. See official documentation.
alphon – Anaddb input variable. See official documentation.
dos_method – Possible choices: “tetra”, “gaussian” or “gaussian:0.001 eV”. In the later case, the value 0.001 eV is used as gaussian broadening
directions – list of 3D directions along which the non analytical contribution will be calculated. If None the three cartesian direction will be used. Used only when dte=True.
anaddb_args – List of tuples (key, value) with Anaddb input variables (default: empty)
anaddb_kwargs – Dictionary with Anaddb input variables (default: empty)
comment – Optional string with a comment that will be placed at the beginning of the file.
- property structure: Structure¶
abipy.core.structure.Structure
object.
- to_string(sortmode=None, mode='text', verbose=0, files_file=False) str [source]¶
String representation.
- Parameters:
sortmode – “a” for alphabetical order, None if no sorting is wanted
mode – Either text or html if HTML output with links is wanted.
- set_qpath(ndivsm, qptbounds=None) dict [source]¶
Set the variables for the computation of the phonon band structure.
- Parameters:
ndivsm – Number of divisions for the smallest segment.
qptbounds – q-points defining the path in k-space. If None, we use the default high-symmetry k-path defined in the pymatgen database.
- set_autoqmesh(nqsmall) dict [source]¶
Set the variable nqpt for the sampling of the BZ.
- Parameters:
nqsmall – Number of divisions used to sample the smallest lattice vector.
- abivalidate(workdir: str = None, manager=None)[source]¶
Run ANADDB in dry-run mode to validate the input file.
- Parameters:
workdir – Working directory of the fake task used to compute the ibz. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- Returns:
retcode: Return code. 0 if OK. output_file: output file of the run. log_file: log file of the Abinit run, use log_file.read() to access its content. stderr_file: stderr file of the Abinit run. use stderr_file.read() to access its content. task: Task object
- Return type:
namedtuple with the following attributes
- class abipy.abio.inputs.OpticInput(**kwargs)[source]¶
Bases:
AbiAbstractInput
,MSONable
Input file for optic executable
- Error¶
alias of
OpticError
- property vars: dict¶
Dictionary with the input variables. Used to implement the dict-like interface.
- classmethod from_dict(d: dict) OpticInput [source]¶
JSON interface used in pymatgen for easier serialization.
- only_independent_chi_components(structure, assume_symmetric_tensor=False, symprec=0.001, angle_tolerance=5)[source]¶
Use the crystal system returned by spglib to find the independent components of the linear susceptibility tensor and set the appropriate variables.
- Parameters:
structure – Crystalline structure
assume_symmetric_tensor – True if tensor can be assumed symmetric. Note that the tensor is symmetric only for a lossless and non-optically active material.
symprec – Parameters passed to spglib.
angle_tolerance – Parameters passed to spglib.
- Returns:
Set internal variables and return list of components to compute.
- abivalidate(workdir=None, manager=None)[source]¶
Run OPTIC in dry-run mode to validate the input file. Note: This method is a stub, it always return retcode 0
- Parameters:
workdir – Working directory of the fake task used to compute the ibz. Use None for temporary dir.
manager –
abipy.flowtk.tasks.TaskManager
of the task. If None, the manager is initialized from the config file.
- Returns:
retcode: Return code. 0 if OK. output_file: output file of the run. log_file: log file of the Abinit run, use log_file.read() to access its content. stderr_file: stderr file of the Abinit run. use stderr_file.read() to access its content. task: Task object
- Return type:
namedtuple with the following attributes
- class abipy.abio.inputs.Cut3DInput(infile_path=None, output_filepath=None, options=None)[source]¶
Bases:
MSONable
This object stores the options to run a single cut3d analysis.
Warning
Converters with nspden > 1 won’t work since cut3d asks for the ispden index.
- classmethod den_to_3d_formatted(density_filepath, output_filepath) Cut3DInput [source]¶
Generates a cut3d input for the conversion to a 3D formatted format.
- Parameters:
density_filepath – absolute or relative path to the input density produced by abinit. Can be None to be defined at a later time.
output_filepath – path to the file that should be produced by cut3D, if required. At this stage it would be safer to use just the file name, as using an absolute or relative path may fail depending on the compiler.
- classmethod den_to_3d_indexed(density_filepath, output_filepath) Cut3DInput [source]¶
Generates a cut3d input for the conversion to a 3D indexed format.
- Parameters:
density_filepath – absolute or relative path to the input density produced by abinit. Can be None to be defined at a later time.
output_filepath – path to the file that should be produced by cut3D, if required. At this stage it would be safer to use just the file name, as using an absolute or relative path may fail depending on the compiler.
- classmethod den_to_molekel(density_filepath, output_filepath) Cut3DInput [source]¶
Generates a cut3d input for the conversion to a Molekel format.
- Parameters:
density_filepath – absolute or relative path to the input density produced by abinit. Can be None to be defined at a later time.
output_filepath – path to the file that should be produced by cut3D, if required. At this stage it would be safer to use just the file name, as using an absolute or relative path may fail depending on the compiler.
- classmethod den_to_tecplot(density_filepath, output_filepath) Cut3DInput [source]¶
Generates a cut3d input for the conversion to a Tecplot format.
- Parameters:
density_filepath – absolute or relative path to the input density produced by abinit. Can be None to be defined at a later time.
output_filepath – path to the file that should be produced by cut3D, if required. At this stage it would be safer to use just the file name, as using an absolute or relative path may fail depending on the compiler.
- classmethod den_to_xsf(density_filepath, output_filepath, shift=None) Cut3DInput [source]¶
Generates a cut3d input for the conversion to an xsf format.
- Parameters:
density_filepath – absolute or relative path to the input density produced by abinit. Can be None to be defined at a later time.
output_filepath – path to the file that should be produced by cut3D, if required. At this stage it would be safer to use just the file name, as using an absolute or relative path may fail depending on the compiler.
shift – a list of three integers defining the shift along the x, y, z axis. None if no shift is required.
- classmethod den_to_cube(density_filepath, output_filepath) Cut3DInput [source]¶
Generates a cut3d input for the conversion to a cube format.
- Parameters:
density_filepath – absolute or relative path to the input density produced by abinit. Can be None to be defined at a later time.
output_filepath – path to the file that should be produced by cut3D, if required. At this stage it would be safer to use just the file name, as using an absolute or relative path may fail depending on the compiler.
- classmethod hirshfeld(density_filepath, all_el_dens_paths) Cut3DInput [source]¶
Generates a cut3d input for the calculation of the Hirshfeld charges from the density.
- Parameters:
density_filepath – absolute or relative path to the input density produced by abinit. Can be None to be defined at a later time.
all_el_dens_paths – a list of paths to the all-electron density files corresponding to the elements defined in the abinit input. See https://www.abinit.org/downloads/all_core_electron for files.
- classmethod hirshfeld_from_fhi_path(density_filepath, structure, fhi_all_el_path) Cut3DInput [source]¶
Generates a cut3d input for the calculation of the Hirshfeld charges from the density. Automatically selects the all-electron density files from a folder containing the fhi all-electron density files: https://www.abinit.org/downloads/all_core_electron
This will work only if the input has been generated with AbinitInput and the Structure object is the same provided to AbinitInput because the order of the pseudos depends on znucl.
- Parameters:
density_filepath – absolute or relative path to the input density produced by abinit. Can be None to be defined at a later time.
structure – the structure used for the ground state calculation. Used to determine the elements
fhi_all_el_path – path to the folder containing the fhi all-electron density files
- classmethod from_dict(d: dict) Cut3DInput [source]¶
JSON interface used in pymatgen for easier serialization.
- abipy.abio.inputs.product_dict(d: dict)[source]¶
This function receives a dictionary d where each key defines a list of items or a simple scalar. It constructs the Cartesian product of the values (equivalent to nested for-loops), and returns a list of dictionaries with the values that would be used inside the loop.
>>> d = OrderedDict([("foo", [2, 4]), ("bar", 1)]) >>> product_dict(d) == [OrderedDict([('foo', 2), ('bar', 1)]), OrderedDict([('foo', 4), ('bar', 1)])] True >>> d = OrderedDict([("bar", [1,2]), ('foo', [3,4])]) >>> product_dict(d) == [{'bar': 1, 'foo': 3}, ... {'bar': 1, 'foo': 4}, ... {'bar': 2, 'foo': 3}, ... {'bar': 2, 'foo': 4}] True
- abipy.abio.inputs.kpoints_from_line_density(structure, line_density, symprec=0.01)[source]¶
Compute an high-symmetry k-path using pymatgen conventions and line_density
- Parameters:
line_density – Number of points in each segment is computed as: int(ceil(distance * line_density)) where distance is the lenght of the segment. This option is the recommended one if the k-path contains two consecutive high symmetry k-points that are very close as ndivsm > 0 may produce a very large number of wavevectors.
symprec – Symmetry precision passed to spglib.
Return: (nkpt, 3) numpy array with k-points in reduced coords.
outputs
Module¶
Objects used to extract and plot results from output files in text format.
- class abipy.abio.outputs.AbinitTextFile(filepath: str)[source]¶
Bases:
TextFile
Base class for the ABINIT main output files and log files.
- get_timer() AbinitTimerParser [source]¶
Timer data.
- class abipy.abio.outputs.AbinitLogFile(filepath: str)[source]¶
Bases:
AbinitTextFile
,NotebookWriter
Class representing the ABINIT log file.
Inheritance Diagram
- class abipy.abio.outputs.AbinitOutputFile(filepath: str)[source]¶
Bases:
AbinitTextFile
,NotebookWriter
Class representing the main Abinit output file.
Inheritance Diagram
- initial_structures()[source]¶
List of initial
abipy.core.structure.Structure
.
- final_structures()[source]¶
List of final
abipy.core.structure.Structure
.
- initial_structure()[source]¶
The
abipy.core.structure.Structure
defined in the output file.If the input file contains multiple datasets AND the datasets have different structures, this property returns None. In this case, one has to access the structure of the individual datasets. For example:
self.initial_structures[0]
gives the structure of the first dataset.
- final_structure()[source]¶
The
abipy.core.structure.Structure
defined in the output file.If the input file contains multiple datasets AND the datasets have different structures, this property returns None. In this case, one has to access the structure of the individual datasets. For example:
self.final_structures[0]
gives the structure of the first dataset.
- diff_datasets(dt_list1, dt_list2, with_params=True, differ='html', dryrun=False)[source]¶
Compare datasets
- get_dims_spginfo_dataframe(verbose: int = 0) DataFrame [source]¶
Parse the section with the dimensions of the calculation. Return Dataframe.
- get_dims_spginfo_dataset(verbose=0) tuple[dict, dict] [source]¶
Parse the section with the dimensions of the calculation. Return dictionaries
- Parameters:
verbose – Verbosity level.
Return: (dims_dataset, spginfo_dataset)
where dims_dataset[i] is an OrderedDict with the dimensions of dataset i spginfo_dataset[i] is a dictionary with space group information.
- next_gs_scf_cycle() GroundStateScfCycle [source]¶
Return the next
GroundStateScfCycle
in the file. None if not found.
- get_all_gs_scf_cycles() list[GroundStateScfCycle] [source]¶
Return list of
GroundStateScfCycle
objects. Empty list if no entry is found.
- next_d2de_scf_cycle() D2DEScfCycle [source]¶
Return
D2DEScfCycle
with information on the DFPT iterations. None if not found.
- get_all_d2de_scf_cycles() list[D2DEScfCycle] [source]¶
Return list of
D2DEScfCycle
objects. Empty list if no entry is found.
- plot(tight_layout=True, with_timer=False, show=True)[source]¶
Plot GS/DFPT SCF cycles and timer data found in the output file.
- Parameters:
with_timer – True if timer section should be plotted
- yield_figs(**kwargs)[source]¶
This function generates a predefined list of matplotlib figures with minimal input from the user.
- yield_plotly_figs(**kwargs)[source]¶
This function generates a predefined list of plotly figures with minimal input from the user.
- compare_gs_scf_cycles(others, show=True) list[Any] [source]¶
Produce and returns a list of matplotlib figure comparing the GS self-consistent cycle in self with the ones in others.
- Parameters:
others – list of
AbinitOutputFile
objects or strings with paths to output files.show – True to diplay plots.
- compare_d2de_scf_cycles(others, show=True) list[Any] [source]¶
Produce and returns a list of matplotlib figure comparing the DFPT self-consistent cycle in self with the ones in others.
- Parameters:
others – list of
AbinitOutputFile
objects or strings with paths to output files.show – True to diplay plots.
- abipy.abio.outputs.validate_output_parser(abitests_dir=None, output_files=None) int [source]¶
Validate/test Abinit output parser.
- Parameters:
dirpath – Abinit tests directory.
output_files – List of Abinit output files.
Return: Exit code.
- class abipy.abio.outputs.AboRobot(*args)[source]¶
Bases:
Robot
This robot analyzes the results contained in multiple Abinit output files. Can compare dimensions, SCF cycles, analyze timers.
Inheritance Diagram
- EXT = 'abo'¶
- get_dims_dataframe(with_time=True, index=None) DataFrame [source]¶
Build and return a
pandas.DataFrame
with the dimensions of the calculation.- Parameters:
with_time – True if walltime and cputime should be added
index – Index of the dataframe. Use relative paths of files if None.
- get_dataframe(with_geo=True, with_dims=True, abspath=False, funcs=None) DataFrame [source]¶
Return a
pandas.DataFrame
with the most important results and the filenames as index.- Parameters:
with_geo – True if structure info should be added to the dataframe
with_dims – True if dimensions should be added
abspath – True if paths in index should be absolute. Default: Relative to getcwd().
funcs – Function or list of functions to execute to add more data to the DataFrame. Each function receives a
abipy.electrons.gsr.GsrFile
object and returns a tuple (key, value) where key is a string with the name of column and value is the value to be inserted.
- get_time_dataframe() DataFrame [source]¶
Return a
pandas.DataFrame
with the wall-time, cpu time in seconds and the filenames as index.
- class abipy.abio.outputs.OutNcFile(filepath: str)[source]¶
Bases:
AbinitNcFile
Class representing the _OUT.nc file containing the dataset results produced at the end of the run. The netcdf variables can be accessed via instance attribute e.g.
outfile.ecut
. Provides integration with ipython.
robots
Module¶
This module defines the Robot BaseClass. Robots operates on multiple files and provide helper functions to plot the data e.g. convergence studies and to build pandas dataframes from the output files.
- class abipy.abio.robots.Robot(*args)[source]¶
Bases:
NotebookWriter
This is the base class from which all Robot subclasses should derive. A Robot supports the with context manager:
Usage example:
with Robot([("label1", "file1"), (label2, "file2")]) as robot: # Do something with robot. files are automatically closed when we exit. for label, abifile in self.items(): print(label)
- start = None¶
- HATCH = '/'¶
- classmethod get_supported_extensions() list[str] [source]¶
List of strings with extensions supported by Robot subclasses.
- classmethod class_for_ext(ext: str)[source]¶
Return the Robot subclass associated to the given extension.
- classmethod from_dir(top: str, walk: bool = True, abspath: bool = False) Robot [source]¶
Build a robot by scanning all files located within directory top. This method should be invoked with a concrete robot class, for example:
robot = GsrRobot.from_dir(“.”)
- Parameters:
top – Root directory
walk – if True, directories inside top are included as well.
abspath – True if paths in index should be absolute. Default: Relative to top.
- classmethod from_dirs(dirpaths: list[str], walk: bool = True, abspath: bool = False) Robot [source]¶
Similar to from_dir but accepts a list of directories instead of a single directory.
- Parameters:
walk – if True, directories inside top are included as well.
abspath – True if paths in index should be absolute. Default: Relative to top.
- classmethod from_dir_glob(pattern: str, walk: bool = True, abspath: bool = False) Robot [source]¶
This class method builds a robot by scanning all files located within the directories matching pattern as implemented by glob.glob This method should be invoked with a concrete robot class, for example:
robot = GsrRobot.from_dir_glob(“flow_dir/w*/outdata/”)
- Parameters:
pattern – Pattern string.
walk – if True, directories inside top are included as well.
abspath – True if paths in index should be absolute. Default: Relative to getcwd().
- classmethod class_handles_filename(filename: str) bool [source]¶
True if robot class handles filename.
- classmethod from_files(filenames, labels=None, abspath=False) Robot [source]¶
Build a Robot from a list of filenames. If labels is None, labels are automatically generated from absolute paths.
- Parameters:
abspath – True if paths in index should be absolute. Default: Relative to top.
- classmethod from_top_and_json_basename(top: str, json_basename: str) Robot [source]¶
Build a robot by scanning all json files located within directory top. and matching json_basename.
Example
gwr_robot = Robot.from_top_and_json_basename(“.”, “gwr_robot.json”)
- classmethod from_json_files(json_paths: list[str])[source]¶
Build Robot from a list of json files. Each json file should have a list of filepaths and @module and @class as required by msonable.
- classmethod from_work(work, outdirs='all', nids=None, ext=None, task_class=None) Robot [source]¶
Build a robot from a
abipy.flowtk.works.Work
object.
- classmethod from_flow(flow, outdirs='all', nids=None, ext=None, task_class=None) Robot [source]¶
Build a robot from a
abipy.flowtk.flows.Flow
object.- Parameters:
flow –
abipy.flowtk.flows.Flow
objectoutdirs – String used to select/ignore the files in the output directory of flow, works and tasks outdirs=”work” selects only the outdir of the Works, outdirs=”flow+task” selects the outdir of the Flow and the outdirs of the tasks outdirs=”-work” excludes the outdir of the Works. Cannot use
+
and-
flags in the same string. Default: all that is equivalent to “flow+work+task”nids – List of node identifiers used to select particular nodes. Not used if None
ext – File extension associated to the robot. Mainly used if method is invoked with the BaseClass
task_class – Task class or string with the class name used to select the tasks in the flow. None implies no filtering.
Usage example:
with abilab.GsrRobot.from_flow(flow) as robot: print(robot)
- Returns:
Robot
subclass.
- add_extfile_of_node(node, nids=None, task_class=None) None [source]¶
Add the file produced by this node to the robot.
- Parameters:
node –
abipy.flowtk.flows.Flow
orabipy.flowtk.works.Work
orabipy.flowtk.tasks.Task
object.nids – List of node identifiers used to select particular nodes. Not used if None
task_class – Task class or string with class name used to select the tasks in the flow. None implies no filtering.
- scan_dir(top: str, walk: bool = True) int [source]¶
Scan directory tree starting from
top
. Add files to the robot instance. Return: Number of files found.- Parameters:
top – Root directory
walk – if True, directories inside
top
are included as well.
- add_file(label, abifile, filter_abifile=None) None [source]¶
Add a file to the robot with the given label.
- Parameters:
label – String used to identify the file (must be unique, ax exceptions is raised if label is already present.
abifile – Specify the file to be added. Accepts strings (filepath) or abipy file-like objects.
filter_abifile – Function that receives an
abifile
object and returns True if the file should be added to the plotter.
- get_pyscript(filepath: str) RobotPythonScript [source]¶
Return RobotPythonScript to br used as context manager.
- static ordered_intersection(list_1, list_2) list [source]¶
Return ordered intersection of two lists. Items must be hashable.
- remove() None [source]¶
Close the file handle, remove the file from disk for each file in the robot.
- change_labels(new_labels: list[str], dryrun: bool = False) dict [source]¶
Change labels of the files.
- Parameters:
new_labels – List of strings (same length as self.abifiles)
dryrun – True to activate dryrun mode.
- Returns:
mapping new_label –> old_label.
- remap_labels(function: Callable, dryrun: bool = False) dict [source]¶
Change labels of the files by executing
function
- Parameters:
function – Callable object e.g. lambda function. The output of function(abifile) is used as new label. Note that the function shall not return duplicated labels when applied to self.abifiles.
dryrun – True to activate dryrun mode.
- Returns:
mapping new_label –> old_label.
- trim_paths(start=None) str [source]¶
Replace absolute filepaths in the robot with relative paths wrt to
start
directory. If start is None, os.getcwd() is used. Setself.start
attribute, returnself.start
.
- property labels: list[str]¶
List of strings used to create labels in matplotlib figures when plotting results taked from multiple files. By default, labels is initialized with the path of the files in the robot. Use change_labels to change the list.
- show_files(stream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>) None [source]¶
Show label –> file path
- getattrs_alleq(*aname_args) list [source]¶
Return list of attribute values for each attribute name in *aname_args.
- getattr_alleq(aname: str)[source]¶
Return the value of attribute aname. Try firs in self then in self.r Raises ValueError if value is not the same across all the files in the robot.
- has_different_structures(rtol=1e-05, atol=1e-08) str [source]¶
Check if structures are equivalent, return string with info about differences (if any).
- is_sortable(aname: str, raise_exc: bool = False) bool [source]¶
Return True if
aname
is an attribute of the netcdf file If raise_exc is True, AttributeError with an explicit message is raised.
- sortby(func_or_string: Callable | str | None, reverse: bool = False, unpack: bool = False) list[tuple] [source]¶
Sort files in the robot by
func_or_string
.- Parameters:
func_or_string – Can be None, string or callable defining the quantity to be used for sorting. If string, it’s assumed that the abifile has an attribute with the same name and getattr is invoked. If callable, the output of func_or_string(abifile) is used. If None, no sorting is performed.
reverse – If set to True, then the list elements are sorted as if each comparison were reversed.
unpack – Return (labels, abifiles, xs) if True
- Return: list of (label, abifile, xs) tuples where xs is obtained via
func_or_string
. or labels, abifiles, xs if
unpack
- group_and_sortby(hue: Callable | str, func_or_string: Callable | str | None) list[HueGroup] [source]¶
Group files by
hue
and, inside each group` sort items byfunc_or_string
.- Parameters:
hue – Variable that defines subsets of the data, which will be drawn on separate lines. Accepts callable or string If string, it’s assumed that the abifile has an attribute with the same name and getattr is invoked. Dot notation is also supported e.g. hue=”structure.formula” –> abifile.structure.formula If callable, the output of hue(abifile) is used.
func_or_string – Either None, string, callable defining the quantity to be used for sorting. If string, it’s assumed that the abifile has an attribute with the same name and getattr is invoked. If callable, the output of func_or_string(abifile) is used. If None, no sorting is performed.
Return: List of
HueGroup
instance.
- static sortby_label(sortby, param) str [source]¶
Return the label to be used when files are sorted with
sortby
.
- get_structure_dataframes(abspath=False, filter_abifile=None, **kwargs)[source]¶
Wrap dataframes_from_structures function.
- Parameters:
abspath – True if paths in index should be absolute. Default: Relative to getcwd().
filter_abifile – Function that receives an
abifile
object and returns True if the file should be added to the plotter.
- get_lattice_dataframe(**kwargs) DataFrame [source]¶
Return
pandas.DataFrame
with lattice parameters.
- get_coords_dataframe(**kwargs) DataFrame [source]¶
Return
pandas.DataFrame
with atomic positions.
- get_params_dataframe(abspath: bool = False) DataFrame [source]¶
Return
pandas.DataFrame
with the most important parameters. that are usually subject to convergence studies.- Parameters:
abspath – True if paths in index should be absolute. Default: Relative to top.
- static plot_xy_with_hue(data: DataFrame, x: str, y: str, hue: str, decimals=None, ax=None, xlims=None, ylims=None, fontsize=8, **kwargs) Any ¶
Plot y = f(x) relation for different values of hue. Useful for convergence tests done wrt two parameters.
- Parameters:
data –
pandas.DataFrame
containing columns x, y, and hue.x – Name of the column used as x-value
y – Name of the column(s) used as y-value
hue – Variable that define subsets of the data, which will be drawn on separate lines
decimals – Number of decimal places to round hue columns. Ignore if None
ax –
matplotlib.axes.Axes
or None if a new figure should be created.xlims – Set the data limits for the x(y)-axis. Accept tuple e.g. (left, right) or scalar e.g. left. If left (right) is None, default values are used
ylims – Set the data limits for the x(y)-axis. Accept tuple e.g. (left, right) or scalar e.g. left. If left (right) is None, default values are used
fontsize – Legend fontsize.
kwargs – Keyword arguments are passed to ax.plot method.
Returns:
matplotlib.figure.Figure
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- plot_convergence(item: str | Callable, sortby=None, hue=None, abs_conv=None, ax=None, fontsize=8, **kwargs) Any [source]¶
Plot the convergence of
item
wrt thesortby
parameter. Values can optionally be grouped byhue
.- Parameters:
item – Define the quantity to plot. Accepts callable or string If string, it’s assumed that the abifile has an attribute with the same name and getattr is invoked. Dot notation is also supported e.g. hue=”structure.formula” –> abifile.structure.formula. If callable, the output of item(abifile) is used.
sortby – Define the convergence parameter, sort files and produce plot labels. Can be None, string or function. If None, no sorting is performed. If string and not empty it’s assumed that the abifile has an attribute with the same name and getattr is invoked. If callable, the output of sortby(abifile) is used.
hue – Variable that define subsets of the data, which will be drawn on separate lines. Accepts callable or string. If string, it’s assumed that the abifile has an attribute with the same name and getattr is invoked. If callable, the output of hue(abifile) is used.
abs_conv – If not None, plot f(x) and abs(f(x) - f(x_inf)) in log scale.
ax –
matplotlib.axes.Axes
or None if a new figure should be created.fontsize – legend and label fontsize.
kwargs – keyword arguments passed to matplotlib plot method.
Returns:
matplotlib.figure.Figure
Example
robot.plot_convergence(“energy”) robot.plot_convergence(“energy”, sortby=”nkpt”) robot.plot_convergence(“pressure”, sortby=”nkpt”, hue=”tsmear”) robot.plot_convergence(“pressure”, sortby=”nkpt”, hue=”tsmear”, abs_conv=1e-3)
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- plot_convergence_items(items: list[str | Callable], sortby=None, hue=None, abs_conv=None, fontsize=8, **kwargs) Any [source]¶
Plot the convergence of a list of
items
wrt to thesortby
parameter. Values can optionally be grouped byhue
.- Parameters:
items – List of attributes (or callables) to be analyzed.
sortby – Define the convergence parameter, sort files and produce plot labels. Can be None, string or function. If None, no sorting is performed. If string and not empty it’s assumed that the abifile has an attribute with the same name and getattr is invoked. If callable, the output of sortby(abifile) is used.
hue – Variable that define subsets of the data, which will be drawn on separate lines. Accepts callable or string If string, it’s assumed that the abifile has an attribute with the same name and getattr is invoked. Dot notation is also supported e.g. hue=”structure.formula” –> abifile.structure.formula If callable, the output of hue(abifile) is used.
abs_conv – If not None, plot f(x) and abs(f(x) - f(x_inf)) in log scale. Since we are plotting multiple quantities, abs_conv is a dict mapping the name of the item to to the convergence.
fontsize – legend and label fontsize.
kwargs – keyword arguments are passed to ax.plot
Returns:
matplotlib.figure.Figure
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- get_convergence_analyzer(xname: str, ytols_dict: dict) ConvergenceAnalyzer [source]¶
The main difference is that ConvergenceAnalyze supports multiple convergence tolerances for a given y-value.
- Parameters:
xname – Name of the x-variable.
ytols_dict – dict mapping the name of the y-variable to the tolerance(s).
Example:
- plot_lattice_convergence(what_list=None, sortby=None, hue=None, fontsize=8, **kwargs) Any [source]¶
Plot the convergence of the lattice parameters (a, b, c, alpha, beta, gamma). wrt the``sortby`` parameter. Values can optionally be grouped by
hue
.- Parameters:
what_list – List of strings with the quantities to plot e.g. [“a”, “alpha”, “beta”]. None means all.
item – Define the quantity to plot. Accepts callable or string If string, it’s assumed that the abifile has an attribute with the same name and getattr is invoked. If callable, the output of item(abifile) is used.
sortby – Define the convergence parameter, sort files and produce plot labels. Can be None, string or function. If None, no sorting is performed. If string and not empty it’s assumed that the abifile has an attribute with the same name and getattr is invoked. If callable, the output of sortby(abifile) is used.
hue – Variable that define subsets of the data, which will be drawn on separate lines. Accepts callable or string If string, it’s assumed that the abifile has an attribute with the same name and getattr is invoked. Dot notation is also supported e.g. hue=”structure.formula” –> abifile.structure.formula If callable, the output of hue(abifile) is used.
ax –
matplotlib.axes.Axes
or None if a new figure should be created.fontsize – legend and label fontsize.
Returns:
matplotlib.figure.Figure
Example
robot.plot_lattice_convergence()
robot.plot_lattice_convergence(sortby=”nkpt”)
robot.plot_lattice_convergence(sortby=”nkpt”, hue=”tsmear”)
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- get_baserobot_code_cells(title=None) list [source]¶
Return list of jupyter cells with calls to methods provided by the base class.
- static get_yvals_item_abifiles(item: Any, abifiles: list) ndarray [source]¶
Extract values for a list of Abinit files.
- static plot_xvals_or_xstr_ax(ax, xs, yvals, fontsize, **kwargs) list [source]¶
Plot xs where xs can contain either numbers or strings.
- static plot_abs_conv(ax1, ax2, xs, yvals, abs_conv, xlabel, fontsize, hatch, **kwargs) None [source]¶
Plot |y - y_xmax| in log scale on ax2 and add hspan to ax1.
- class abipy.abio.robots.HueGroup(hvalue, xvalues, abifiles, labels)[source]¶
Bases:
object
This small object is used by
group_and_sortby
to store information about the group.
- class abipy.abio.robots.RobotPythonScript(robot: Robot, filepath_py: str)[source]¶
Bases:
object
Small object used to generate a python script that reconstructs the robot from a json file containing the list of files. Client code can then add additional logic to the script and write it to disk.
This object is typically used in the on_all_ok method of Works to generate ready-to-use python scripts to post-process/visualize the results.
Example
- with gsr_robot.get_pyscript(work.outdir.path_in(“gsr_robot.py”)) as script:
script.add_text(“a = 1”)
timer
Module¶
This module provides objects for extracting timing data from the ABINIT output files It also provides tools to analyze and to visualize the parallel efficiency.
- class abipy.abio.timer.AbinitTimerParser[source]¶
Bases:
AbinitTimerParser
,NotebookWriter
variable
Module¶
Support for Abinit input variables.