.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "generated_examples/1-advanced/09-symplectic-flashmd.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code. .. rst-class:: sphx-glr-example-title .. _sphx_glr_generated_examples_1-advanced_09-symplectic-flashmd.py: Training a Symplectic FlashMD Model =================================== This tutorial demonstrates how to train a symplectic FlashMD model using the FlashMD framework. Symplectic integrators are designed to preserve the geometric properties of Hamiltonian dynamics, making them particularly suitable for long molecular dynamics simulations. By leveraging symplectic integrators, we can achieve more accurate and stable simulations over extended periods. .. GENERATED FROM PYTHON SOURCE LINES 13-27 .. code-block:: Python import copy import subprocess import ase import ase.build import ase.io import ase.units from ase.calculators.emt import EMT from ase.md import VelocityVerlet from ase.md.langevin import Langevin from ase.md.velocitydistribution import MaxwellBoltzmannDistribution .. GENERATED FROM PYTHON SOURCE LINES 28-35 Dataset Creation ---------------- We will create a dataset of molecular dynamics trajectories using ASE and its built-in EMT potential. The dataset will consist of atomic configurations, forces, and energies obtained from NVE simulations. In reality, you might want to use a more accurate baseline such as ab initio MD or a machine-learned interatomic potential (MLIP). .. GENERATED FROM PYTHON SOURCE LINES 36-60 .. code-block:: Python # We start by creating a simple system (a small box of aluminum). atoms = ase.build.bulk("Al", "fcc", cubic=True) * (2, 2, 2) # We first equilibrate the system at 300K using a Langevin thermostat. MaxwellBoltzmannDistribution(atoms, temperature_K=300) atoms.calc = EMT() dyn = Langevin( atoms, 2 * ase.units.fs, temperature_K=300, friction=1 / (100 * ase.units.fs) ) dyn.run(1000) # 2 ps equilibration (around 10 ps is better in practice) # Then, we run a production simulation in the NVE ensemble. trajectory = [] def store_trajectory(): trajectory.append(copy.deepcopy(atoms)) dyn = VelocityVerlet(atoms, 1 * ase.units.fs) dyn.attach(store_trajectory, interval=1) dyn.run(2000) # 2 ps NVE run .. GENERATED FROM PYTHON SOURCE LINES 61-69 Data Preparation ---------------- Note that the data preparation process is similar to the one in the `04-flashmd.py` example, with one key difference. Instead of storing a phase-space point and its future state after one time step, we store the input to the symplectic fixed-point solver. The input is a midpoint that is mapped to the difference in positions and momenta after one time step. .. GENERATED FROM PYTHON SOURCE LINES 70-101 .. code-block:: Python time_lag = 32 spacing = 200 def get_structure_for_dataset(frame_now, frame_ahead): s = copy.deepcopy(frame_now) s.arrays["delta_positions"] = ( frame_ahead.get_positions() - frame_now.get_positions() ) s.arrays["delta_momenta"] = frame_ahead.get_momenta() - frame_now.get_momenta() s.set_positions(0.5 * (frame_now.get_positions() + frame_ahead.get_positions())) s.set_momenta(0.5 * (frame_now.get_momenta() + frame_ahead.get_momenta())) return s structures_for_dataset = [] for i in range(0, len(trajectory) - time_lag, spacing): frame_now = trajectory[i] frame_ahead = trajectory[i + time_lag] s = get_structure_for_dataset(frame_now, frame_ahead) structures_for_dataset.append(s) frame_now_trev = copy.deepcopy(frame_now) frame_ahead_trev = copy.deepcopy(frame_ahead) frame_now_trev.set_momenta(-frame_now_trev.get_momenta()) frame_ahead_trev.set_momenta(-frame_ahead_trev.get_momenta()) s = get_structure_for_dataset(frame_ahead_trev, frame_now_trev) structures_for_dataset.append(s) ase.io.write("midpoint-to-delta.xyz", structures_for_dataset) .. GENERATED FROM PYTHON SOURCE LINES 102-111 Model Training -------------- We can now train a symplectic FlashMD model using the prepared dataset. For example, you can use the following options file: .. literalinclude:: options-flashmd-symplectic.yaml :language: yaml .. GENERATED FROM PYTHON SOURCE LINES 112-114 .. code-block:: Python subprocess.run(["mtt", "train", "options-flashmd-symplectic.yaml"], check=True) .. rst-class:: sphx-glr-script-out .. code-block:: none CompletedProcess(args=['mtt', 'train', 'options-flashmd-symplectic.yaml'], returncode=0) .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 34.150 seconds) .. _sphx_glr_download_generated_examples_1-advanced_09-symplectic-flashmd.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: 09-symplectic-flashmd.ipynb <09-symplectic-flashmd.ipynb>` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: 09-symplectic-flashmd.py <09-symplectic-flashmd.py>` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: 09-symplectic-flashmd.zip <09-symplectic-flashmd.zip>` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_