wfl.calculators.orca package#
Submodules#
wfl.calculators.orca.basinhopping module#
- class wfl.calculators.orca.basinhopping.BasinHoppingORCA(atoms=None, n_hop=10, n_run=3, scratchdir=None, n_missing_tolerance=1, energy_tol=0.001, forces_tol=0.05, seed='orca', n_orb=10, max_angle=60.0, smearing=5000.0, maxiter=500, chained_hops=True, orcasimpleinput='UHF PBE def2-SVP tightscf', orcablocks=None, orca_command='orca', rng=None, keep_files=False, uhf=True, **kwargs)#
Bases:
Calculator
ORCA calculator with basin hopping in wavefunction space for smooth PES of radicals
Call n_runs (3 <= recommended) instances of calculation on each frame. If all agree in energy within a given margin, then it is a successful.
calculations:
no rotation, only smearing calculation
n_hop - 1 number of calculations with random rotations
- Parameters
atoms (ase.Atoms) – molecule for calculation
n_hop (int) – number of hopping steps
n_run (int) – number of independent runs to perform, resultant energy/force is checked to be the same
scratchdir (path_like) – put temporary directories here for the calculations, one per task very efficient if you have an SSD scratch disk on nodes that perform the tasks
n_missing_tolerance (int) – number of files allowed to be missing
energy_tol (float) – energy tolerance in eV/atom
force_tol (float) – force tolerance in eV/Angstrom per component
seed (str) – file name seed
n_orb (int) – number of orbital pairs to rotate
max_angle (float) – maximal angle to rotate orbitals to
smearing (float, default 5000.) – smearing temperature in K, recommended is 5000 K
maxiter (int, default 500) – maximum number of SCF iterations to do
chained_hops (bool, default True) – chain hops, ie. to take the previous wavefunction result or the very initial one
orcasimpleinput (str) – What you’d put after the “!” in an orca input file.
orcablocks (str) – What you’d put in the “% … end”-blocks.
orca_command (str) – command of ORCA calculator as in path
rng (np.random.Generator) – random number generator
keep_files (bool) – to keep the resultant files from each calculation
uhf (bool, default True) – to use Unrestricted HF/KS method. This is advised to be used at all times, implemented for development purposes.
kwargs –
- calculate(atoms=None, properties=None, system_changes=['positions', 'numbers', 'cell', 'pbc', 'initial_charges', 'initial_magmoms'], keep_raw_results=False, verbose=False)#
Do the calculation.
- properties: list of str
List of what needs to be calculated. Can be any combination of ‘energy’, ‘forces’, ‘stress’, ‘dipole’, ‘charges’, ‘magmom’ and ‘magmoms’.
- system_changes: list of str
List of what has changed since last calculation. Can be any combination of these six: ‘positions’, ‘numbers’, ‘cell’, ‘pbc’, ‘initial_charges’ and ‘initial_magmoms’.
Subclasses need to implement this, but can ignore properties and system_changes if they want. Calculated properties should be inserted into results dictionary like shown in this dummy example:
self.results = {'energy': 0.0, 'forces': np.zeros((len(atoms), 3)), 'stress': np.zeros(6), 'dipole': np.zeros(3), 'charges': np.zeros(len(atoms)), 'magmom': 0.0, 'magmoms': np.zeros(len(atoms))}
The subclass implementation should first call this implementation to set the atoms attribute and create any missing directories.
- get_homo()#
Get the index of the HOMO.
ORCA counts form 0 and fills the alpha (indexed 0) spin more than the beta (indexed 1) for odd n_elec.
- Returns
i_homo_alpha, i_homo_beta
- Return type
int
- get_multiplicity()#
Gets the multiplicity of the atoms object.
ORCA format: singlet = 1, doublet = 2, etc.
- Returns
multiplicity
- Return type
int
- implemented_properties: List[str] = ['energy', 'forces']#
Properties calculator can handle (energy, forces, …)
- process_results(energy_array, force_array)#
Compares the resultant energies and forces
Maximal difference between energy minima per atom and corresponding force components are checked against the set thresholds.
- Parameters
energy_array (array, shape(n_runs, n_hop)) –
force_array (array, shape(n_runs, n_hop, len_atoms, 3)) –
- Returns
energy
forces
- wfl.calculators.orca.basinhopping.evaluate_basin_hopping(*args, **kwargs)#
Evaluate with BasinHoppingORCA calculator
- Parameters
inputs (list(Atoms) / ConfigSet) – input atomic configs, needs to be iterable
outputs (OutputSpec) – output atomic configs
workdir_root (path-like, default os.getcwd()) – directory to put calculation directories into
rundir_prefix (str, default 'ORCA_') – directory name prefix for calculations
keep_files ("default" / bool) –
what kind of files to keep from the run
”default : .out, .inp, .ase, .engrad is kept – ase can read the results again
True : all files kept
False: none kept
orca_kwargs (dict) – kwargs for BasinHoppingORCA calculator
output_prefix (str, default None) – prefix for keys in the
- Returns
results – ConfigSet(outputs)
- Return type
Module contents#
- class wfl.calculators.orca.ORCA(keep_files='default', rundir_prefix='ORCA_', scratchdir=None, workdir=None, calculator_exec=None, post_process=None, **kwargs)#
Bases:
WFLFileIOCalculator
,ORCA
Extension of ASE’s ORCA calculator that can be used by wfl.calculators.generic
Notes
“directory” argument cannot be present. Use rundir_prefix and workdir instead.
Unless specified, multiplicity is set to singlet/doublet for closed-/open-shell structures
Results from additional optimisation task (set to “opt” or “copt”) are stored in extra_results dictionary.
- Parameters
keep_files (bool / None / "default" / list(str), default "default") –
- what kind of files to keep from the run
True : everything kept
None, False : nothing kept, unless calculation fails
”default” : only ones needed for NOMAD uploads (‘*.pwo’)
list(str) : list of file globs to save
rundir_prefix (str / Path, default 'ORCA_') – Run directory name prefix
workdir (str / Path, default . at calculate time) – Path in which rundir will be created.
scratchdir (str / Path, default None) – temporary directory to execute calculations in and delete or copy back results (set by “keep_files”) if needed. For example, directory on a local disk with fast file I/O.
calculator_exec (str, default "orca") – command for ORCA, without any prefix or redirection set. for example: “/path/to/orca” mutually exclusive with “command”
post_process (function that takes the current instance of the calculator and is to be) – executed after reading back results from file, but before all the files are deleted. For example, a localisation scheme that uses the wavefunction files and saves local charges to ORCA.extra_results.
**kwargs (arguments for ase.calculators.espresso.Espresso) –
- calculate(atoms=None, properties=['energy', 'forces'], system_changes=['positions', 'numbers', 'cell', 'pbc', 'initial_charges', 'initial_magmoms'])#
Does the calculation. Handles the working directories in addition to regular ASE calculation operations (writing input, executing, reading_results)
- cleanup()#
Clean all (empty) directories that could not have been removed immediately after the calculation, for example, because other parallel process might be using them.
- default_parameters = {'charge': 0, 'mult': None, 'orcablocks': '%scf maxiter 200 end', 'orcasimpleinput': 'tightscf PBE def2-SVP', 'task': 'engrad'}#
- static get_default_multiplicity(atoms, charge=0)#
Gets the multiplicity of the atoms object.
ORCA format: singlet = 1, doublet = 2, etc.
- Parameters
atoms (ase.Atoms) –
charge (int) –
- Returns
multiplicity
- Return type
int
- implemented_properties: List[str] = ['energy', 'forces', 'dipole']#
Properties calculator can handle (energy, forces, …)
- is_converged()#
checks for warnings about SCF/wavefunction not converging.
- Returns
None
if “FINAL SINGLE POINT ENERGY” is not in the outputFalse
if “Wavefunction not fully converged” is in the fileTrue
otherwise
Based on ase.calculators.orca.read_energy().
- pick_task()#
- read_dipole()#
Note
Dipole is calculated by default, though only written up to 0.00001 so the rest of the digits are not meaningful in the output.
- read_frequencies()#
Reads frequencies (dynamical matrix eigenvalues) and normal modes (eigenvectors) from output. Currently broken.
- read_opt_atoms()#
Reads the result of the geometry optimisation
- read_results()#
Reads all results
- read_trajectory()#
Reads the trajectory of the geometry optimisation
Notes
Each comment line in output xyz has “Coordinates from ORCA-job orca E -154.812399026326”; # TODO parse out forces as well
- wfl_generic_default_autopara_info = {'num_inputs_per_python_subprocess': 1}#
- write_input(atoms, properties=None, system_changes=None)#
Writes orca.inp, based on the wfl ORCA calculator parameters
- wfl.calculators.orca.natural_population_analysis(janpa_home_dir, orca_calc)#
- wfl.calculators.orca.parse_npa_output(fname)#