On-the-fly training

ASE_OTF is the on-the-fly training module for ASE, WITHOUT molecular dynamics engine. It needs to be used adjointly with ASE MD engine.

class flare.ase.otf.ASE_OTF(atoms, timestep, number_of_steps, dft_calc, md_engine, md_kwargs, calculator=None, trajectory=None, **otf_kwargs)

On-the-fly training module using ASE MD engine, a subclass of OTF.

Parameters:
  • atoms (ASE Atoms) – the ASE Atoms object for the on-the-fly MD run.
  • calculator – ASE calculator. Must have “get_uncertainties” method implemented.
  • timestep – the timestep in MD. Please use ASE units, e.g. if the timestep is 1 fs, then set timestep = 1 * units.fs
  • number_of_steps (int) – the total number of steps for MD.
  • dft_calc (ASE Calculator) – any ASE calculator is supported, e.g. Espresso, VASP etc.
  • md_engine (str) – the name of MD thermostat, only VelocityVerlet, NVTBerendsen, NPTBerendsen, NPT and Langevin, NoseHoover are supported.
  • md_kwargs (dict) – Specify the args for MD as a dictionary, the args are as required by the ASE MD modules consistent with the md_engine.
  • trajectory (ASE Trajectory) – default None, not recommended, currently in experiment.

The following arguments are for on-the-fly training, the user can also refer to flare.otf.OTF

Parameters:
  • prev_pos_init ([type], optional) – Previous positions. Defaults to None.
  • rescale_steps (List[int], optional) – List of frames for which the velocities of the atoms are rescaled. Defaults to [].
  • rescale_temps (List[int], optional) – List of rescaled temperatures. Defaults to [].
  • calculate_energy (bool, optional) – If True, the energy of each frame is calculated with the GP. Defaults to False.
  • write_model (int, optional) – If 0, write never. If 1, write at end of run. If 2, write after each training and end of run. If 3, write after each time atoms are added and end of run. If 4, write after each training and end of run, and back up after each write.
  • std_tolerance_factor (float, optional) – Threshold that determines when DFT is called. Specifies a multiple of the current noise hyperparameter. If the epistemic uncertainty on a force component exceeds this value, DFT is called. Defaults to 1.
  • skip (int, optional) – Number of frames that are skipped when dumping to the output file. Defaults to 0.
  • init_atoms (List[int], optional) – List of atoms from the input structure whose local environments and force components are used to train the initial GP model. If None is specified, all atoms are used to train the initial GP. Defaults to None.
  • output_name (str, optional) – Name of the output file. Defaults to ‘otf_run’.
  • max_atoms_added (int, optional) – Number of atoms added each time DFT is called. Defaults to 1.
  • freeze_hyps (int, optional) – Specifies the number of times the hyperparameters of the GP are optimized. After this many updates to the GP, the hyperparameters are frozen. Defaults to 10.
  • n_cpus (int, optional) – Number of cpus used during training. Defaults to 1.
compute_properties()
Compute energies, forces, stresses, and their uncertainties with
the FLARE ASE calcuator, and write the results to the OTF structure object.
md_step()

Get new position in molecular dynamics based on the forces predicted by FLARE_Calculator or DFT calculator

rescale_temperature(new_pos)

Change the previous positions to update the temperature

Parameters:new_pos (np.ndarray) – Positions of atoms in the next MD frame.
update_gp(train_atoms, dft_frcs, dft_energy=None, dft_stress=None)

Updates the current GP model.

Parameters:
  • train_atoms (List[int]) – List of atoms whose local environments will be added to the training set.
  • dft_frcs (np.ndarray) – DFT forces on all atoms in the structure.
update_temperature()

Updates the instantaneous temperatures of the system.

Parameters:new_pos (np.ndarray) – Positions of atoms in the next MD frame.