pymatgen.io.aims package

IO interface for FHI-aims.

Subpackages

Submodules

pymatgen.io.aims.inputs module

Classes for reading/manipulating/writing FHI-aims input files.

Works for aims cube objects, geometry.in and control.in

class AimsControlIn(_parameters: dict[str, ~typing.Any] = <factory>)[source]

Bases: MSONable

Class representing and FHI-aims control.in file

_parameters[source]

The parameters dictionary containing all input flags (key) and values for the control.in file

Type:

dict[str, Any]

as_dict() dict[str, Any][source]

Get a dictionary representation of the geometry.in file.

classmethod from_dict(dct: dict[str, Any]) AimsControlIn[source]

Initialize from dictionary.

Parameters:

dct (dict[str, Any]) – The MontyEncoded dictionary

Returns:

The AimsControlIn for dct

get_aims_control_parameter_str(key: str, value: Any, fmt: str) str[source]

Get the string needed to add a parameter to the control.in file

Parameters:
  • key (str) – The name of the input flag

  • value (Any) – The value to be set for the flag

  • fmt (str) – The format string to apply to the value

Returns:

The line to add to the control.in file

Return type:

str

get_content(structure: Structure | Molecule, verbose_header: bool = False, directory: str | Path | None = None) str[source]

Get the content of the file

Parameters:
  • structure (Structure or Molecule) – The structure to write the input file for

  • verbose_header (bool) – If True print the input option dictionary

  • directory – str | Path | None = The directory for the calculation,

Returns:

The content of the file for a given structure

Return type:

str

get_species_block(structure: Structure | Molecule, species_dir: str | Path) str[source]

Get the basis set information for a structure

Parameters:
  • structure (Molecule or Structure) – The structure to get the basis set information for

  • Pat (species_dir (str or) – ): The directory to find the species files in

Returns:

The block to add to the control.in file for the species

Raises:

ValueError – If a file for the species is not found

property parameters: dict[str, Any][source]

The dictionary of input parameters for control.in

write_file(structure: Structure | Molecule, directory: str | Path | None = None, verbose_header: bool = False, overwrite: bool = False) None[source]

Writes the control.in file

Parameters:
  • structure (Structure or Molecule) – The structure to write the input file for

  • directory (str or Path) – The directory to write the control.in file. If None use cwd

  • verbose_header (bool) – If True print the input option dictionary

  • overwrite (bool) – If True allow to overwrite existing files

Raises:
  • ValueError – If a file must be overwritten and overwrite is False

  • ValueError – If k-grid is not provided for the periodic structures

class AimsCube(type: str = <factory>, origin: Sequence[float] | tuple[float, float, float] = <factory>, edges: Sequence[Sequence[float]] = <factory>, points: Sequence[int] | tuple[int, int, int] = <factory>, format: str = 'cube', spin_state: int | None = None, kpoint: int | None = None, filename: str | None = None, elf_type: int | None = None)[source]

Bases: MSONable

Class representing the FHI-aims cubes

type[source]

The value to be outputted as a cube file

Type:

str

origin[source]

The origin of the cube

Type:

Sequence[float] or tuple[float, float, float]

edges[source]

Specifies the edges of a cube: dx, dy, dz dx (float): The length of the step in the x direction dy (float): The length of the step in the y direction dx (float): The length of the step in the x direction

Type:

Sequence[Sequence[float]]

points[source]

The number of points along each edge

Type:

Sequence[int] or tuple[int, int, int]

spin_state[source]

The spin-channel to use either 1 or 2

Type:

int

kpoint[source]

The k-point to use (the index of the list printed from output k_point_list)

Type:

int

filename[source]

The filename to use

Type:

str

format[source]

The format to output the cube file in: cube, gOpenMol, or xsf

Type:

str

elf_type[source]

The type of electron localization function to use ( see FHI-aims manual)

Type:

int

as_dict() dict[str, Any][source]

Get a dictionary representation of the geometry.in file.

property control_block: str[source]

Get the block of text for the control.in file of the Cube

edges: Sequence[Sequence[float]][source]
elf_type: int | None = None[source]
filename: str | None = None[source]
format: str = 'cube'[source]
classmethod from_dict(dct: dict[str, Any]) AimsCube[source]

Initialize from dictionary.

Parameters:

dct (dict[str, Any]) – The MontyEncoded dictionary

Returns:

AimsCube

kpoint: int | None = None[source]
origin: Sequence[float] | tuple[float, float, float][source]
points: Sequence[int] | tuple[int, int, int][source]
spin_state: int | None = None[source]
type: str[source]
class AimsGeometryIn(_content: str, _structure: Structure | Molecule)[source]

Bases: MSONable

Representation of an aims geometry.in file

_content[source]

The content of the input file

Type:

str

_structure[source]

The structure or molecule representation of the file

Type:

Structure or Molecule

as_dict() dict[str, Any][source]

Get a dictionary representation of the geometry.in file.

property content: str[source]

Access the contents of the file

classmethod from_dict(dct: dict[str, Any]) AimsGeometryIn[source]

Initialize from dictionary.

Parameters:

dct (dict[str, Any]) – The MontyEncoded dictionary of the AimsGeometryIn object

Returns:

The input object represented by the dict

classmethod from_file(filepath: str | Path) AimsGeometryIn[source]

Create an AimsGeometryIn from an input file.

Parameters:

filepath (str | Path) – The path to the input file (either plain text of gzipped)

Returns:

The input object represented in the file

Return type:

AimsGeometryIn

classmethod from_str(contents: str) AimsGeometryIn[source]

Create an input from the content of an input file

Parameters:

contents (str) – The content of the file

Returns:

The AimsGeometryIn file for the string contents

classmethod from_structure(structure: Structure | Molecule) AimsGeometryIn[source]

Construct an input file from an input structure.

Parameters:

structure (Structure or Molecule) – The structure for the file

Returns:

The input object for the structure

Return type:

AimsGeometryIn

property structure: Structure | Molecule[source]

Access structure for the file

write_file(directory: str | Path | None = None, overwrite: bool = False) None[source]

Write the geometry.in file

Parameters:
  • directory (str | Path | None) – The directory to write the geometry.in file

  • overwrite (bool) – If True allow to overwrite existing files

pymatgen.io.aims.outputs module

A representation of FHI-aims output (based on ASE output parser).

class AimsOutput(results: Molecule | Structure | Sequence[Molecule | Structure], metadata: dict[str, Any], structure_summary: dict[str, Any])[source]

Bases: MSONable

The main output file for FHI-aims.

AimsOutput object constructor.

Parameters:
  • results (Molecule or Structure or Sequence[Molecule or Structure]) – A list of all images in an output file

  • metadata (Dict[str, Any]) – The metadata of the executable used to perform the calculation

  • structure_summary (Dict[str, Any]) – The summary of the starting atomic structure

property aims_version: str[source]

The version of FHI-aims used for the calculation.

property all_forces: list[list[Vector3D]][source]

The forces for all images in the calculation.

as_dict() dict[str, Any][source]

Create a dict representation of the outputs for MSONable.

property band_gap: float[source]

The band gap for the final structure in the calculation.

property cbm: float[source]

The LUMO level for the final structure in the calculation.

property completed: bool[source]

Did the calculation complete.

property direct_band_gap: float[source]

The direct band gap for the final structure in the calculation.

property fermi_energy: float[source]

The Fermi energy for the final structure in the calculation.

property final_energy: float[source]

The total energy for the final structure in the calculation.

property final_structure: Structure | Molecule[source]

The final structure for the calculation.

property forces: Sequence[Vector3D] | None[source]

The forces for the final image of the calculation.

classmethod from_dict(dct: dict[str, Any]) AimsOutput[source]

Construct an AimsOutput from a dictionary.

Parameters:

dct (dict[str, Any]) – The dictionary used to create AimsOutput

Returns:

AimsOutput

classmethod from_outfile(outfile: str | Path) AimsOutput[source]

Construct an AimsOutput from an output file.

Parameters:

outfile – str | Path: The aims.out file to parse

Returns:

The AimsOutput object for the output file

classmethod from_str(content: str) AimsOutput[source]

Construct an AimsOutput from an output file.

Parameters:

content (str) – The content of the aims.out file

Returns:

The AimsOutput for the output file content

get_results_for_image(image_ind: int) Structure | Molecule[source]

Get the results dictionary for a particular image or slice of images.

Parameters:

image_ind (int) – The index of the image to get the results for

Returns:

The results of the image with index images_ind

property initial_structure: Structure | Molecule[source]

The initial structure for the calculations.

property metadata: dict[str, Any][source]

The system metadata.

property n_images: int[source]

The number of images in results.

property stress: Matrix3D[source]

The stress for the final image of the calculation.

property stresses: Sequence[Matrix3D] | None[source]

The atomic virial stresses for the final image of the calculation.

property structure_summary: dict[str, Any][source]

The summary of the material/molecule that the calculations represent.

property structures: Sequence[Structure | Molecule][source]

All images in the output file.

property vbm: float[source]

The HOMO level for the final structure in the calculation.

pymatgen.io.aims.parsers module

AIMS output parser, taken from ASE with modifications.

class AimsOutCalcChunk(lines: list[str], header: AimsOutHeaderChunk)[source]

Bases: AimsOutChunk

A part of the aims.out file corresponding to a single structure.

Construct the AimsOutCalcChunk.

Parameters:
  • lines (list[str]) – The lines used for the structure

  • header (.AimsOutHeaderChunk) – A summary of the relevant information from the aims.out header

property E_f: float | None[source]

The Fermi energy

property cbm: float[source]

The conduction band minimnum

property converged: bool[source]

True if the calculation is converged

property coords: list[Vector3D][source]

The cartesian coordinates of the atoms

property dielectric_tensor: Matrix3D | None[source]

The dielectric tensor from the aims.out file.

property dipole: Vector3D | None[source]

The electric dipole moment from the aims.out file.

property direct_gap: float[source]

The direct bandgap

property electronic_temperature: float[source]

The electronic temperature for the chunk.

property energy: float[source]

The energy from the aims.out file.

property forces: np.array[Vector3D] | None[source]

The forces from the aims.out file.

property free_energy: float | None[source]

The free energy of the calculation

property gap: float[source]

The band gap

property hirshfeld_atomic_dipoles: Sequence[Vector3D] | None[source]

The Hirshfeld atomic dipoles of the system

property hirshfeld_charges: Sequence[float] | None[source]

The Hirshfeld charges of the system

property hirshfeld_dipole: None | Vector3D[source]

The Hirshfeld dipole of the system

property hirshfeld_volumes: Sequence[float] | None[source]

The Hirshfeld atomic dipoles of the system

property initial_lattice: Lattice | None[source]

The initial Lattice of the structure

property initial_structure: Structure | Molecule[source]

The initial structure for the calculation

property is_metallic: bool[source]

Is the system is metallic.

property k_point_weights: Sequence[float][source]

The k-point weights for the calculation.

property k_points: Sequence[Vector3D][source]

All k-points listed in the calculation.

property lattice: Lattice[source]

The Lattice object for the structure

property magmom: float | None[source]

The magnetic moment of the structure

property n_atoms: int[source]

The number of atoms in the structure

property n_bands: int[source]

The number of Kohn-Sham states for the chunk.

property n_electrons: int[source]

The number of electrons for the chunk.

property n_iter: int | None[source]

The number of steps needed to converge the SCF cycle for the chunk.

property n_k_points: int[source]

The number of k_ppoints for the calculation.

property n_spins: int[source]

The number of spin channels for the chunk.

property polarization: Vector3D | None[source]

The polarization vector from the aims.out file.

property results: dict[str, Any][source]

Convert an AimsOutChunk to a Results Dictionary.

property species: list[str][source]

The list of atomic symbols for all atoms in the structure

property stress: Matrix3D | None[source]

The stress from the aims.out file and convert to kbar.

property stresses: np.array[Matrix3D] | None[source]

The stresses from the aims.out file and convert to kbar.

property structure: Structure | Molecule[source]

The pytmagen SiteCollection of the chunk.

property vbm: float[source]

The valance band maximum

property velocities: list[Vector3D][source]

The velocities of the atoms

class AimsOutChunk(lines: list[str] = <factory>)[source]

Bases: object

Base class for AimsOutChunks.

lines[source]

The list of all lines in the chunk

Type:

list[str]

lines: list[str][source]
parse_scalar(property: str) float | None[source]

Parse a scalar property from the chunk.

Parameters:

property (str) – The property key to parse

Returns:

The scalar value of the property or None if not found

reverse_search_for(keys: list[str], line_start: int = 0) int[source]

Find the last time one of the keys appears in self.lines.

Parameters:
  • keys (list[str]) – The key strings to search for in self.lines

  • line_start (int) – The lowest index to search for in self.lines

Returns:

The last time one of the keys appears in self.lines

search_for_all(key: str, line_start: int = 0, line_end: int = -1) list[int][source]

Find the all times the key appears in self.lines.

Parameters:
  • key (str) – The key string to search for in self.lines

  • line_start (int) – The first line to start the search from

  • line_end (int) – The last line to end the search at

Returns:

All times the key appears in the lines

class AimsOutHeaderChunk(lines: list[str] = <factory>, _cache: dict[str, ~typing.Any] = <factory>)[source]

Bases: AimsOutChunk

The header of the aims.out file containing general information.

property aims_uuid: str[source]

The aims-uuid for the calculation.

property build_type: list[str][source]

The optional build flags passed to cmake.

property c_compiler: str | None[source]

The C compiler used to make FHI-aims.

property c_compiler_flags: str | None[source]

The C compiler flags used to make FHI-aims.

property commit_hash: str[source]

The commit hash for the FHI-aims version.

property electronic_temperature: float[source]

The electronic temperature for the chunk.

property fortran_compiler: str | None[source]

The fortran compiler used to make FHI-aims.

property fortran_compiler_flags: str | None[source]

The fortran compiler flags used to make FHI-aims.

property header_summary: dict[str, Any][source]

Dictionary summarizing the information inside the header.

property initial_charges: Sequence[float][source]

The initial charges for the structure

property initial_lattice: Lattice | None[source]

The initial lattice vectors from the aims.out file.

property initial_magnetic_moments: Sequence[float][source]

The initial magnetic Moments

property initial_structure: Structure | Molecule[source]

The initial structure

Using the FHI-aims output file recreate the initial structure for the calculation.

property is_md: bool[source]

Is the output for a molecular dynamics calculation?

property is_relaxation: bool[source]

Is the output for a relaxation?

property k_point_weights: Sequence[float][source]

The k-point weights for the calculation.

property k_points: Sequence[Vector3D][source]

All k-points listed in the calculation.

lines: list[str][source]
property linked_against: list[str][source]

Get all libraries used to link the FHI-aims executable.

property metadata_summary: dict[str, list[str] | str | None][source]

Dictionary containing all metadata for FHI-aims build.

property n_atoms: int[source]

The number of atoms for the material.

property n_bands: int | None[source]

The number of Kohn-Sham states for the chunk.

property n_electrons: int | None[source]

The number of electrons for the chunk.

property n_k_points: int | None[source]

The number of k_ppoints for the calculation.

property n_spins: int | None[source]

The number of spin channels for the chunk.

property version_number: str[source]

The commit hash for the FHI-aims version.

exception AimsParseError(message: str)[source]

Bases: Exception

Exception raised if an error occurs when parsing an Aims output file.

Initialize the error with the message, message

exception ParseError[source]

Bases: Exception

Parse error during reading of a file

check_convergence(chunks: list[AimsOutCalcChunk], non_convergence_ok: bool = False) bool[source]

Check if the aims output file is for a converged calculation.

Parameters:
  • chunks (list[.AimsOutCalcChunk]) – The list of chunks for the aims calculations

  • non_convergence_ok (bool) – True if it is okay for the calculation to not be converged

  • chunks – list[AimsOutCalcChunk]:

  • non_convergence_ok – bool: (Default value = False)

Returns:

True if the calculation is converged

get_aims_out_chunks(content: str | TextIOWrapper, header_chunk: AimsOutHeaderChunk) Generator[source]

Yield unprocessed chunks (header, lines) for each AimsOutChunk image.

Parameters:
  • content (str or TextIOWrapper) – the content to parse

  • header_chunk (AimsOutHeaderChunk) – The AimsOutHeader for the calculation

Yields:

The next AimsOutChunk

get_header_chunk(content: str | TextIOWrapper) AimsOutHeaderChunk[source]

Get the header chunk for an output

Parameters:

content (str or TextIOWrapper) – the content to parse

Returns:

The AimsHeaderChunk of the file

get_lines(content: str | TextIOWrapper) list[str][source]

Get a list of lines from a str or file of content

Parameters:

content – the content of the file to parse

Returns:

The list of lines

read_aims_header_info(filename: str | Path) tuple[dict[str, None | list[str] | str], dict[str, Any]][source]

Read the FHI-aims header information.

Parameters:

filename (str or Path) – The file to read

Returns:

The metadata for the header of the aims calculation

read_aims_header_info_from_content(content: str) tuple[dict[str, list[str] | None | str], dict[str, Any]][source]

Read the FHI-aims header information.

Parameters:

content (str) – The content of the output file to check

Returns:

The metadata for the header of the aims calculation

read_aims_output(filename: str | Path, index: int | slice = -1, non_convergence_ok: bool = False) Structure | Molecule | Sequence[Structure | Molecule][source]

Import FHI-aims output files with all data available.

Includes all structures for relaxations and MD runs with FHI-aims

Parameters:
  • filename (str or Path) – The file to read

  • index (int or slice) – The index of the images to read

  • non_convergence_ok (bool) – True if the calculations do not have to be converged

Returns:

The list of images to get

read_aims_output_from_content(content: str, index: int | slice = -1, non_convergence_ok: bool = False) Structure | Molecule | Sequence[Structure | Molecule][source]

Read and aims output file from the content of a file

Parameters:
  • content (str) – The content of the file to read

  • index – int | slice: (Default value = -1)

  • non_convergence_ok – bool: (Default value = False)

Returns:

The list of images to get