pymatgen.io.abinit package

This package implements basic input and output capabilities for Abinit.

Submodules

pymatgen.io.abinit.abiobjects module

Low-level objects providing an abstraction for the objects involved in the calculation.

class AbivarAble[source]

Bases: ABC

An AbivarAble object provides a method to_abivars that returns a dictionary with the abinit variables.

abstract to_abivars()[source]

Returns a dictionary with the abinit variables.

class Constraints[source]

Bases: AbivarAble

This object defines the constraints for structural relaxation.

to_abivars()[source]

Dictionary with Abinit variables.

class Electrons(spin_mode='polarized', smearing='fermi_dirac:0.1 eV', algorithm=None, nband=None, fband=None, charge=0.0, comment=None)[source]

Bases: AbivarAble, MSONable

The electronic degrees of freedom.

Constructor for Electrons object.

Parameters:
  • comment – String comment for Electrons

  • charge – Total charge of the system. Default is 0.

as_dict()[source]

JSON friendly dict representation.

classmethod from_dict(dct: dict) Self[source]

Build object from dictionary.

property nspden[source]

Number of independent density components.

property nspinor[source]

Number of independent spinor components.

property nsppol[source]

Number of independent spin polarizations.

to_abivars()[source]

Return dictionary with Abinit variables.

class ElectronsAlgorithm(*args, **kwargs)[source]

Bases: dict, AbivarAble, MSONable

Variables controlling the SCF/NSCF algorithm.

Initialize object.

as_dict()[source]

Convert object to dict.

classmethod from_dict(dct: dict) Self[source]

Build object from dict.

to_abivars()[source]

Dictionary with Abinit input variables.

class ExcHamiltonian(bs_loband, nband, mbpt_sciss, coulomb_mode, ecuteps, spin_mode='polarized', mdf_epsinf=None, exc_type='TDA', algo='haydock', with_lf=True, bs_freq_mesh=None, zcut=None, **kwargs)[source]

Bases: AbivarAble

This object contains the parameters for the solution of the Bethe-Salpeter equation.

Parameters:
  • bs_loband – Lowest band index (Fortran convention) used in the e-h basis set. Can be scalar or array of shape (nsppol,). Must be >= 1 and <= nband

  • nband – Max band index used in the e-h basis set.

  • mbpt_sciss – Scissors energy in Hartree.

  • coulomb_mode – Treatment of the Coulomb term.

  • ecuteps – Cutoff energy for W in Hartree.

  • mdf_epsinf – Macroscopic dielectric function \(\\epsilon_\\inf\) used in the model dielectric function.

  • exc_type – Approximation used for the BSE Hamiltonian

  • with_lf – True if local field effects are included <==> exchange term is included

  • bs_freq_mesh – Frequency mesh for the macroscopic dielectric function (start, stop, step) in Ha.

  • zcut – Broadening parameter in Ha.

  • **kwargs – Extra keywords.

property inclvkb[source]

Treatment of the dipole matrix element (NC pseudos, default is 2).

to_abivars()[source]

Returns a dictionary with the abinit variables.

property use_cg[source]

True if we are using the conjugate gradient method.

property use_direct_diago[source]

True if we are performing the direct diagonalization of the BSE Hamiltonian.

property use_haydock[source]

True if we are using the Haydock iterative technique.

class HilbertTransform(nomegasf, domegasf=None, spmeth=1, nfreqre=None, freqremax=None, nfreqim=None, freqremin=None)[source]

Bases: AbivarAble

Parameters for the Hilbert-transform method (Screening code) i.e. the parameters defining the frequency mesh used for the spectral function and the frequency mesh used for the polarizability.

Parameters:
  • nomegasf – Number of points for sampling the spectral function along the real axis.

  • domegasf – Step in Ha for the linear mesh used for the spectral function.

  • spmeth – Algorithm for the representation of the delta function.

  • nfreqre – Number of points along the real axis (linear mesh).

  • freqremax – Maximum frequency for W along the real axis (in hartree).

  • nfreqim – Number of point along the imaginary axis (Gauss-Legendre mesh).

  • freqremin – Minimum frequency for W along the real axis (in hartree).

to_abivars()[source]

Returns a dictionary with the abinit variables.

class KSampling(mode=KSamplingModes.monkhorst, num_kpts=0, kpts=((1, 1, 1),), kpt_shifts=(0.5, 0.5, 0.5), kpts_weights=None, use_symmetries=True, use_time_reversal=True, chksymbreak=None, comment=None)[source]

Bases: AbivarAble, MSONable

Input variables defining the K-point sampling.

Highly flexible constructor for KSampling objects. The flexibility comes at the cost of usability and in general, it is recommended that you use the default constructor only if you know exactly what you are doing and requires the flexibility. For most usage cases, the object be constructed far more easily using the convenience static constructors:

  1. gamma_only

  2. gamma_centered

  3. monkhorst

  4. monkhorst_automatic

  5. path

and it is recommended that you use those.

Parameters:
  • mode – Mode for generating k-poits. Use one of the KSamplingModes enum types.

  • num_kpts – Number of kpoints if mode is “automatic” Number of division for the sampling of the smallest segment if mode is “path”. Not used for the other modes

  • kpts – Number of divisions. Even when only a single specification is required, e.g. in the automatic scheme, the kpts should still be specified as a 2D array. e.g., [[20]] or [[2,2,2]].

  • kpt_shifts – Shifts for Kpoints.

  • use_symmetries – False if spatial symmetries should not be used to reduce the number of independent k-points.

  • use_time_reversal – False if time-reversal symmetry should not be used to reduce the number of independent k-points.

  • kpts_weights – Optional weights for kpoints. For explicit kpoints.

  • chksymbreak – Abinit input variable: check whether the BZ sampling preserves the symmetry of the crystal.

  • comment – String comment for Kpoints

Note

The default behavior of the constructor is monkhorst.

as_dict()[source]

Convert object to dict.

classmethod automatic_density(structure, kppa, chksymbreak=None, use_symmetries=True, use_time_reversal=True, shifts=(0.5, 0.5, 0.5))[source]

Returns an automatic Kpoint object based on a structure and a kpoint density. Uses Gamma centered meshes for hexagonal cells and Monkhorst-Pack grids otherwise.

Algorithm:

Uses a simple approach scaling the number of divisions along each reciprocal lattice vector proportional to its length.

Parameters:
  • structure – Input structure

  • kppa – Grid density

classmethod explicit_path(ndivsm, kpath_bounds)[source]

See _path for the meaning of the variables.

classmethod from_dict(dct: dict) Self[source]

Build object from dict.

classmethod gamma_centered(kpts=(1, 1, 1), use_symmetries=True, use_time_reversal=True)[source]

Convenient static constructor for an automatic Gamma centered Kpoint grid.

Parameters:
  • kpts – Subdivisions N_1, N_2 and N_3 along reciprocal lattice vectors.

  • use_symmetries – False if spatial symmetries should not be used to reduce the number of independent k-points.

  • use_time_reversal – False if time-reversal symmetry should not be used to reduce the number of independent k-points.

Returns:

KSampling object.

classmethod gamma_only()[source]

Gamma-only sampling.

property is_homogeneous: bool[source]

Homogeneous sampling.

classmethod monkhorst(ngkpt, shiftk=(0.5, 0.5, 0.5), chksymbreak=None, use_symmetries=True, use_time_reversal=True, comment=None)[source]

Convenient static constructor for a Monkhorst-Pack mesh.

Parameters:
  • ngkpt – Subdivisions N_1, N_2 and N_3 along reciprocal lattice vectors.

  • shiftk – Shift to be applied to the kpoints.

  • use_symmetries – Use spatial symmetries to reduce the number of k-points.

  • use_time_reversal – Use time-reversal symmetry to reduce the number of k-points.

Returns:

KSampling object.

classmethod monkhorst_automatic(structure, ngkpt, use_symmetries=True, use_time_reversal=True, chksymbreak=None, comment=None)[source]

Convenient static constructor for an automatic Monkhorst-Pack mesh.

Parameters:
  • structure – Structure object.

  • ngkpt – Subdivisions N_1, N_2 and N_3 along reciprocal lattice vectors.

  • use_symmetries – Use spatial symmetries to reduce the number of k-points.

  • use_time_reversal – Use time-reversal symmetry to reduce the number of k-points.

Returns:

KSampling object.

classmethod path_from_structure(ndivsm, structure)[source]

See _path for the meaning of the variables.

to_abivars()[source]

Dictionary with Abinit variables.

class KSamplingModes(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]

Bases: Enum

Enum if the different samplings of the BZ.

automatic = 3[source]
monkhorst = 1[source]
path = 2[source]
class ModelDielectricFunction(mdf_epsinf)[source]

Bases: AbivarAble

Model dielectric function used for BSE calculation.

Parameters:

mdf_epsinf – Value of epsilon_infinity.

to_abivars()[source]

Return dictionary with abinit variables.

class PPModel(mode='godby', plasmon_freq=None)[source]

Bases: AbivarAble, MSONable

Parameters defining the plasmon-pole technique. The common way to instantiate a PPModel object is via the class method PPModel.as_ppmodel(string).

Parameters:
  • mode – ppmodel type

  • plasmon_freq – Plasmon frequency in Ha.

as_dict()[source]

Convert object to dictionary.

classmethod as_ppmodel(obj)[source]

Constructs an instance of PPModel from obj.

Accepts obj in the form:
  • PPmodel instance

  • string. e.g “godby:12.3 eV”, “linden”.

classmethod from_dict(dct: dict) Self[source]

Build object from dictionary.

classmethod get_noppmodel()[source]

Calculation without plasmon-pole model.

to_abivars()[source]

Return dictionary with Abinit variables.

class PPModelModes(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]

Bases: Enum

Different kind of plasmon-pole models.

farid = 4[source]
godby = 1[source]
hybersten = 2[source]
linden = 3[source]
noppmodel = 0[source]
class RelaxationMethod(*args, **kwargs)[source]

Bases: AbivarAble, MSONable

This object stores the variables for the (constrained) structural optimization ionmov and optcell specify the type of relaxation. The other variables are optional and their use depend on ionmov and optcell. A None value indicates that we use abinit default. Default values can be modified by passing them to the constructor. The set of variables are constructed in to_abivars depending on ionmov and optcell.

Initialize object.

IONMOV_DEFAULT = 3[source]
OPTCELL_DEFAULT = 2[source]
as_dict()[source]

Convert object to dict.

classmethod atoms_and_cell(atoms_constraints=None)[source]

Relax atomic positions as well as unit cell.

classmethod atoms_only(atoms_constraints=None)[source]

Relax atomic positions, keep unit cell fixed.

classmethod from_dict(dct: dict) Self[source]

Build object from dictionary.

property move_atoms[source]

True if atoms must be moved.

property move_cell[source]

True if lattice parameters must be optimized.

to_abivars()[source]

Returns a dictionary with the abinit variables.

class Screening(ecuteps, nband, w_type='RPA', sc_mode='one_shot', hilbert=None, ecutwfn=None, inclvkb=2)[source]

Bases: AbivarAble

This object defines the parameters used for the computation of the screening function.

Parameters:
  • ecuteps – Cutoff energy for the screening (Ha units).

  • function (nband Number of bands for the Green's) –

  • w_type – Screening type

  • sc_mode – Self-consistency mode.

  • hilbert – Instance of HilbertTransform defining the parameters for the Hilber transform method.

  • ecutwfn – Cutoff energy for the wavefunctions (Default: ecutwfn == ecut).

  • inclvkb – Option for the treatment of the dipole matrix elements (NC pseudos).

to_abivars()[source]

Returns a dictionary with the abinit variables.

property use_hilbert[source]

True if we are using the Hilbert transform method.

class SelfEnergy(se_type, sc_mode, nband, ecutsigx, screening, gw_qprange=1, ppmodel=None, ecuteps=None, ecutwfn=None, gwpara=2)[source]

Bases: AbivarAble

This object defines the parameters used for the computation of the self-energy.

Parameters:
  • se_type – Type of self-energy (str)

  • sc_mode – Self-consistency mode.

  • nband – Number of bands for the Green’s function

  • ecutsigx – Cutoff energy for the exchange part of the self-energy (Ha units).

  • screening – Screening instance.

  • 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 computie the QP energies for all the point in the IBZ and one band above and one band below the Fermi level.

  • ppmodel – PPModel instance with the parameters used for the plasmon-pole technique.

  • ecuteps – Cutoff energy for the screening (Ha units).

  • ecutwfn – Cutoff energy for the wavefunctions (Default: ecutwfn == ecut).

property gwcalctyp[source]

Returns the value of the gwcalctyp input variable.

property symsigma[source]

1 if symmetries can be used to reduce the number of q-points.

to_abivars()[source]

Returns a dictionary with the abinit variables.

property use_ppmodel[source]

True if we are using the plasmon-pole approximation.

class Smearing(occopt, tsmear)[source]

Bases: AbivarAble, MSONable

Variables defining the smearing technique. The preferred way to instantiate a Smearing object is via the class method Smearing.as_smearing(string).

Build object with occopt and tsmear.

as_dict()[source]

JSON-friendly dict representation of Smearing.

classmethod as_smearing(obj)[source]

Constructs an instance of Smearing from obj. Accepts obj in the form:

  • Smearing instance

  • “name:tsmear” e.g. “gaussian:0.004” (Hartree units)

  • “name:tsmear units” e.g. “gaussian:0.1 eV”

  • None –> no smearing

classmethod from_dict(dct: dict) Self[source]

Build object from dict.

property mode[source]

String with smearing technique.

static nosmearing()[source]

Build object for calculations without smearing.

to_abivars()[source]

Return dictionary with Abinit variables.

class SpinMode(mode, nsppol, nspinor, nspden)[source]

Bases: SpinMode, AbivarAble, MSONable

Different configurations of the electron density as implemented in abinit: One can use as_spinmode to construct the object via SpinMode.as_spinmode (string) where string can assume the values:

  • polarized

  • unpolarized

  • afm (anti-ferromagnetic)

  • spinor (non-collinear magnetism)

  • spinor_nomag (non-collinear, no magnetism)

Create new instance of SpinMode(mode, nsppol, nspinor, nspden)

as_dict()[source]

Convert object to dict.

classmethod as_spinmode(obj)[source]

Converts obj into a SpinMode instance.

classmethod from_dict(dct: dict) Self[source]

Build object from dict.

to_abivars()[source]

Dictionary with Abinit input variables.

contract(string)[source]

assert contract(“1 1 1 2 2 3”) == “3*1 2*2 1*3” assert contract(“1 1 3 2 3”) == “2*1 1*3 1*2 1*3”

lattice_from_abivars(cls=None, *args, **kwargs)[source]

Returns a Lattice object from a dictionary with the Abinit variables acell and either rprim in Bohr or angdeg If acell is not given, the Abinit default is used i.e. [1,1,1] Bohr.

Parameters:

cls – Lattice class to be instantiated. Defaults to pymatgen.core.Lattice.

Example

lattice_from_abivars(acell=3*[10], rprim=np.eye(3))

species_by_znucl(structure: Structure) list[Species][source]

Return list of unique specie found in structure ordered according to sites.

Example

Site0: 0.5 0 0 O Site1: 0 0 0 Si

produces [Specie_O, Specie_Si] and not set([Specie_O, Specie_Si]) as in types_of_specie

structure_from_abivars(cls=None, *args, **kwargs) Structure[source]

Build a Structure object from a dictionary with ABINIT variables.

Parameters:

cls – Structure class to be instantiated. Defaults to Structure.

Example

al_structure = structure_from_abivars(

acell=3*[7.5], rprim=[0.0, 0.5, 0.5, 0.5, 0.0, 0.5, 0.5, 0.5, 0.0], typat=1, xred=[0.0, 0.0, 0.0], ntypat=1, znucl=13,

)

xred can be replaced with xcart or xangst.

structure_to_abivars(structure: Structure, enforce_znucl: list | None = None, enforce_typat: list | None = None, **kwargs)[source]

Receives a structure and returns a dictionary with ABINIT variables.

Parameters:
  • enforce_znucl (list) – ntypat entries with the value of Z for each type of atom. Used to change the default ordering. Defaults to None.

  • enforce_typat (list) – natom entries with the type index. Fortran conventions: start to count from 1. Used to change the default ordering.

pymatgen.io.abinit.abitimer 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 AbinitTimer(sections, info, cpu_time, wall_time)[source]

Bases: object

Container class storing the timing results.

Parameters:
  • sections – List of sections

  • info – Dictionary with extra info.

  • cpu_time – Cpu-time in seconds.

  • wall_time – Wall-time in seconds.

cpuwall_histogram(ax: Axes = None, **kwargs)[source]

Plot histogram with cpu- and wall-time on axis ax.

Parameters:

ax – matplotlib Axes or None if a new figure should be created.

Returns:

matplotlib figure

Return type:

plt.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.

get_dataframe(sort_key='wall_time', **kwargs)[source]

Return a pandas DataFrame with entries sorted according to sort_key.

get_section(section_name)[source]

Return section associated to section_name.

get_values(keys)[source]

Return a list of values associated to a particular list of keys.

names_and_values(key, minval=None, minfract=None, sorted=True)[source]

Select the entries whose value[key] is >= minval or whose fraction[key] is >= minfract Return the names of the sections and the corresponding values.

property ncpus[source]

Total number of CPUs employed.

order_sections(key, reverse=True)[source]

Sort sections according to the value of key.

pie(key='wall_time', minfract=0.05, ax: Axes = None, **kwargs)[source]

Plot pie chart for this timer.

Parameters:
  • key – Keyword used to extract data from the timer.

  • minfract – Don’t show sections whose relative weight is less that minfract.

  • ax – matplotlib Axes or None if a new figure should be created.

Returns:

matplotlib figure

Return type:

plt.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.

scatter_hist(ax: Axes = None, **kwargs)[source]

Scatter plot + histogram.

Parameters:

ax – matplotlib Axes or None if a new figure should be created.

Returns:

matplotlib figure

Return type:

plt.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.

sum_sections(keys)[source]

Sum value of keys.

to_csv(fileobj=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>)[source]

Write data on file fileobj using CSV format.

to_table(sort_key='wall_time', stop=None)[source]

Return a table (list of lists) with timer data.

totable(sort_key='wall_time', stop=None)[source]

Return a table (list of lists) with timer data.

exception AbinitTimerParseError[source]

Bases: ParseError

Errors raised by AbinitTimerParser.

class AbinitTimerParser[source]

Bases: Iterable

Responsible for parsing a list of output files, extracting the timing results and analyzing the results. Assume the Abinit output files have been produced with timopt -1.

Example

parser = AbinitTimerParser() parser.parse(list_of_files)

To analyze all *.abo files within top, use:

parser, paths, okfiles = AbinitTimerParser.walk(top=”.”, ext=”.abo”)

Initialize object.

BEGIN_TAG = '-<BEGIN_TIMER'[source]
END_TAG = '-<END_TIMER>'[source]
Error[source]

alias of AbinitTimerParseError

property filenames[source]

List of files that have been parsed successfully.

get_sections(section_name)[source]

Return the list of sections stored in self.timers() given section_name A fake section is returned if the timer does not have section_name.

parse(filenames)[source]

Read and parse a filename or a list of filenames. Files that cannot be opened are ignored. A single filename may also be given.

Returns:

list of successfully read files.

pefficiency()[source]

Analyze the parallel efficiency.

Returns:

ParallelEfficiency object.

plot_all(show=True, **kwargs)[source]

Call all plot methods provided by the parser.

plot_efficiency(key='wall_time', what='good+bad', nmax=5, ax: Axes = None, **kwargs)[source]

Plot the parallel efficiency.

Parameters:
  • key – Parallel efficiency is computed using the wall_time.

  • what – Specifies what to plot: good for sections with good parallel efficiency. bad for sections with bad efficiency. Options can be concatenated with +.

  • nmax – Maximum number of entries in plot

  • ax – matplotlib Axes or None if a new figure should be created.

kwargs

Meaning

linewidth

matplotlib linewidth. Default: 2.0

markersize

matplotlib markersize. Default: 10

Returns:

matplotlib 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.

plot_pie(key='wall_time', minfract=0.05, **kwargs)[source]

Plot pie charts of the different timers.

Parameters:
  • key – Keyword used to extract data from timers.

  • minfract – Don’t show sections whose relative weight is less that minfract.

Returns:

matplotlib 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.

plot_stacked_hist(key='wall_time', nmax=5, ax: Axes = None, **kwargs)[source]

Plot stacked histogram of the different timers.

Parameters:
  • key – Keyword used to extract data from the timers. Only the first nmax sections with largest value are show.

  • nmax – Maximum number of sections to show. Other entries are grouped together in the others section.

  • ax – matplotlib Axes or None if a new figure should be created.

Returns:

matplotlib 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.

section_names(ordkey='wall_time')[source]

Return the names of sections ordered by ordkey. For the time being, the values are taken from the first timer.

summarize(**kwargs)[source]

Return pandas DataFrame with the most important results stored in the timers.

timers(filename=None, mpi_rank='0')[source]

Return the list of timers associated to the given filename and MPI rank mpi_rank.

classmethod walk(top='.', ext='.abo')[source]

Scan directory tree starting from top, look for files with extension ext and parse timing data.

Returns:

the new object paths: the list of files found okfiles: list of files that have been parsed successfully.

(okfiles == paths) if all files have been parsed.

Return type:

parser

class AbinitTimerSection(name, cpu_time, cpu_fract, wall_time, wall_fract, ncalls, gflops)[source]

Bases: object

Record with the timing results associated to a section of code.

Parameters:
  • name – Name of the sections.

  • cpu_time – CPU time in seconds.

  • cpu_fract – Percentage of CPU time.

  • wall_time – Wall-time in seconds.

  • wall_fract – Percentage of wall-time.

  • ncalls – Number of calls

  • gflops – Gigaflops.

FIELDS = ('name', 'wall_time', 'wall_fract', 'cpu_time', 'cpu_fract', 'ncalls', 'gflops')[source]
NUMERIC_FIELDS = ('wall_time', 'wall_fract', 'cpu_time', 'cpu_fract', 'ncalls', 'gflops')[source]
STR_FIELDS = ('name',)[source]
classmethod fake()[source]

Return a fake section. Mainly used to fill missing entries if needed.

to_csvline(with_header=False)[source]

Return a string with data in CSV format. Add header if with_header.

to_dict()[source]

Convert object to dictionary.

to_tuple()[source]

Convert object to tuple.

class ParallelEfficiency(filenames, ref_idx, *args, **kwargs)[source]

Bases: dict

Store results concerning the parallel efficiency of the job.

Parameters:
  • filenames – List of filenames

  • ref_idx – Index of the Reference time (calculation done with the smallest number of cpus).

bad_sections(key='wall_time', criterion='mean', nmax=5)[source]

Return first nmax sections with worst value of key key using criterion criterion.

good_sections(key='wall_time', criterion='mean', nmax=5)[source]

Return first nmax sections with best value of key key using criterion criterion.

totable(stop=None, reverse=True)[source]

Return table (list of lists) with timing results.

Parameters:
  • stop – Include results up to stop. None for all

  • reverse – Put items with highest wall_time in first positions if True.

alternate(*iterables)[source]

[a[0], b[0], … , a[1], b[1], …, a[n], b[n] …] >>> alternate([1,4], [2,5], [3,6]) [1, 2, 3, 4, 5, 6].

pymatgen.io.abinit.inputs module

This module defines a simplified interface for generating ABINIT input files. Note that not all the features of Abinit are supported by BasicAbinitInput. For a more comprehensive implementation, use the AbinitInput object provided by AbiPy.

class AbstractInput[source]

Bases: MutableMapping, ABC

Abstract class defining the methods that must be implemented by Input objects.

deepcopy()[source]

Deep copy of the input.

pop_vars(keys)[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: Sequence[str], strict: bool = True) dict[str, InputVariable][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.

set_vars(*args, **kwargs)[source]

Set the value of the variables. Return dict with the variables added to the input.

Example

input.set_vars(ecut=10, ionmov=3)

set_vars_ifnotin(*args, **kwargs)[source]

Set the value of the variables but 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)

abstract to_str()[source]

Returns a string with the input.

abstract property vars[source]

Dictionary with the input variables. Used to implement dict-like interface.

write(filepath='run.abi')[source]

Write the input file to file to filepath.

class BasicAbinitInput(structure, pseudos: str | list[str] | list[Pseudo] | PseudoTable, pseudo_dir=None, comment=None, abi_args=None, abi_kwargs=None)[source]

Bases: AbstractInput, MSONable

This object stores the ABINIT variables for a single dataset.

Parameters:
  • structure (file with) – Parameters defining the crystalline structure. Accepts |Structure| object

  • structure

  • pseudos – Pseudopotentials to be used for the calculation. Accepts: string or list of strings with the name of the pseudopotential files, list of |Pseudo| objects or |PseudoTable| object.

  • pseudo_dir – Name of the directory where the pseudopotential files are located.

  • ndtset – Number of datasets.

  • comment – Optional string with a comment that will be placed at the beginning of the file.

  • abi_args – list of tuples (key, value) with the initial set of variables. Default: Empty

  • abi_kwargs – Dictionary with the initial set of variables. Default: Empty.

Error[source]

alias of BasicAbinitInputError

add_abiobjects(*abi_objects)[source]

This function receive a list of AbiVarable objects and add the corresponding variables to the input.

as_dict()[source]

JSON interface used in pymatgen for easier serialization.

property comment[source]

Optional string with comment. None if comment is not set.

classmethod from_dict(dct: dict) Self[source]

JSON interface used in pymatgen for easier serialization.

property isnc[source]

True if norm-conserving calculation.

property ispaw[source]

True if PAW calculation.

new_with_vars(*args, **kwargs)[source]

Return a new input with the given variables.

Example

new = input.new_with_vars(ecut=20)

pop_irdvars()[source]

Remove all the ird* variables present in self. Return dictionary with the variables that have been removed.

pop_tolerances()[source]

Remove all the tolerance variables present in self. Return dictionary with the variables that have been removed.

property pseudos[source]

List of |Pseudo| objects.

set_comment(comment)[source]

Set a comment to be included at the top of the file.

set_gamma_sampling()[source]

Gamma-only sampling of the BZ.

set_kmesh(ngkpt, shiftk, kptopt=1)[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_kpath(ndivsm, kptbounds=None, iscf=-2)[source]

Set the variables for the computation of the electronic band structure.

Parameters:
  • ndivsm – Number of divisions for the smallest segment.

  • 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_spin_mode(spin_mode)[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_structure(structure: Structure)[source]

Set structure.

property structure[source]

The |Structure| object associated to this input.

to_str(post=None, with_structure=True, with_pseudos=True, exclude=None)[source]

String representation.

Parameters:
  • 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.

property vars[source]

Dictionary with variables.

exception BasicAbinitInputError[source]

Bases: Exception

Base error class for exceptions raised by BasicAbinitInput.

class BasicMultiDataset(structure: Structure, pseudos, pseudo_dir='', ndtset=1)[source]

Bases: object

This object is essentially a list of BasicAbinitInput objects. that provides an easy-to-use interface to apply global changes to the the inputs stored in the objects.

Let’s assume for example that multi contains two BasicAbinitInput objects and we want to set ecut to 1 in both dictionaries. The direct 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)

BasicMultiDataset provides its own implementation of __getattr__ so that one can simply use:

multi.set_vars(ecut=1)

multi.get(“ecut”) returns a list of values. It’s 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

BasicMultiDataset 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 BasicMultiDataset must have the same set of pseudopotentials.

Parameters:
  • structure – file with the structure, |Structure| object or dictionary with ABINIT geo variable Accepts also list of objects that can be converted to Structure object. In this case, however, ndtset must be equal to the length of the list.

  • pseudos – String or list of string with the name of the pseudopotential files.

  • pseudo_dir – Name of the directory where the pseudopotential files are located.

  • ndtset – Number of datasets.

Error[source]

alias of BasicAbinitInputError

addnew_from(dtindex)[source]

Add a new entry in the multidataset by copying the input with index dtindex.

append(abinit_input)[source]

Add a |BasicAbinitInput| to the list.

deepcopy()[source]

Deep copy of the BasicMultiDataset.

extend(abinit_inputs)[source]

Extends self with a list of |BasicAbinitInput| objects.

classmethod from_inputs(inputs)[source]

Build object from a list of BasicAbinitInput objects.

property has_same_structures[source]

True if all inputs in BasicMultiDataset are equal.

property isnc[source]

True if norm-conserving calculation.

property ispaw[source]

True if PAW calculation.

property ndtset[source]

Number of inputs in self.

property pseudos[source]

Pseudopotential objects.

classmethod replicate_input(input, ndtset)[source]

Construct a multidataset with ndtset from the BasicAbinitInput input.

split_datasets()[source]

Return list of |BasicAbinitInput| objects..

to_str(with_pseudos=True)[source]

String representation i.e. the input file read by Abinit.

Parameters:

with_pseudos – False if JSON section with pseudo data should not be added.

write(filepath='run.abi')[source]

Write ndset input files to disk. The name of the file is constructed from the dataset index e.g. run0.abi.

class ShiftMode(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]

Bases: Enum

Class defining the mode to be used for the shifts. G: Gamma centered M: Monkhorst-Pack ((0.5, 0.5, 0.5)) S: Symmetric. Respects the chksymbreak with multiple shifts O: OneSymmetric. Respects the chksymbreak with a single shift (as in ‘S’ if a single shift is given, gamma

centered otherwise.

GammaCentered = 'G'[source]
MonkhorstPack = 'M'[source]
OneSymmetric = 'O'[source]
Symmetric = 'S'[source]
classmethod from_object(obj)[source]

Returns an instance of ShiftMode based on the type of object passed. Converts strings to ShiftMode depending on the initial letter of the string. G for GammaCentered, M for MonkhorstPack, S for Symmetric, O for OneSymmetric. Case insensitive.

as_structure(obj)[source]

Convert obj into a Structure. Accepts:

  • Structure object.

  • Filename

  • Dictionaries (MSONable format or dictionaries with abinit variables).

calc_shiftk(structure, symprec: float = 0.01, angle_tolerance=5)[source]

Find the values of shiftk and nshiftk appropriated for the sampling of the Brillouin zone.

When the primitive vectors of the lattice do NOT form a FCC or a BCC lattice, the usual (shifted) Monkhorst-Pack grids are formed by using nshiftk=1 and shiftk 0.5 0.5 0.5 . This is often the preferred k point sampling. For a non-shifted Monkhorst-Pack grid, use nshiftk=1 and shiftk 0.0 0.0 0.0, but there is little reason to do that.

When the primitive vectors of the lattice form a FCC lattice, with rprim:

0.0 0.5 0.5 0.5 0.0 0.5 0.5 0.5 0.0

the (very efficient) usual Monkhorst-Pack sampling will be generated by using nshiftk= 4 and shiftk:

0.5 0.5 0.5 0.5 0.0 0.0 0.0 0.5 0.0 0.0 0.0 0.5

When the primitive vectors of the lattice form a BCC lattice, with rprim:

-0.5 0.5 0.5

0.5 -0.5 0.5 0.5 0.5 -0.5

the usual Monkhorst-Pack sampling will be generated by using nshiftk= 2 and shiftk:

0.25 0.25 0.25

-0.25 -0.25 -0.25

However, the simple sampling nshiftk=1 and shiftk 0.5 0.5 0.5 is excellent.

For hexagonal lattices with hexagonal axes, e.g. rprim:

1.0 0.0 0.0

-0.5 sqrt(3)/2 0.0

0.0 0.0 1.0

one can use nshiftk= 1 and shiftk 0.0 0.0 0.5 In rhombohedral axes, e.g. using angdeg 3*60., this corresponds to shiftk 0.5 0.5 0.5, to keep the shift along the symmetry axis.

Returns:

Suggested value of shiftk.

ebands_input(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)[source]

Returns a |BasicMultiDataset| object for band structure calculations.

Parameters:
  • structure|Structure| object.

  • pseudos – List of filenames or list of |Pseudo| objects or |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 of 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).

gs_input(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)[source]

Returns a |BasicAbinitInput| for ground-state calculation.

Parameters:
  • structure|Structure| object.

  • pseudos – List of filenames or list of |Pseudo| objects or |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.

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')[source]

Returns a |BasicMultiDataset| for a structural relaxation. The first dataset optimizes the atomic positions at fixed unit cell. The second datasets optimizes both ions and unit cell parameters.

Parameters:
  • structure|Structure| object.

  • pseudos – List of filenames or list of |Pseudo| objects or |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.

num_valence_electrons(structure, pseudos)[source]

Returns the number of valence electrons.

Parameters:

pseudos – List of |Pseudo| objects or list of filenames.

pymatgen.io.abinit.netcdf module

Wrapper for netCDF readers.

class AbinitHeader(*args, **kwargs)[source]

Bases: AttrDict

Stores the values reported in the Abinit header.

Parameters:
  • args – Passthrough arguments for standard dict.

  • kwargs – Passthrough keyword arguments for standard dict.

to_str(verbose=0, title=None, **kwargs)[source]

String representation. kwargs are passed to pprint.pformat.

Parameters:
  • verbose – Verbosity level

  • title – Title string.

to_string(verbose=0, title=None, **kwargs)[source]

String representation. kwargs are passed to pprint.pformat.

Parameters:
  • verbose – Verbosity level

  • title – Title string.

class EtsfReader(path)[source]

Bases: NetcdfReader

This object reads data from a file written according to the ETSF-IO specifications.

We assume that the netcdf file contains at least the crystallographic section.

Open the Netcdf file specified by path (read mode).

chemical_symbols()[source]

Chemical symbols char [number of atom species][symbol length].

read_abinit_hdr()[source]

Read the variables associated to the Abinit header.

Return AbinitHeader

read_abinit_xcfunc()[source]

Read ixc from an Abinit file. Return XcFunc object.

read_structure(cls=<class 'pymatgen.core.structure.Structure'>)[source]

Returns the crystalline structure stored in the rootgrp.

type_idx_from_symbol(symbol)[source]

Returns the type index from the chemical symbol. Note python convention.

class NO_DEFAULT[source]

Bases: object

Signal that read_value should raise an Error.

class NetcdfReader(path)[source]

Bases: object

Wraps and extends netCDF4.Dataset. Read only mode. Supports with statements.

Additional documentation available at:

http://netcdf4-python.googlecode.com/svn/trunk/docs/netCDF4-module.html

Open the Netcdf file specified by path (read mode).

Error[source]

alias of NetcdfReaderError

close()[source]

Close the file.

print_tree()[source]

Print all the groups in the file.

read_dimvalue(dimname, path='/', default=<class 'pymatgen.io.abinit.netcdf.NO_DEFAULT'>)[source]

Returns the value of a dimension.

Parameters:
  • dimname – Name of the variable

  • path – path to the group.

  • default – return default if dimname is not present and default is not NO_DEFAULT else raise self.Error.

read_keys(keys, dict_cls=<class 'monty.collections.AttrDict'>, path='/')[source]

Read a list of variables/dimensions from file. If a key is not present the corresponding entry in the output dictionary is set to None.

read_value(varname, path='/', cmode=None, default=<class 'pymatgen.io.abinit.netcdf.NO_DEFAULT'>)[source]

Returns the values of variable with name varname in the group specified by path.

Parameters:
  • varname – Name of the variable

  • path – path to the group.

  • cmode – if cmode==”c”, a complex ndarrays is constructed and returned (netcdf does not provide native support from complex datatype).

  • default – returns default if varname is not present. self.Error is raised if default is set to NO_DEFAULT

Returns:

numpy array if varname represents an array, scalar otherwise.

read_variable(varname, path='/')[source]

Returns the variable with name varname in the group specified by path.

read_varnames(path='/')[source]

List of variable names stored in the group specified by path.

walk_tree(top=None)[source]

Navigate all the groups in the file starting from top. If top is None, the root group is used.

exception NetcdfReaderError[source]

Bases: Exception

Base error class for NetcdfReader.

as_etsfreader(file)[source]

Return an EtsfReader. Accepts filename or EtsfReader.

as_ncreader(file)[source]

Convert file into a NetcdfReader instance. Returns reader, close_it where close_it is set to True if we have to close the file before leaving the procedure.

structure_from_ncdata(ncdata, site_properties=None, cls=<class 'pymatgen.core.structure.Structure'>)[source]

Reads and returns a pymatgen structure from a NetCDF file containing crystallographic data in the ETSF-IO format.

Parameters:
  • ncdata – filename or NetcdfReader instance.

  • site_properties – Dictionary with site properties.

  • cls – The Structure class to instantiate.

pymatgen.io.abinit.pseudos module

This module provides objects describing the basic parameters of the pseudopotentials used in Abinit, and a parser to instantiate pseudopotential objects..

class AbinitHeader[source]

Bases: dict

Dictionary whose keys can be also accessed as attributes.

class AbinitPseudo(path, header)[source]

Bases: Pseudo

An AbinitPseudo is a pseudopotential whose file contains an abinit header.

Parameters:
  • path – Filename.

  • header – AbinitHeader instance.

property Z[source]

The atomic number of the atom.

property Z_val[source]

Valence charge.

property l_local[source]

Angular momentum used for the local part.

property l_max[source]

Maximum angular momentum.

property summary[source]

Summary line reported in the ABINIT header.

property supports_soc[source]

True if the pseudo can be used in a calculation with spin-orbit coupling. Base classes should provide a concrete implementation that computes this value.

class Hint(ecut, pawecutdg=None)[source]

Bases: object

Suggested value for the cutoff energy [Hartree units] and the cutoff energy for the dense grid (only for PAW pseudos).

as_dict()[source]

Return dictionary for MSONable protocol.

classmethod from_dict(dct: dict) Self[source]

Build instance from dictionary (MSONable protocol).

class NcAbinitHeader(summary, **kwargs)[source]

Bases: AbinitHeader

The abinit header found in the NC pseudopotential files.

static fhi_header(filename, ppdesc)[source]

Parse the FHI abinit header. Example:

Troullier-Martins psp for element Sc Thu Oct 27 17:33:22 EDT 1994

21.00000 3.00000 940714 zatom, zion, pspdat 1 1 2 0 2001 .00000 pspcod,pspxc,lmax,lloc,mmax,r2well 1.80626423934776 .22824404341771 1.17378968127746 rchrg,fchrg,qchrg

static gth_header(filename, ppdesc)[source]

Parse the GTH abinit header. Example:

Goedecker-Teter-Hutter Wed May 8 14:27:44 EDT 1996 1 1 960508 zatom,zion,pspdat 2 1 0 0 2001 0. pspcod,pspxc,lmax,lloc,mmax,r2well 0.2000000 -4.0663326 0.6778322 0 0 rloc, c1, c2, c3, c4 0 0 0 rs, h1s, h2s 0 0 rp, h1p

1.36 .2 0.6 rcutoff, rloc

static hgh_header(filename, ppdesc)[source]

Parse the HGH abinit header. Example:

Hartwigsen-Goedecker-Hutter psp for Ne, from PRB58, 3641 (1998)

10 8 010605 zatom,zion,pspdat 3 1 1 0 2001 0 pspcod,pspxc,lmax,lloc,mmax,r2well

static oncvpsp_header(filename, ppdesc)[source]

Parse the ONCVPSP abinit header. Example:

Li ONCVPSP r_core= 2.01 3.02

3.0000 3.0000 140504 zatom,zion,pspd

8 2 1 4 600 0 pspcod,pspxc,lmax,lloc,mmax,r2well

5.99000000 0.00000000 0.00000000 rchrg fchrg qchrg

2 2 0 0 0 nproj 0 extension_switch

0 -2.5000025868368D+00 -1.2006906995331D+00

1 0.0000000000000D+00 0.0000000000000D+00 0.0000000000000D+00 2 1.0000000000000D-02 4.4140499497377D-02 1.9909081701712D-02

static tm_header(filename, ppdesc)[source]

Parse the TM abinit header. Example:

Troullier-Martins psp for element Fm Thu Oct 27 17:28:39 EDT 1994 100.00000 14.00000 940714 zatom, zion, pspdat

1 1 3 0 2001 .00000 pspcod,pspxc,lmax,lloc,mmax,r2well 0 4.085 6.246 0 2.8786493 l,e99.0,e99.9,nproj,rcpsp .00000000 .0000000000 .0000000000 .00000000 rms,ekb1,ekb2,epsatm 1 3.116 4.632 1 3.4291849 l,e99.0,e99.9,nproj,rcpsp .00000000 .0000000000 .0000000000 .00000000 rms,ekb1,ekb2,epsatm 2 4.557 6.308 1 2.1865358 l,e99.0,e99.9,nproj,rcpsp .00000000 .0000000000 .0000000000 .00000000 rms,ekb1,ekb2,epsatm 3 23.251 29.387 1 2.4776730 l,e99.0,e99.9,nproj,rcpsp .00000000 .0000000000 .0000000000 .00000000 rms,ekb1,ekb2,epsatm 3.62474762267880 .07409391739104 3.07937699839200 rchrg,fchrg,qchrg

class NcAbinitPseudo(path, header)[source]

Bases: NcPseudo, AbinitPseudo

Norm-conserving pseudopotential in the Abinit format.

Parameters:
  • path – Filename.

  • header – AbinitHeader instance.

property Z[source]

The atomic number of the atom.

property Z_val[source]

Number of valence electrons.

property l_local[source]

Angular momentum used for the local part.

property l_max[source]

Maximum angular momentum.

property nlcc_radius[source]

Radius at which the core charge vanish (i.e. cut-off in a.u.). Returns 0.0 if nlcc is not used.

property summary[source]

Summary line reported in the ABINIT header.

class NcPseudo[source]

Bases: ABC

Abstract class defining the methods that must be implemented by the concrete classes representing norm-conserving pseudopotentials.

property has_nlcc[source]

True if the pseudo is generated with non-linear core correction.

abstract property nlcc_radius[source]

Radius at which the core charge vanish (i.e. cut-off in a.u.). Returns 0.0 if nlcc is not used.

property rcore[source]

Radius of the pseudization sphere in a.u.

class PawAbinitHeader(summary, **kwargs)[source]

Bases: AbinitHeader

The abinit header found in the PAW pseudopotential files.

static paw_header(filename, ppdesc)[source]

Parse the PAW abinit header. Examples:

Paw atomic data for element Ni - Generated by AtomPAW (N. Holzwarth) + AtomPAW2Abinit v3.0.5

28.000 18.000 20061204 : zatom,zion,pspdat 7 7 2 0 350 0. : pspcod,pspxc,lmax,lloc,mmax,r2well

paw3 1305pspfmt,creatorID

5 13 : basis_size,lmn_size

0 0 1 1 2 : orbitals 3 : number_of_meshes 1 3 350 1.1803778368E-05 3.5000000000E-02 : mesh 1, type,size,rad_step[,log_step] 2 1 921 2.500000000000E-03 : mesh 2, type,size,rad_step[,log_step] 3 3 391 1.1803778368E-05 3.5000000000E-02 : mesh 3, type,size,rad_step[,log_step]

2.3000000000 : r_cut(SPH)

2 0.

Another format:

C (US d-loc) - PAW data extracted from US-psp (D.Vanderbilt) - generated by USpp2Abinit v2.3.0

6.000 4.000 20090106 : zatom,zion,pspdat

7 11 1 0 560 0. : pspcod,pspxc,lmax,lloc,mmax,r2well

paw4 2230pspfmt,creatorID

4 8 : basis_size,lmn_size

0 0 1 1 : orbitals 5 : number_of_meshes 1 2 560 1.5198032759E-04 1.6666666667E-02 : mesh 1, type,size,rad_step[,log_step] 2 2 556 1.5198032759E-04 1.6666666667E-02 : mesh 2, type,size,rad_step[,log_step] 3 2 576 1.5198032759E-04 1.6666666667E-02 : mesh 3, type,size,rad_step[,log_step] 4 2 666 1.5198032759E-04 1.6666666667E-02 : mesh 4, type,size,rad_step[,log_step] 5 2 673 1.5198032759E-04 1.6666666667E-02 : mesh 5, type,size,rad_step[,log_step]

1.5550009124 : r_cut(PAW)

3 0. : shape_type,rshape

Yet nnother one:

Paw atomic data for element Si - Generated by atompaw v3.0.1.3 & AtomPAW2Abinit v3.3.1

14.000 4.000 20120814 : zatom,zion,pspdat 7 11 1 0 663 0. : pspcod,pspxc,lmax,lloc,mmax,r2well

paw5 1331pspfmt,creatorID

4 8 : basis_size,lmn_size

0 0 1 1 : orbitals 5 : number_of_meshes 1 2 663 8.2129718540404674E-04 1.1498160595656655E-02 : mesh 1, type,size,rad_step[,log_step] 2 2 658 8.2129718540404674E-04 1.1498160595656655E-02 : mesh 2, type,size,rad_step[,log_step] 3 2 740 8.2129718540404674E-04 1.1498160595656655E-02 : mesh 3, type,size,rad_step[,log_step] 4 2 819 8.2129718540404674E-04 1.1498160595656655E-02 : mesh 4, type,size,rad_step[,log_step] 5 2 870 8.2129718540404674E-04 1.1498160595656655E-02 : mesh 5, type,size,rad_step[,log_step]

1.5669671236 : r_cut(PAW)

2 0. : shape_type,rshape

class PawAbinitPseudo(path, header)[source]

Bases: PawPseudo, AbinitPseudo

Paw pseudopotential in the Abinit format.

Parameters:
  • path – Filename.

  • header – AbinitHeader instance.

property paw_radius[source]

Radius of the PAW sphere in a.u.

property supports_soc[source]

True if the pseudo can be used in a calculation with spin-orbit coupling. Base classes should provide a concrete implementation that computes this value.

class PawPseudo[source]

Bases: ABC

Abstract class that defines the methods that must be implemented by the concrete classes representing PAW pseudopotentials.

abstract property paw_radius[source]

Radius of the PAW sphere in a.u.

property rcore[source]

Alias of paw_radius.

class PawXmlSetup(filepath)[source]

Bases: Pseudo, PawPseudo

Setup class for PawXml.

Parameters:

filepath (str) – Path to the XML file.

property Z[source]

The atomic number of the atom.

property Z_val[source]

Number of valence electrons.

ae_core_density()[source]

The all-electron radial density.

ae_partial_waves()[source]

Dictionary with the AE partial waves indexed by state.

property l_local[source]

Angular momentum used for the local part.

property l_max[source]

Maximum angular momentum.

property paw_radius[source]

Radius of the PAW sphere in a.u.

plot_densities(ax: plt.Axes = None, **kwargs)[source]

Plot the PAW densities.

Parameters:

ax – matplotlib Axes or None if a new figure should be created.

Returns:

matplotlib 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.

plot_projectors(ax: plt.Axes = None, fontsize=12, **kwargs)[source]

Plot the PAW projectors.

Parameters:

ax – matplotlib Axes or None if a new figure should be created.

Returns:

matplotlib figure

Return type:

plt.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.

plot_waves(ax: plt.Axes = None, fontsize=12, **kwargs)[source]

Plot the AE and the pseudo partial waves.

Parameters:
  • ax – matplotlib Axes or None if a new figure should be created.

  • fontsize – fontsize for legends and titles

Returns:

matplotlib figure

Return type:

plt.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.

projector_functions()[source]

Dictionary with the PAW projectors indexed by state.

pseudo_core_density()[source]

The pseudized radial density.

property pseudo_partial_waves[source]

Dictionary with the pseudo partial waves indexed by state.

root()[source]

Root tree of XML.

property summary[source]

String summarizing the most important properties.

property supports_soc[source]

Here I assume that the ab-initio code can treat the SOC within the on-site approximation.

yield_figs(**kwargs)[source]

This function generates a predefined list of matplotlib figures with minimal input from the user.

class Pseudo[source]

Bases: MSONable, ABC

Abstract base class defining the methods that must be implemented by the concrete pseudo-potential sub-classes.

abstract property Z: int[source]

The atomic number of the atom.

abstract property Z_val: int[source]

Valence charge.

as_dict(**kwargs)[source]

Return dictionary for MSONable protocol.

classmethod as_pseudo(obj)[source]

Convert obj into a pseudo. Accepts:

  • Pseudo object.

  • string defining a valid path.

as_tmpfile(tmpdir=None)[source]

Copy the pseudopotential to a temporary a file and returns a new pseudopotential object. Useful for unit tests in which we have to change the content of the file.

Parameters:

tmpdir – If None, a new temporary directory is created and files are copied here else tmpdir is used.

property basename: str[source]

File basename.

compute_md5()[source]

Compute and return MD5 hash value.

property djrepo_path[source]

The path of the djrepo file. None if file does not exist.

property element: Element[source]

Pymatgen Element.

property filepath: str[source]

Absolute path to pseudopotential file.

classmethod from_dict(dct: dict) Self[source]

Build instance from dictionary (MSONable protocol).

static from_file(filename) Pseudo[source]

Build an instance of a concrete Pseudo subclass from filename. Note: the parser knows the concrete class that should be instantiated Client code should rely on the abstract interface provided by Pseudo.

property has_dojo_report[source]

True if the pseudo has an associated DOJO_REPORT section.

property has_hints[source]

True if self provides hints on the cutoff energy.

hint_for_accuracy(accuracy='normal')[source]

Returns a Hint object with the suggested value of ecut [Ha] and pawecutdg [Ha] for the given accuracy. ecut and pawecutdg are set to zero if no hint is available.

Parameters:

accuracy – [“low”, “normal”, “high”]

property isnc: bool[source]

True if norm-conserving pseudopotential.

property ispaw: bool[source]

True if PAW pseudopotential.

abstract property l_local: int[source]

Angular momentum used for the local part.

abstract property l_max: int[source]

Maximum angular momentum.

md5()[source]

MD5 hash value.

open_pspsfile(ecut=20, pawecutdg=None)[source]

Calls Abinit to compute the internal tables for the application of the pseudopotential part. Returns PspsFile object providing methods to plot and analyze the data or None if file is not found or it’s not readable.

Parameters:
  • ecut – Cutoff energy in Hartree.

  • pawecutdg – Cutoff energy for the PAW double grid.

abstract property summary: str[source]

String summarizing the most important properties.

abstract property supports_soc[source]

True if the pseudo can be used in a calculation with spin-orbit coupling. Base classes should provide a concrete implementation that computes this value.

property symbol: str[source]

Element symbol.

to_str(verbose=0) str[source]

String representation.

property type: str[source]

Type of pseudo.

exception PseudoParseError[source]

Bases: ParseError

Base Error class for the exceptions raised by PseudoParser.

class PseudoParser[source]

Bases: object

Responsible for parsing pseudopotential files and returning pseudopotential objects.

Usage:

pseudo = PseudoParser().parse(“filename”)

Error[source]

alias of PseudoParseError

parse(filename)[source]

Read and parse a pseudopotential file. Main entry point for client code.

Returns:

pseudopotential object or None if filename is not a valid pseudopotential file.

read_ppdesc(filename)[source]

Read the pseudopotential descriptor from filename.

Returns:

Pseudopotential descriptor. None if filename is not a valid pseudopotential file.

Raises:

PseudoParseError

scan_directory(dirname, exclude_exts=(), exclude_fnames=())[source]

Analyze the files contained in directory dirname.

Parameters:
  • dirname – directory path

  • exclude_exts – list of file extensions that should be skipped.

  • exclude_fnames – list of file names that should be skipped.

Returns:

List of pseudopotential objects.

class PseudoTable(pseudos: Sequence[Pseudo])[source]

Bases: Sequence, MSONable

Define the pseudopotentials from the element table. Individidual elements are accessed by name, symbol or atomic number.

For example, the following all retrieve iron:

print elements[26] Fe print elements.Fe Fe print elements.symbol(‘Fe’) Fe print elements.name(‘iron’) Fe print elements.isotope(‘Fe’) Fe

Parameters:

pseudos – List of pseudopotentials or filepaths.

all_combinations_for_elements(element_symbols)[source]

Return a list with all the possible combination of pseudos for the given list of element_symbols. Each item is a list of pseudopotential objects.

Example

table.all_combinations_for_elements([“Li”, “F”])

property allnc: bool[source]

True if all pseudos are norm-conserving.

property allpaw[source]

True if all pseudos are PAW.

as_dict(**kwargs)[source]

Return dictionary for MSONable protocol.

classmethod as_table(items)[source]

Return an instance of PseudoTable from the iterable items.

classmethod from_dict(dct: dict) Self[source]

Build instance from dictionary (MSONable protocol).

classmethod from_dir(top, exts=None, exclude_dirs='_*')[source]

Find all pseudos in the directory tree starting from top.

Parameters:
  • top – Top of the directory tree

  • exts – List of files extensions. if exts == “all_files” we try to open all files in top

  • exclude_dirs – Wildcard used to exclude directories.

Returns:

PseudoTable sorted by atomic number Z.

get_pseudos_for_structure(structure: Structure)[source]

Return the list of Pseudo objects to be used for this Structure.

Parameters:

structure – pymatgen Structure.

Raises:
  • ValueError

  • multiple occurrences are present in the table.

is_complete(zmax=118) bool[source]

True if table is complete i.e. all elements with Z < zmax have at least on pseudopotential.

print_table(stream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>, filter_function=None)[source]

A pretty ASCII printer for the periodic table, based on some filter_function.

Parameters:
  • stream – file-like object

  • filter_function – A filtering function that take a Pseudo as input and returns a boolean. For example, setting filter_function = lambda p: p.Z_val > 2 will print a periodic table containing only pseudos with Z_val > 2.

pseudo_with_symbol(symbol, allow_multi=False)[source]

Return the pseudo with the given chemical symbol.

Parameters:
  • symbols – String with the chemical symbol of the element

  • allow_multi – By default, the method raises ValueError if multiple occurrences are found. Use allow_multi to prevent this.

Raises:

ValueError if symbol is not found or multiple occurrences are present and not allow_multi

pseudos_with_symbols(symbols)[source]

Return the pseudos with the given chemical symbols.

Raises:

ValueError if one of the symbols is not found or multiple occurrences are present.

select(condition) PseudoTable[source]

Select only those pseudopotentials for which condition is True.

Parameters:

condition – Function that accepts a Pseudo object and returns True or False.

Returns:

New PseudoTable instance with pseudos for which condition is True.

Return type:

PseudoTable

select_family(family)[source]

Return PseudoTable with element belonging to the specified family, e.g. family=”alkaline”.

select_rows(rows)[source]

Return new class:PseudoTable object with pseudos in the given rows of the periodic table. rows can be either a int or a list of integers.

select_symbols(symbols, ret_list=False)[source]

Return a PseudoTable with the pseudopotentials with the given list of chemical symbols.

Parameters:
  • symbols – str or list of symbols Prepend the symbol string with “-”, to exclude pseudos.

  • ret_list – if True a list of pseudos is returned instead of a PseudoTable

sort_by_z()[source]

Return a new PseudoTable with pseudos sorted by Z.

sorted(attrname, reverse=False)[source]

Sort the table according to the value of attribute attrname.

Returns:

PseudoTable object

Return type:

New class

to_table(filter_function=None)[source]

Return string with data in tabular form.

with_dojo_report()[source]

Select pseudos containing the DOJO_REPORT section. Return new class:PseudoTable object.

property zlist[source]

Ordered list with the atomic numbers available in the table.

class RadialFunction(mesh, values)[source]

Bases: RadialFunction

Radial Function class.

Create new instance of RadialFunction(mesh, values)

l2str(l_ang_mom)[source]

Convert the angular momentum l (int) to string.

str2l(s)[source]

Convert a string to the angular momentum l (int).

straceback()[source]

Returns a string with the traceback.

pymatgen.io.abinit.variable module

Support for Abinit input variables.

class InputVariable(name: str, value, units: str = '', valperline: int = 3)[source]

Bases: object

An Abinit input variable.

Parameters:
  • name – Name of the variable.

  • value – Value of the variable.

  • units – String specifying one of the units supported by Abinit. Default: atomic units.

  • valperline – Number of items printed per line.

property basename[source]

Return the name trimmed of any dataset index.

property dataset[source]

Return the dataset index in string form.

format_list(values, float_decimal=0)[source]

Format a list of values into a string. The result might be spread among several lines.

static format_list2d(values, float_decimal=0)[source]

Format a list of lists.

static format_scalar(val, float_decimal=0)[source]

Format a single numerical value into a string with the appropriate number of decimal.

get_value()[source]

Return the value.

property name[source]

Name of the variable.

property units[source]

Return the units.

flatten(iterable)[source]

Make an iterable flat, i.e. a 1d iterable object.