maxwelllink.mxl_drivers.python.models.ase_model module¶
- class maxwelllink.mxl_drivers.python.models.ase_model.ASEModel[source]¶
Bases:
DummyModelGeneral BOMD (Born-Oppenheimer MD) driver using ASE.
This model provides the two key functionalities like other Models supported by MaxwellLink:
MD coupled to E-field: Injects \(F_i^{\\mathrm{ext}} = q_i \\mathbf{E}\) (uniform field) to the molecular forces in MD simulations, where \(q_i\) are per-atom charges (constant user-supplied or calculator-reported each step).
Return source amplitude: \(\\dot{\\boldsymbol{\\mu}} = \\sum_i q_i \\mathbf{v}_i\) (converted to atomic units).
- __init__(atoms, calculator='psi4', calc_kwargs='', charges=None, recompute_charges=False, temperature_K=0.0, verbose=False, checkpoint=False, restart=False, **extra)[source]¶
- Parameters:
atoms (ase.Atoms or str) – Either an ASE
Atomsobject or a path to a structure file (e.g.,.xyz) readable by ASE.calculator (str, default: 'psi4') – Name of ASE calculator (
'psi4','dftb','orca', …).calc_kwargs (str, optional) – String of kwargs passed to the calculator constructor.
charges (str or sequence of float or None, optional) – A string like
"[-1.0 1.0]"representing an array of per-atom charges (in \(\lvert e \rvert\)), separated by space (not comma). IfNone, setrecompute_charges=True.recompute_charges (bool, default: False) – If
True, query charges each step (e.g., Mulliken).temperature_K (float, default: 0.0) – Initial temperature in Kelvin for the Maxwell–Boltzmann distribution.
verbose (bool, default: False) – Whether to print verbose output.
checkpoint (bool, default: False) – Whether to enable checkpointing.
restart (bool, default: False) – Whether to restart from a checkpoint if available.
**extra – Additional keyword arguments forwarded to the calculator constructor.
- append_additional_data()[source]¶
Append additional data to be sent back to MaxwellLink.
The data can be retrieved by the user via the Python interface:
maxwelllink.SocketMolecule.additional_data_history, whereadditional_data_historyis a list of dictionaries.- Returns:
A dictionary containing additional data.
- Return type:
dict
- calc_amp_vector()[source]¶
Return the amplitude vector \(\\mathrm{d}\\boldsymbol{\\mu}/\\mathrm{d}t\) for the current time step in atomic units.
In classical MD: \(\\displaystyle \\frac{\\mathrm{d}\\boldsymbol{\\mu}}{\\mathrm{d}t} = \\sum_i q_i \\mathbf{v}_i\).
- commit_step()¶
Commit the previewed step and return the staged amplitude.
This method applies the changes from the staged step to the internal state and returns the calculated amplitude vector.
Notes
This method should not be overridden by subclasses.
- Returns:
Amplitude vector in the form \([\mathrm{d}\mu_x/\mathrm{d}t,\ \mathrm{d}\mu_y/\mathrm{d}t,\ \mathrm{d}\mu_z/\mathrm{d}t]\).
- Return type:
numpy.ndarray of float, shape (3,)
- have_result()¶
Check if a staged step is ready to be committed.
Notes
This method should not be overridden by subclasses.
- Returns:
Whether a staged step is ready.
- Return type:
bool
- initialize(dt_new, molecule_id)[source]¶
Initialize the model with the new time step and molecule ID.
- Parameters:
dt_new (float) – The new time step in atomic units (a.u.).
molecule_id (int) – The ID of the molecule.
- propagate(effective_efield_vec)[source]¶
Propagate the BO molecular dynamics given the effective electric field vector.
- Parameters:
effective_efield_vec (array-like of float, shape (3,)) – Effective electric field vector in the form
[E_x, E_y, E_z].
- stage_step(E_vec)¶
Stage a propagation step with the given effective electric field vector.
This method performs the propagation and calculates the amplitude vector, but does not commit the changes to the internal state. The result can be committed later using the
self.commit_stepmethod.Notes
This method should not be overridden by subclasses.
- Parameters:
E_vec (array-like of float, shape (3,)) – Effective electric field vector in the form
[E_x, E_y, E_z].
- class maxwelllink.mxl_drivers.python.models.ase_model.ForceAugmenter[source]¶
Bases:
CalculatorASE Calculator wrapper that adds an external uniform E-field force \(F_i^{\mathrm{ext}} = q_i \mathbf{E}\) to each atom \(i\), where \(q_i\) are per-atom charges (either fixed or recomputed each step).
- __init__(base, charges=None, recompute_charges=False, verbose=False)[source]¶
- Parameters:
base (ase.calculators.calculator.Calculator) – An ASE Calculator instance to wrap.
charges (array-like or None, optional) – Per-atom charges in \(\lvert e \rvert\) units. If
None, setrecompute_charges=True.recompute_charges (bool, default: False) – If
True, query charges each step (e.g., Mulliken).verbose (bool, default: False) – Whether to print verbose output.
- static __new__(cls, *args, **kwargs)¶
- Parameters:
args (Any)
kwargs (Any)
- Return type:
Any
- calculate(atoms=None, properties=('energy',), system_changes=all_changes)[source]¶
Calculate the requested properties for the given atomic configuration.
- Parameters:
atoms (ase.Atoms, optional) – An ASE
Atomsobject.properties (tuple, default: ('energy',)) – Properties to be calculated (e.g.,
'energy','forces').system_changes (list) – List of changes in the atomic configuration.
- calculate_external_force(atoms)[source]¶
Calculate the external force on each atom due to the uniform E-field.
- Parameters:
atoms (ase.Atoms) – An ASE
Atomsobject.- Returns:
External forces on each atom in eV/Å (ASE internal force units).
- Return type:
numpy.ndarray of float, shape (N, 3)
- calculation_required(atoms, properties)[source]¶
Determine whether a recalculation is required based on changes in the atomic configuration.
- Parameters:
atoms (ase.Atoms) – An ASE
Atomsobject.properties (tuple) – Properties to be calculated (e.g.,
'energy','forces').
- Returns:
Whether a recalculation is required.
- Return type:
bool
- implemented_properties = ('energy', 'forces')¶