deMon-Nano#

deMon-Nano is a density-functional based tight-binding (DFTB) code using atom centered orbitals. This interface makes it possible to use deMon-Nano as a calculator in ASE. You need Slater-Koster files for the combination of atom types of your system. These can be obtained at dftb.org.

Environment variables#

Set environment variables in your configuration file (what is the directory for the Slater-Koster files and what is the name of the executable):

  • bash:

    $ DEMONNANO_BASIS_PATH="/path/to/basis/"  (an example)
    $ ASE_DEMONNANO_COMMAND="/path/to/bin/deMon.username.x (an example)
    

deMon-Nano Calculator (a FileIOCalculator)#

class ase.calculators.demonnano.DemonNano(**kwargs)[source]#

Calculator interface to the deMon-nano code.

ASE interface to the deMon-nano code.

The deMon-nano code can be obtained from http://demon-nano.ups-tlse.fr/

The ASE_DEMONNANO_COMMAND environment variable must be set to run the executable, in bash it would be set along the lines of export ASE_DEMONNANO_COMMAND=”pathway-to-deMon-binary/deMon.username.x”

Parameters:

labelstr

relative path to the run directory

atomsAtoms object

the atoms object

commandstr

Command to run deMon. If not present, the environment variable ASE_DEMONNANO_COMMAND is used

basis_pathstr

Relative path to the directory containing DFTB-SCC or DFTB-0 parameters If not present, the environment variable DEMONNANO_BASIS_PATH is used

restart_pathstr

Relative path to the deMon restart dir

titlestr

Title in the deMon input file.

forcesbool

If True a force calculation is enforced

print_outstr | list

Options for the printing in deMon

input_argumentsdict

Explicitly given input arguments. The key is the input keyword and the value is either a str, a list of str (will be written on the same line as the keyword), or a list of lists of str (first list is written on the first line, the others on following lines.)

Example: Geometry Optimization with ASE#

import numpy as np

from ase import Atoms
from ase.calculators.demonnano import DemonNano
from ase.io import write
from ase.optimize import BFGS

d = 0.9775
t = np.pi / 180 * 110.51
mol = Atoms(
    'H2O', positions=[(d, 0, 0), (d * np.cos(t), d * np.sin(t), 0), (0, 0, 0)]
)

input_arguments = {'DFTB': 'SCC', 'CHARGE': '0.0', 'PARAM': 'PTYPE=BIO'}

calc = DemonNano(label='rundir/', input_arguments=input_arguments)
mol.calc = calc

# optimize geometry
dyn = BFGS(mol, trajectory='test.traj')
dyn.run(fmax=0.01)

write('h2o_optimized.xyz', mol)