Utilities

Methods for structure manipulations

class scine_chemoton.utilities.reactive_complexes.ReactiveComplexes[source]

The base class for all reactive complex generators.

class Options[source]

Options attribute to be implemented by child classes

set_options(option_dict)[source]

Sets the options for the ReactiveComplexes from a dictionary. Generates a warning if an option is unknown.

Return type:

None

Parameters:
option_dictDict[str, Any]

Dictionary with options to be used for generating reactive complexes.

class scine_chemoton.utilities.reactive_complexes.inter_reactive_complexes.InterReactiveComplexes[source]

Class to generate reactive complexes from two structures.

class Options[source]

The options for the InterReactiveComplexes

multiple_attack_points
bool

Whether to consider multiple attack points for each active centers involved in intermolecular reactive pairs or just one. (default: True)

number_rotamers
int

The number of rotamers to be generated for reactive complexes with at least one active center being an atom. (default: 2)

number_rotamers_two_on_two
int

The number of rotamers to be generated for reactive complexes with both reactive centers being diatomic. (default: 1)

generate_reactive_complexes(structure1, structure2, reactive_inter_coords)[source]

Generates a set of reactive complexes for two given structures arising from the given intermolecular reactive pairs.

Return type:

Generator[Tuple[List[Tuple[int, int]], ndarray, ndarray, float, float], None, None]

Parameters:
structure1, structure2scine_database.Structure (Scine::Database::Structure)

The two structures for which a set of reactive complexes is to be generated. The structures have to be linked to a collection.

reactive_inter_coordsList[List[Tuple[int, int]]]

A list of intermolecular reactive atom pairs corresponding to one trial reaction coordinate. Each reactive pair tuple has to be ordered such that its first element belongs to structure1 and the second to structure2. The indices are expected to refer to be on structure level, i.e. the first atom of structure2 has index 0 and not index n_atoms(structure1).

Yields:
inter_coordTuple[Tuple[Tuple[int]]

Tuple of Tuples of one or two atom pairs composing the reactive atoms of the interstructural component of this reactive complex reaction. First atom per pair belongs to structure1, second to structure2.

align1, align2np.array

Rotation matrices aligning the two sites along the x-axis (rotations assume that the geometric mean of the reactive atoms of each structure is translated into the origin)

xrotfloat

Angle of rotation around the x-axis

spreadfloat

Spread to be applied along the x-axis between the two structures.

options: Options
set_options(option_dict)

Sets the options for the ReactiveComplexes from a dictionary. Generates a warning if an option is unknown.

Return type:

None

Parameters:
option_dictDict[str, Any]

Dictionary with options to be used for generating reactive complexes.

scine_chemoton.utilities.reactive_complexes.inter_reactive_complexes.assemble_reactive_complex(atoms1, atoms2, lhs_list, rhs_list, x_alignment_0=None, x_alignment_1=None, x_rotation=0.0, x_spread=2.0, displacement=0.0)[source]

Assembles a reactive complex from the parameters generated by the InterReactiveComplexes class.

Return type:

Tuple[AtomCollection, List[int], List[int]]

Parameters:
atoms1, atoms2utils.AtomCollection

The atoms of both structures, that are to be combined in the reactive complex. atoms1 refers to the LHS and atoms2 to the RHS.

lhs_list, rhs_listList[int]

Indices of the reactive sites within the reactive complex. The lhs_list should correspond to atoms1 and rhs_list to atoms2.

x_alignment_0List[float], length=9

In case of two structures building the reactive complex, this option describes a rotation of the first structure (index 0) that aligns the reaction coordinate along the x-axis (pointing towards +x). The rotation assumes that the geometric mean position of all atoms in the reactive site (lhs_list) is shifted into the origin.

x_alignment_1List[float], length=9

In case of two structures building the reactive complex, this option describes a rotation of the second structure (index 1) that aligns the reaction coordinate along the x-axis (pointing towards -x). The rotation assumes that the geometric mean position of all atoms in the reactive site (rhs_list) is shifted into the origin.

x_rotationfloat

In case of two structures building the reactive complex, this option describes a rotation angle around the x-axis of one of the two structures after x_alignment_0 and x_alignment_1 have been applied.

x_spreadfloat

In case of two structures building the reactive complex, this option gives the distance by which the two structures are moved apart along the x-axis after x_alignment_0, x_alignment_1, and x_rotation have been applied.

displacementfloat

In case of two structures building the reactive complex, this option adds a random displacement to all atoms (random direction, random length). The maximum length of this displacement (per atom) is set to be the value of this option.

Returns:
utils.AtomCollection

The reactive complex structure.

List[int], List[int]

The LHS and RHS lists with indices adapted to the reactive complex structure.

class scine_chemoton.utilities.reactive_complexes.lebedev_sphere.LebedevSphere[source]

A 5810 point Lebedev unit sphere and the nearest neighbors for each point.

property nearest_neighbors: List[List[int]]

List[List[int]] : The nearest neighbors of all points on the sphere.

property points: ndarray

np.ndarray : The points of the Lebedev sphere.

class scine_chemoton.utilities.reactive_complexes.unit_circle.UnitCircle[source]

A class holding a unit circle of 100 points.

property nearest_neighbors: List[List[int]]

List[List[int]] : The nearest neighbors of all points on the circle.

property points: ndarray

np.ndarray : The points of the circle.

scine_chemoton.utilities.insert_initial_structure.insert_initial_structure(database, molecule_path, charge, multiplicity, model, label=<Label.USER_GUESS: 1>, job=<scine_database.Job object>, settings=scine_utilities.ValueCollection({}), start_concentration=None)[source]

Insert a structure to the database and set up a calculation working on it.

Parameters:
databasedb.Manager

Database to use.

molecule_pathUnion[str, utils.AtomCollection]

Atom collection or path to the xyz file with the structure to be inserted.

chargeint

Charge of the structure.

multiplicityint

Multiplicity of the structure.

modeldb.Model

Model to be used for the calculation.

labeldb.Label, optional

Label of the inserted structure, by default db.Label.MINIMUM_GUESS.

jobdb.Job, optional

Job to be performed on the initial structure, by default db.Job(‘scine_geometry_optimization’).

settingsutils.ValueCollection, optional

Job settings, by default none.

start_concentrationfloat

The start concentratoin of the compound that will be generated from this structure.

Returns:
db.Structure, db.Calculation

The inserted structure and the calculation generated for it

class scine_chemoton.utilities.masm.ComponentDistanceTuple(mol_idx, components, distance)

Create new instance of ComponentDistanceTuple(mol_idx, components, distance)

components

Alias for field number 1

count(value, /)

Return number of occurrences of value.

distance

Alias for field number 2

index(value, start=0, stop=sys.maxsize, /)

Return first index of value.

Raises ValueError if the value is not present.

mol_idx

Alias for field number 0

scine_chemoton.utilities.masm.deserialize_molecules(structure)[source]

Retrieves all molecules stored for a structure

Return type:

List[Molecule]

Parameters:
structuredb.Structure

The structure whose contained molecules to deserialize

Returns:
moleculesList[masm.Molecule]

A list of all the molecules contained in the database structure

scine_chemoton.utilities.masm.distinct_atoms(mol, h_only)[source]

Generates a list of distinct atom indices

Return type:

List[int]

Parameters:
molmasm.Molecule

A molecule whose atoms to list distinct atoms for

h_onlybool

Whether to only apply ranking deduplication to hydrogen atoms

Returns:
componentsList[int]

A list of ranking-distinct atoms

scine_chemoton.utilities.masm.distinct_components(mol, h_only)[source]

Generates a flat map of atom index to component identifier

Return type:

List[int]

Parameters:
molmasm.Molecule

A molecule whose atoms to generate distinct components for

h_onlybool

Whether to only apply ranking deduplication to hydrogen atoms

Returns:
componentsList[int]

A flat per-atom index mapping to a component index. Contains only sequential numbers starting from zero.

scine_chemoton.utilities.masm.distinguish_components(components, map_unary)[source]

Splits components by the result of a unary mapping function

Return type:

List[int]

Parameters:
componentsList[int]

A per-index mapping to a component index. Must contain only sequential numbers starting from zero.

map_unaryCallable[[int], Any]

A unary callable that is called with an index, not a component index, yielding some comparable type. Components of indices are then split by matching results of invocations of this callable.

Returns:
componentsList[int]

A per-index mapping to a component index. Contains only sequential numbers starting from zero.

scine_chemoton.utilities.masm.get_atom_pairs(structure, distance_bounds, prune='None', superset=None)[source]

Gets a list of all atom pairs whose graph distance is smaller or equal to max_graph_distance and larger or equal to min_graph_distance on the basis of the interpreted graph representation.

Return type:

Set[Tuple[int, int]]

Parameters:
structuredb.Structure

The structure that is investigated

distance_boundsTuple[int, int]

The minimum and maximum distance between two points that is allowed so that they are considered a valid atom pair.

prunestr

Whether to prune atom pairings by Molassembler’s ranking distinct atoms descriptor. Allowed values: ‘None’, ‘Hydrogen’, ‘All’

supersetOptional[Set[Tuple[int, int]]]

Optional superset of pairs to filter. If set, will filter the passed set. Otherwise, generates atom pairings from all possible pairs in the molecule.

Returns:
pairsSet[Tuple[int, int]]

The indices of valid atom pairs.

scine_chemoton.utilities.masm.make_sorted_pair(a, b)[source]
Return type:

Tuple[int, int]

scine_chemoton.utilities.masm.mol_from_cbor(cbor_str)[source]

Convert base-64 encoded CBOR to a molassembler Molecule

Converts a single base-64 encoded CBOR string (no ‘;’ separator as stored in the database) into a molecule

Return type:

Molecule

Parameters:
cbor_strstr

String to deserialize into a Molecule

Returns:
moleculemasm.Molecule

The deserialized molecule

scine_chemoton.utilities.masm.mol_to_cbor(mol)[source]

Convert a molecule into a base-64 encoded CBOR serialization

Return type:

str

Parameters:
molmasm.Molecule

Molecule to serialize

Returns:
serializationstr

The string-serialized molecule representation

scine_chemoton.utilities.masm.mols_from_properties(structure, properties)[source]

Generate all molecules based on atomic positions in a structure and the bond orders stored in attache properties.

Return type:

Optional[List[Molecule]]

Parameters:
structuredb.Structure

The structure whose contained molecule(s) to analyze.

propertiesdb.Collection

The collection holding all properties.

Returns:
moleculesList[masm.Molecule]

A list of all the molecules contained in the database structure.

scine_chemoton.utilities.masm.pruned_atom_pairs(molecules, idx_map, distance_bounds, prune)[source]
Return type:

Set[Tuple[int, int]]

scine_chemoton.utilities.masm.unpruned_atom_pairs(molecules, idx_map, distance_bounds)[source]

Helper function to generate the set of unpruned atom pairs

Return type:

Set[Tuple[int, int]]

Database Object Wrappers

class scine_chemoton.utilities.db_object_wrappers.aggregate_cache.AggregateCache(manager, electronic_model, hessian_model, only_electronic=False)[source]

A cache for aggregate objects.

Parameters:
manager

The database manager.

electronic_model

The electronic structure model from which the electronic energy is taken.

hessian_model

The electronic structure model with which the geometrical hessian was calculated.

only_electronic

If true, only the electronic energies are used to determine the thermodynamics.

get_electronic_model()[source]
get_or_produce(aggregate_id)[source]

Get an instance of the aggregate wrapper corresponding the reaction ID. The instance may be newly constructed or retrieved from the cache.

Return type:

Aggregate

Parameters:
aggregate_iddb.ID

The aggregate database ID.

Returns:
Aggregate

The aggregate wrapper instance.

class scine_chemoton.utilities.db_object_wrappers.aggregate_wrapper.Aggregate(aggregate_id, manager, electronic_model, hessian_model, only_electronic=False)[source]

A class that wraps around db.Compound and db.Flask to provide easy access to all thermodynamic contributions. Free energies are calculated with the harmonic oscillator, rigid rotor, particle in a box model according to the given thermodynamic reference state (allows caching or works on the fly if necessary).

Parameters:
aggregate_iddb.ID

The database ID.

manager

The database manager.

electronic_model

The electronic structure model from which the electronic energy is taken.

hessian_model

The electronic structure model with which the geometrical hessian was calculated.

only_electronic

If true, only the electronic energies are used to determine the thermodynamics.

analyze()

Returns true if the database object is set to analyze. Returns false otherwise.

complete()

Returns true if a free energy approximation is available for the ensemble. Returns false otherwise.

explore()

Returns true if the database object is set to explore. Returns false otherwise.

get_aggregate_type()[source]

Getter for the aggregate type.

Return type:

CompoundOrFlask

get_concentration_flux()[source]

Getter for the aggregates vertex flux.

Return type:

float

get_db_id()

Returns the database ID.

Return type:

ID

get_db_object()

Return the database object.

get_electronic_model()

Getter for the electronic structure mode with which the electronic energies are evaluated.

get_element_count()[source]

Getter for a dictionary with the element and its count in the aggregate.

Return type:

Dict[str, int]

get_enthalpy(reference_state)[source]

Getter for the aggregate’s enthalpy.

Return type:

Optional[float]

get_entropy(reference_state)[source]

Getter for the aggregate’s entropy.

Return type:

Optional[float]

get_free_energy(reference_state)[source]

Getter for the aggregate’s free energy.

Return type:

Optional[float]

get_lowest_n_structures(n, energy_cut_off, reference_state)

Getter for the lowest n structures according to their energy.

Return type:

List[ID]

Parameters:
nint

The number of structures requested.

energy_cut_offfloat

The energy cutoff after which structures are no longer considered for the list.

reference_stateReferenceState

The thermodynamic reference state for the free energy approximation.

Returns:
List[Tuple[db.ID, float]]

The lowest n structures according to their free energy approximation.

get_manager()

Returns the database manager.

Return type:

Manager

get_model_combination()

Getter for the model combination of electronic and Hessian model.

get_sorted_structure_list(reference_state)[source]

Getter a list of structures and their energies sorted by energy (ascending).

Return type:

List[Tuple[ID, float]]

Parameters:
reference_stateReferenceState

The thermodynamic reference state for the free energy approximation.

get_starting_concentration()[source]

Getter for the aggregates starting concentration.

Return type:

float

class scine_chemoton.utilities.db_object_wrappers.ensemble_wrapper.Ensemble(db_id, manager, electronic_model, hessian_model, only_electronic=False)[source]

Base class of the wrapper for an ensemble of structures that is derived from either a database compound, flask, or reaction.

Parameters:
db_iddb.ID

The database ID.

managerdb.Manager

The database manager.

electronic_modeldb.Model

The electronic structure model for the electronic energy contribution to the free energy.

hessian_modeldb.Model

The electronic structure model with which the vibrational (Hessian) and other free energy corrections are calculated.

only_electronicbool

If true, only the electronic energy is considered for the free energies.

analyze()[source]

Returns true if the database object is set to analyze. Returns false otherwise.

complete()[source]

Returns true if a free energy approximation is available for the ensemble. Returns false otherwise.

explore()[source]

Returns true if the database object is set to explore. Returns false otherwise.

get_db_id()[source]

Returns the database ID.

Return type:

ID

get_db_object()[source]

Return the database object.

get_electronic_model()[source]

Getter for the electronic structure mode with which the electronic energies are evaluated.

abstract get_element_count()[source]

Getter for a dictionary with the element and its count in the aggregate.

Return type:

Dict[str, int]

get_lowest_n_structures(n, energy_cut_off, reference_state)[source]

Getter for the lowest n structures according to their energy.

Return type:

List[ID]

Parameters:
nint

The number of structures requested.

energy_cut_offfloat

The energy cutoff after which structures are no longer considered for the list.

reference_stateReferenceState

The thermodynamic reference state for the free energy approximation.

Returns:
List[Tuple[db.ID, float]]

The lowest n structures according to their free energy approximation.

get_manager()[source]

Returns the database manager.

Return type:

Manager

get_model_combination()[source]

Getter for the model combination of electronic and Hessian model.

abstract get_sorted_structure_list(_)[source]

Getter a list of structures and their energies sorted by energy (ascending).

Return type:

List[Tuple[ID, float]]

Parameters:
_ReferenceState

The thermodynamic reference state for the free energy approximation.

class scine_chemoton.utilities.db_object_wrappers.reaction_cache.ReactionCache(manager, electronic_model, hessian_model, aggregate_cache=None, only_electronic=False)[source]

A cache for reaction objects.

Parameters:
manager

The database manager.

electronic_model

The electronic structure model from which the electronic energy is taken.

hessian_model

The electronic structure model with which the geometrical hessian was calculated.

aggregate_cacheOptional[AggregateCache]

An optional aggregate cache with which the aggregates taking part in the reaction are constructed.

only_electronic

If true, only the electronic energies are used to determine the thermodynamics.

get_or_produce(reaction_id)[source]

Get an instance of the reaction wrapper corresponding the reaction ID. The instance may be newly constructed or retrieved from the cache.

Return type:

Reaction

Parameters:
reaction_iddb.ID

The reaction database ID.

Returns:
Reaction

The reaction wrapper instance.

class scine_chemoton.utilities.db_object_wrappers.reaction_wrapper.Reaction(reaction_id, manager, electronic_model, hessian_model, aggregate_cache=None, only_electronic=False)[source]

A class that wraps around db.Reaction to provide easy access to all thermodynamic contributions (barriers, reaction energies etc.). Free energies are calculated with the harmonic oscillator, rigid rotor, particle in a box model according to the given thermodynamic reference state (caching + on the fly if necessary).

Parameters:
reaction_iddb.ID

The database ID.

manager

The database manager.

electronic_model

The electronic structure model from which the electronic energy is taken.

hessian_model

The electronic structure model with which the geometrical hessian was calculated.

aggregate_cache

A cache of already existing aggregates (optional).

only_electronic

If true, only the electronic energies are used to determine the thermodynamics (optional).

analyze()[source]

Returns true if the database object is set to analyze. Returns false otherwise.

Return type:

bool

barrierless(reference_state)[source]

Checks if the reaction has a valid barrier-less elementary step.

Return type:

bool

Parameters:
reference_stateReferenceState

The reference state.

Returns:
bool

True if the reaction has a barrier-less elementary step.

circle_reaction()[source]
complete()[source]

Returns true if a free energy approximation is available for the ensemble. Returns false otherwise.

Return type:

bool

explore()

Returns true if the database object is set to explore. Returns false otherwise.

static get_arrhenius_prefactor(reference_state)[source]

Getter for the factor kBT/h from Eyring’s transition state theory.

Return type:

float

Parameters:
reference_stateReferenceState

The reference state for the temperature.

Returns:
float

The factor kBT/h.

get_db_id()

Returns the database ID.

Return type:

ID

get_db_object()

Return the database object.

get_electronic_model()

Getter for the electronic structure mode with which the electronic energies are evaluated.

get_element_count()[source]

Getter for a dictionary with the element and its count in the aggregate.

Return type:

Dict[str, int]

get_free_energy_of_activation(reference_state, in_j_per_mol=False)[source]

Getter for the free energy of activation/barriers as a tuple for lhs and rhs.

Return type:

Tuple[Optional[float], Optional[float]]

Parameters:
reference_stateReferenceState

The reference state (temperature, and pressure)

in_j_per_molbool, optional

If true, the barriers are returned in J/mol (NOT kJ/mol), by default False

Returns:
Tuple[Optional[float], Optional[float]]

A tuple for the lhs and rhs barriers. Returns None if the energies are incomplete. For barrier-less reactions, one barrier will be the reaction energy, the other zero.

get_lhs_aggregates()[source]

Getter for the aggregate wrappers of the reaction’s LHS.

Return type:

List[Aggregate]

Returns:
List[Aggregate]

The LHS aggregates.

get_lhs_free_energy(reference_state)[source]

Getter for the total free energy of the LHS.

Return type:

Optional[float]

Parameters:
reference_stateReferenceState

The reference state.

Returns:
Optional[float]

The total free energy if available. Otherwise, None.

get_lowest_n_steps(n, energy_cut_off, reference_state)[source]

Getter for the n elementary steps with the lowest energy transition state.

Return type:

List[ID]

Parameters:
nint

The number of elementary steps.

energy_cut_offfloat

An energy cutoff after which no additional elementary steps are considered for the list.

reference_stateReferenceState

The thermodynamic reference state used to calculate the free energies.

Returns:
List[db.ID]

The n elementary step IDs with the lowest energy transition state.

get_lowest_n_structures(n, energy_cut_off, reference_state)

Getter for the lowest n structures according to their energy.

Return type:

List[ID]

Parameters:
nint

The number of structures requested.

energy_cut_offfloat

The energy cutoff after which structures are no longer considered for the list.

reference_stateReferenceState

The thermodynamic reference state for the free energy approximation.

Returns:
List[Tuple[db.ID, float]]

The lowest n structures according to their free energy approximation.

get_manager()

Returns the database manager.

Return type:

Manager

get_model_combination()

Getter for the model combination of electronic and Hessian model.

get_reaction_free_energy(reference_state, in_j_per_mol=False)[source]

Getter for the reaction free energy.

Return type:

Optional[float]

Parameters:
reference_stateReferenceState

The reference state (temperature, and pressure)

in_j_per_molbool, optional

If true, the barriers are returned in J/mol (NOT kJ/mol), by default False

Returns:
Optional[float]

The reaction free energy. May return None, if the energies are incomplete.

get_rhs_aggregates()[source]

Getter for the aggregate wrappers of the reaction’s RHS.

Return type:

List[Aggregate]

Returns:
List[Aggregate]

The RHS aggregates.

get_rhs_free_energy(reference_state)[source]

Getter for the total free energy of the RHS.

Return type:

Optional[float]

Parameters:
reference_stateReferenceState

The reference state.

Returns:
Optional[float]

The total free energy if available. Otherwise, None.

get_sorted_structure_list(reference_state)[source]

Getter a list of structures and their energies sorted by energy (ascending).

Return type:

List[Tuple[ID, float]]

Parameters:
_ReferenceState

The thermodynamic reference state for the free energy approximation.

get_sorted_structure_step_list(reference_state)[source]
Return type:

List[Tuple[ID, ID, float]]

get_transition_state_free_energy(reference_state)[source]

Getter fo the free energy of the transition state ensemble in Hartree.

Return type:

Optional[float]

Parameters:
reference_state

The reference state (temperature, and pressure)

Returns:
The free energy of the transition in Hartree.
get_ts_theory_rate_constants(reference_state)[source]

Getter for the transition state theory based reaction rate constants.

Return type:

Tuple[Optional[float], Optional[float]]

Parameters:
reference_stateReferenceState

The reference state (temperature, and pressure)

Returns:
Tuple[Optional[float], Optional[float]]

The transition state theory based reaction rate constants as a tuple for lhs and rhs. May return Tuple[None, None] if the energies are incomplete.

class scine_chemoton.utilities.db_object_wrappers.thermodynamic_properties.PlaceHolderReferenceState[source]

Place-holder reference state. This can be used as a replacement for None default arguments that must be replaced at a later point.

class scine_chemoton.utilities.db_object_wrappers.thermodynamic_properties.ReferenceState(t, p, sym=1)[source]

Reference state for thermodynamic calculations (harmonic vib., rigid rotor, particle in a box). Requires pressure, temperature, and optionally the symmetry of the molecule.

Parameters:
tfloat

The temperature in K.

pfloat

The pressure in Pa.

symint

The symmetry number of the molecule.

class scine_chemoton.utilities.db_object_wrappers.thermodynamic_properties.ThermodynamicProperties(hessian_property_id, energy_property_id, property_collection, structure_collection, structure_id)[source]

A class that provides easy access to free energy and electronic energy contributions of a structure.

Parameters:
hessian_property_id

The ID of the hessian property.

energy_property_id

The ID of the electronic energy property.

property_collection

The property collection.

structure_collection

The structure collection.

get_electronic_energy()[source]

Getter for the electronic energy.

Return type:

float

get_enthalpy(reference_state)[source]

Getter for the enthalpy (may update the reference state).

Return type:

float

Parameters:
reference_state

The thermodynamic reference state (temperature/pressure).

get_entropy(reference_state)[source]

Getter for the entropy (may update the reference state).

Return type:

float

Parameters:
reference_state

The thermodynamic reference state (temperature/pressure).

get_gibbs_free_energy(reference_state)[source]

Getter for the gibbs free energy (may update the reference state).

Return type:

float

Parameters:
reference_state

The thermodynamic reference state (temperature/pressure).

get_gibbs_free_energy_correction(reference_state)[source]

Getter for the gibbs free energy correction (may update the reference state).

Return type:

float

Parameters:
reference_state

The thermodynamic reference state (temperature/pressure).

get_reference_state()[source]

Getter for the last thermodynamic reference state.

Return type:

ReferenceState

get_thermochemistry_calculator()[source]

Getter for the ThermochemistryCalculator.

Return type:

ThermochemistryCalculator

get_zero_point_energy()[source]

Getter for the vibrational zero point energy.

Return type:

float

get_zero_point_energy_correction()[source]

Getter for the vibrational zero point energy correction.

Return type:

float

id()[source]

Getter for the structure’s database ID.

Return type:

ID

class scine_chemoton.utilities.db_object_wrappers.thermodynamic_properties.ThermodynamicPropertiesCache(structures, properties, electronic_model, hessian_model, only_electronic=False)[source]

The purpose of this class is to provide access to all ThermodynamicProperties of an aggregate or set of transition states which provide access to the thermodynamic properties of some structure. Furthermore, it provides access to the minimum Gibbs free energy, enthalpy, entropy etc. since it effectively represents a structure ensemble.

Parameters:
structures

The structure collection.

properties

The property collection

electronic_model

The model for the electronic energies.

hessian_model

The model for the hessian correction. May be None.

only_electronic

If true, only the electronic energies are used to characterize the structures.

TODO: Add conformational free energy contributions after flask deduplication.
clear()[source]

Clear and reset cache.

get_ensemble_enthalpy(reference_state)[source]

Getter for the enthalpy of the structure with the lowest Gibb’s free energy approximation (no conformational contribution).

Return type:

Optional[float]

Returns:
The enthalpy of the structure with the lowest Gibb’s free energy approximation. May return None, if no
enthalpy is available.
get_ensemble_entropy(reference_state)[source]

Getter for the entropy of the structure with the lowest Gibb’s free energy approximation (no conformational contribution).

Return type:

Optional[float]

Returns:
The entropy of the structure with the lowest Gibb’s free energy approximation. May return None, if no
entropy is available.
get_ensemble_gibbs_free_energy(reference_state)[source]

Getter for the Gibb’s free energy of the underlying structure ensemble.

Return type:

Optional[float]

Returns:
At the moment only the minimum value is returned. May return None, if no free energy is available for any
structure.
get_n_cached()[source]

Getter for the number of cached elements.

Return type:

int

get_or_produce(structure_id)[source]

Getter for the thermodynamic properties of a structure.

Return type:

Optional[ThermodynamicProperties]

Parameters:
structure_id

The structure ID.

Returns:
The corresponding ThermodynamicProperties object. May return None, if the structure lacks some free energy
contribution.
get_sorted_structure_list(reference_state)[source]
Return type:

List[Tuple[ID, float]]

minimum_values_need_update(reference_state, n_structures=None)[source]

Check if the cached minimum values must be updated.

Return type:

bool

Parameters:
reference_state

The reference state (temperature/pressure).

n_structures

(optional) The number of elements the cache should contain.

Returns:
Returns true, if the cache is out of date or incomplete.
provides_values()[source]

Returns True if values other than None can be provided by this cache.

Return type:

bool

class scine_chemoton.utilities.db_object_wrappers.wrapper_caches.MultiModelAggregateCache(manager, model_combinations, only_electronic=False)[source]

Multi model cache for aggregate wrappers.

Parameters:
managerManager

The database manager.

model_combinationsList[ModelCombination]

The model combinations.

only_electronicbool

If true, free energies are approximated by their electronic energy contribution only.

get_cache(model_combination)

Getter for a cache corresponding to the given model combination.

Return type:

Union[AggregateCache, ReactionCache, None]

Parameters:
model_combinationModelCombination

The model combination.

Returns:
Optional[Union[AggregateCache, ReactionCache]]

The cache.

get_caches()

Getter for all caches.

Return type:

List[Union[AggregateCache, ReactionCache]]

get_model_combinations()

Getter for all model combinations associated to the caches.

get_or_produce(db_id)

Get an aggregate or reaction from the caches.

Return type:

Union[Aggregate, Reaction]

has_cache(model_combination)

Checks if there is a cache with the given model combination.

Return type:

bool

Parameters:
model_combinationModelCombination

The model combination.

Returns:
Returns true if there is a cache with the given model combination, otherwise false.
class scine_chemoton.utilities.db_object_wrappers.wrapper_caches.MultiModelCacheFactory[source]

Singleton to create multi model aggregate/reaction caches.

clear()[source]

Clear caches.

Return type:

None

get_aggregates_cache(only_electronic, model_combination, manager)[source]

Get an aggregate cache corresponding to the given parameters

Return type:

AggregateCache

Parameters:
only_electronicbool

If true, only the electronic energy is used to approximate free energies.

model_combinationModelCombination

The model combination with which aggregate wrappers will be constructed.

managerdb.Manager

The database manager.

get_reaction_cache(only_electronic, model_combination, manager)[source]

Get a reaction cache corresponding to the given parameters

Return type:

ReactionCache

Parameters:
only_electronicbool

If true, only the electronic energy is used to approximate free energies.

model_combinationModelCombination

The model combination with which reaction wrappers will be constructed.

managerdb.Manager

The database manager.

class scine_chemoton.utilities.db_object_wrappers.wrapper_caches.MultiModelReactionCache(manager, model_combinations, only_electronic=False)[source]

Multi model cache for reaction wrappers.

Parameters:
managerManager

The database manager.

model_combinationsList[ModelCombination]

The model combinations.

only_electronicbool

If true, free energies are approximated by their electronic energy contribution only.

get_cache(model_combination)

Getter for a cache corresponding to the given model combination.

Return type:

Union[AggregateCache, ReactionCache, None]

Parameters:
model_combinationModelCombination

The model combination.

Returns:
Optional[Union[AggregateCache, ReactionCache]]

The cache.

get_caches()

Getter for all caches.

Return type:

List[Union[AggregateCache, ReactionCache]]

get_model_combinations()

Getter for all model combinations associated to the caches.

get_or_produce(db_id)

Get an aggregate or reaction from the caches.

Return type:

Union[Aggregate, Reaction]

has_cache(model_combination)

Checks if there is a cache with the given model combination.

Return type:

bool

Parameters:
model_combinationModelCombination

The model combination.

Returns:
Returns true if there is a cache with the given model combination, otherwise false.
class scine_chemoton.utilities.db_object_wrappers.wrapper_caches.MultipleModelCache[source]

Base class for hierarchical ordered models, aggregate and reaction sets, i.e., the model order is relevant. It is assumed that the model importance/accuracy decreases with its index in the given model combination list.

get_cache(model_combination)[source]

Getter for a cache corresponding to the given model combination.

Return type:

Union[AggregateCache, ReactionCache, None]

Parameters:
model_combinationModelCombination

The model combination.

Returns:
Optional[Union[AggregateCache, ReactionCache]]

The cache.

get_caches()[source]

Getter for all caches.

Return type:

List[Union[AggregateCache, ReactionCache]]

get_model_combinations()[source]

Getter for all model combinations associated to the caches.

get_or_produce(db_id)[source]

Get an aggregate or reaction from the caches.

Return type:

Union[Aggregate, Reaction]

has_cache(model_combination)[source]

Checks if there is a cache with the given model combination.

Return type:

bool

Parameters:
model_combinationModelCombination

The model combination.

Returns:
Returns true if there is a cache with the given model combination, otherwise false.