pymatgen.io.aims package
IO interface for FHI-aims.
Subpackages
- pymatgen.io.aims.sets package
- Submodules
- pymatgen.io.aims.sets.base module
AimsInputGenerator
AimsInputGenerator.user_params
AimsInputGenerator.user_kpoints_settings
AimsInputGenerator.use_structure_charge
AimsInputGenerator.d2k()
AimsInputGenerator.d2k_recip_cell()
AimsInputGenerator.get_input_set()
AimsInputGenerator.get_parameter_updates()
AimsInputGenerator.k2d()
AimsInputGenerator.use_structure_charge
AimsInputGenerator.user_kpoints_settings
AimsInputGenerator.user_params
AimsInputSet
recursive_update()
- pymatgen.io.aims.sets.bs module
- pymatgen.io.aims.sets.core module
- pymatgen.io.aims.sets.magnetism module
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, Any] = <factory>)[source]
Bases:
MSONable
An 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]
- classmethod from_dict(dct: dict[str, Any]) Self [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.
- get_species_block(structure: Structure | Molecule, basis_set: str | dict[str, str], species_dir: str | Path | None = None) str [source]
Get the basis set information for a structure
- Parameters:
structure (Molecule or Structure) – The structure to get the basis set information for
basis_set (str | dict[str, str]) – a name of a basis set (light, tight…) or a mapping from site labels to basis set names. The name of a basis set can either correspond to the subfolder in defaults_2020 folder or be a full path from the FHI-aims/species_defaults directory.
species_dir (str | Path | None) – The base species directory
- Returns:
The block to add to the control.in file for the species
- Raises:
ValueError – If a file for the species is not found
- write_file(structure: Structure | Molecule, directory: str | Path | None = None, verbose_header: bool = False, overwrite: bool = False) None [source]
Write the control.in file.
- Parameters:
- 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] | Tuple3Floats = <factory>, edges: Sequence[Sequence[float]] = <factory>, points: Sequence[int] | Tuple3Ints = <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
The FHI-aims cubes.
- 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]]
- kpoint[source]
The k-point to use (the index of the list printed from output k_point_list)
- Type:
int
- elf_type[source]
The type of electron localization function to use ( see FHI-aims manual)
- Type:
int
- class AimsGeometryIn(_content: str, _structure: Structure | Molecule)[source]
Bases:
MSONable
Representation of an aims geometry.in file.
- classmethod from_dict(dct: dict[str, Any]) Self [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) Self [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:
- classmethod from_str(contents: str) Self [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) Self [source]
Construct an input file from an input structure.
- class AimsSpeciesFile(data: str = '', label: str | None = None)[source]
Bases:
object
An FHI-aims single species’ defaults file.
- classmethod from_dict(dct: dict[str, Any]) Self [source]
Deserialization of the AimsSpeciesFile object
- classmethod from_element_and_basis_name(element: str, basis: str, *, species_dir: str | Path | None = None, label: str | None = None) Self [source]
Initialize from element and basis names.
- Parameters:
element (str) – the element name (not to confuse with the species)
basis (str) – the directory in which the species’ defaults file is located relative to the root species_defaults (or species_defaults/defaults_2020) directory.`.
label (str) – A string representing the name of species. If not specified, then equal to element
- Returns:
AimsSpeciesFile
- class SpeciesDefaults(labels: Sequence[str], basis_set: str | dict[str, str], *, species_dir: str | Path | None = None, elements: dict[str, str] | None = None)[source]
Bases:
list
,MSONable
A list containing a set of species’ defaults objects with methods to read and write them to files
- Parameters:
labels (list[str]) – a list of labels, for which to build species’ defaults
basis_set (str | dict[str, str]) – a name of a basis set (light, tight…) or a mapping from site labels to basis set names. The name of a basis set can either correspond to the subfolder in defaults_2020 folder or be a full path from the FHI-aims/species_defaults directory.
species_dir (str | Path | None) – The base species directory
elements (dict[str, str] | None) – a mapping from site labels to elements. If some label is not in this mapping, it coincides with an element.
- classmethod from_dict(dct: dict[str, Any]) SpeciesDefaults [source]
Deserialization of the SpeciesDefaults object
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.
- Parameters:
- property direct_band_gap: float[source]
The direct band gap for the final structure in the calculation.
- property forces: Sequence[Vector3D] | None[source]
The forces for the final image of the calculation.
- classmethod from_dict(dct: dict[str, Any]) Self [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) Self [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) Self [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 stresses: Sequence[Matrix3D] | None[source]
The atomic virial stresses for the final image of 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 hirshfeld_atomic_dipoles: Sequence[Vector3D] | None[source]
The Hirshfeld atomic dipoles of the system.
- property hirshfeld_volumes: Sequence[float] | None[source]
The Hirshfeld atomic dipoles of the system.
- property initial_structure: Structure | Molecule[source]
The initial structure for the calculation.
- class AimsOutChunk(lines: list[str] = <factory>)[source]
Bases:
object
Base class for AimsOutChunks.
- 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, Any] = <factory>)[source]
Bases:
AimsOutChunk
The header of the aims.out file containing general information.
- 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_lattice: Lattice | None[source]
The initial lattice vectors from the aims.out file.
- property initial_structure: Structure | Molecule[source]
The initial structure.
Using the FHI-aims output file recreate the initial structure for the calculation.
- 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.
- 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