Esempio n. 1
0
    def __iter__(self):
        if type(self.filename) == type(''):
            f = open(self.filename)
            opened = True

        # First eight lines are header
        header = [f.readline() for i in range(8)]

        # Rest of file is data array
        data = np.loadtxt(f)

        lattice = fzeros((3, 3))
        for i in [1, 2, 3]:
            lattice[:, i] = [float(x) for x in header[i + 1].split()[1:]]
        if self.vacuum is not None:
            lattice += np.diag(self.vacuum)
        at = Atoms(n=len(data), lattice=lattice)
        at.pos[...] = data[:, 3:6].T
        at.set_atoms(self.z)
        at.add_property('tag', data[:, 1].astype(int))
        at.add_property('mass', data[:, 2])
        at.add_property('velo', data[:, 6:9].T)
        at.add_property('epot', data[:, 9])
        if self.fix_tags is not None:
            at.add_property('move_mask', 1)
            for tag in self.fix_tags:
                at.move_mask[at.tag == tag] = 0
        if opened:
            f.close()

        yield at
Esempio n. 2
0
def dict2atoms(row):
    from quippy.elasticity import stress_matrix

    at = row.toatoms(add_additional_information=True)

    f = None
    try:
        f = at.get_forces()
    except RuntimeError:
        pass

    e = None
    try:
        e = at.get_potential_energy()
    except RuntimeError:
        pass

    v = None
    try:
        v = -stress_matrix(at.get_stress()) * at.get_volume()
    except RuntimeError:
        pass

    at = Atoms(at)  # convert to quippy Atoms
    if f is not None:
        at.add_property('force', f.T)
    if e is not None:
        at.params['energy'] = e
    if v is not None:
        at.params['virial'] = v

    if 'key_value_pairs' in atoms.params:
        atoms.params.update(atoms.params['key_value_pairs'])
        del atoms.params['key_value_pairs']

    if 'data' in atoms.params:
        for (key, value) in atoms.params['data'].items():
            key = str(key)  # avoid unicode strings
            value = np.array(value)
            if value.dtype.kind == 'U':
                value = value.astype(str)
            if value.dtype.kind != 'S':
                value = value.T
            try:
                at.add_property(key, value)
            except (TypeError, RuntimeError):
                at.params[key] = value

    return at
Esempio n. 3
0
def ASEReader(source, format=None):
    """
    Helper routine to load from ASE trajectories
    """
    from ase.io import read
    from ase.atoms import Atoms as AseAtoms
    from quippy.elasticity import stress_matrix

    format_converter = {
        'POSCAR': 'vasp',
        'CONTCAR': 'vasp',
        'OUTCAR': 'vasp-out',
        'vasprun.xml': 'vasp-xml'
    }
    format = format_converter.get(format, format)
    images = read(source, index=slice(None, None, None), format=format)
    if isinstance(images, AseAtoms):
        images = [images]

    for at in images:
        f = None
        try:
            f = at.get_forces()
        except RuntimeError:
            pass

        e = None
        try:
            e = at.get_potential_energy(force_consistent=True)
        except RuntimeError:
            pass

        v = None
        try:
            v = -stress_matrix(at.get_stress()) * at.get_volume()
        except RuntimeError:
            pass

        at = Atoms(at)  # convert to quippy Atoms
        if f is not None:
            at.add_property('force', f.T)
        if e is not None:
            at.params['energy'] = e
        if v is not None:
            at.params['virial'] = v

        yield at
Esempio n. 4
0
def ASEReader(source):
    """
    Helper routine to load from ASE trajectories
    """
    from ase.io import read
    from ase.atoms import Atoms as AseAtoms
    from quippy.elasticity import stress_matrix

    images = read(source, index=slice(None, None, None))
    if isinstance(images, AseAtoms):
        images = [images]

    for at in images:
        f = None
        try:
            f = at.get_forces()
        except RuntimeError:
            pass

        e = None
        try:
            e = at.get_potential_energy()
        except RuntimeError:
            pass

        v = None
        try:
            v = -stress_matrix(at.get_stress()) * at.get_volume()
        except RuntimeError:
            pass

        at = Atoms(at)  # convert to quippy Atoms
        if f is not None:
            at.add_property('force', f.T)
        if e is not None:
            at.params['energy'] = e
        if v is not None:
            at.params['virial'] = v

        yield at
Esempio n. 5
0
def CP2KOutputReader(fh,
                     module=None,
                     type_map=None,
                     kind_map=None,
                     format=None):

    # mapping from run type to (default module index, list of available module)
    run_types = {
        'QS': ['QUICKSTEP'],
        'QMMM': ['FIST', 'QM/MM', 'QUICKSTEP'],
        'MM': ['FIST']
    }

    filename, lines = read_text_file(fh)
    run_type = cp2k_run_type(cp2k_output=lines)

    if type_map is None:
        type_map = {}
    if kind_map is None:
        kind_map = {}

    try:
        available_modules = run_types[run_type]
    except KeyError:
        raise ValueError('Unknown CP2K run type %s' % run_type)

    if module is None:
        module = available_modules[0]

    try:
        cell_index = available_modules.index(module)
    except ValueError:
        raise ValueError("Don't know how to read module %s from file %s" %
                         (module, filename))

    cell_lines = [
        i for i, line in enumerate(lines) if line.startswith(" CELL| Vector a")
    ]
    if cell_lines == []:
        raise ValueError("Cannot find cell in file %s" % filename)

    try:
        cell_line = cell_lines[cell_index]
    except IndexError:
        raise ValueError(
            "Cannot find cell with index %d in file %s for module %s" %
            (cell_index, filename, module))

    lattice = fzeros((3, 3))
    for i in [0, 1, 2]:
        lattice[:,
                i + 1] = [float(c) for c in lines[cell_line + i].split()[4:7]]

    try:
        start_line = lines.index(
            " MODULE %s:  ATOMIC COORDINATES IN angstrom\n" % module)
    except ValueError:
        raise ValueError(
            "Cannot find atomic positions for module %s in file %s" %
            (module, filename))

    kinds = []
    species = []
    Zs = []
    pos = []
    masses = []
    Zeffs = []
    types = []
    qeffs = []
    for line in lines[start_line + 4:]:
        if line.strip() == '':
            break
        if module == 'FIST':
            atom, kind, typ, x, y, z, qeff, mass = line.split()
            types.append(typ)
            Z = type_map.get(typ, 0)
            kind = int(kind)
            if Z == 0:
                Z = kind_map.get(kind, 0)
            Zs.append(Z)
            qeffs.append(float(qeff))
        else:
            atom, kind, sp, Z, x, y, z, Zeff, mass = line.split()
            species.append(sp)
            Zs.append(int(Z))
            Zeffs.append(float(Zeff))
        kinds.append(int(kind))
        pos.append([float(x), float(y), float(z)])
        masses.append(float(mass))

    at = Atoms(n=len(kinds), lattice=lattice)
    at.pos[...] = farray(pos).T
    at.set_atoms(Zs)
    at.add_property('mass', farray(masses) * MASSCONVERT)
    at.add_property('kind', kinds)
    if module == 'FIST':
        at.add_property('type', ' ' * TABLE_STRING_LENGTH)
        at.add_property('qm', False)
        at.qm[:] = (at.type.stripstrings()
                    == '_QM_') | (at.type.stripstrings() == '_LNK')
        at.type[...] = s2a(types, TABLE_STRING_LENGTH)
        at.add_property('qeff', qeffs)
    else:
        at.species[...] = s2a(species, TABLE_STRING_LENGTH)
        at.add_property('zeff', Zeffs)

    yield at
Esempio n. 6
0
def read_xml_output(xmlfile,energy_from=None, extract_forces=False, extract_dipole=False, datafile=None, cluster=None):
    #parse an xml output file and return cluster with updated info
    # datafile tells which energies, forces to look for, cluster Atoms object which gets returned, this is echoed in the xml file so can be left out
    # If extract_forces is not given and the FORCE keyword is found in datafile, the default is to set extract_forces=True

    log = logging.getLogger('molpro_driver')
    
    if datafile is None:
        datafile=MolproDatafile(xml=xmlfile)
        if 'FORCE' in datafile:
            extract_forces=True

    energy_names = OrderedDict()
    energy_names['CCSD(T)-F12'] = ["total energy"]
    energy_names['CCSD(T)'] = ["total energy"]
    energy_names['MP2'] = ["total energy"]
    energy_names['DF-MP2'] = ["total energy"]
    energy_names['DF-RMP2'] = ["energy"]
    energy_names['RKS'] = ["Energy"]
    energy_names['RHF'] = ["Energy"]
    energy_names['DF-RHF'] = ["Energy"]
    energy_names['HF'] = ["Energy"]
    energy_names['DF-HF'] = ["Energy"]
    #etc
    
    gradient_names = OrderedDict()
    gradient_names['CCSD(T)'] =[""]
    gradient_names['RKS'] =['RKS GRADIENT']
    gradient_names['MP2'] =['MP2 GRADIENT']

    all_methods=OrderedDict()
    all_methods['HF']=["RHF"]
    all_methods['DF-HF']=["RHF"]
    all_methods['RHF']=["RHF"]
    all_methods['DF-RHF']=["RHF"]
    all_methods['MP2']=["MP2"]
    all_methods['DF-MP2']=["MP2"]
    all_methods['DF-RMP2']=["DF-RMP2"]
    all_methods['RKS']=["RKS"]
    all_methods['CCSD(T)-F12']=["CCSD(T)-F12a","CCSD(T)-F12b"]
    all_methods['CCSD(T)']=["CCSD(T)"]

    if energy_from is None:
        log.critical("don't know which energy to extract, use keyword energy_from with options "+str([all_methods[k] for k in iter(all_methods)]).replace('[','').replace(']',''))

    #loop through datafile to look for methods.
    calcs=[] #holds the keys for getting correct method, energy_name, gradient_name
    data_keys_upper = [key.upper() for key in datafile._keys]
    for key in all_methods._keys:
       if key in data_keys_upper:
           calcs.append(key)
    dom = minidom.parse(xmlfile)    

    elements=[]
    position_matrix=[]
    cml = dom.documentElement.getElementsByTagName('cml:atomArray')

    for l in cml[0].childNodes:
        if l.nodeType== 1:
            element=l.attributes['elementType'].value.encode('ascii','ignore')
            elements.append(atomic_number(element))
            posx = l.attributes['x3'].value.encode('ascii','ignore')
            posy = l.attributes['y3'].value.encode('ascii','ignore')
            posz = l.attributes['z3'].value.encode('ascii','ignore')
            position_matrix.append([float(posx),float(posy),float(posz)])
    if cluster is None:
        cluster = Atoms(n=len(elements))
        cluster.set_atoms(elements)
        position_matrix=farray(position_matrix).T
        if not 'ANGSTROM' in datafile._keys and not 'angstrom' in datafile._keys:
            position_matrix = position_matrix * (1.0/0.529177249)
        cluster.pos[:,:]=position_matrix
        #note this leaves the lattice undefined

    #now look for each of these energies in xml file
    energy_found=False
    props = dom.documentElement.getElementsByTagName('property')
    for prop in props:
        prop_name = prop.attributes['name'].value.encode('ascii','ignore')
        prop_method = prop.attributes['method'].value.encode('ascii','ignore')
        for calc in calcs:
            if prop_name in energy_names[calc] and prop_method in all_methods[calc]:
                energy_param_name="_".join([prop_method,prop_name])
                energy_param_name=energy_param_name.replace(" ","_")
                #log.info("found "+energy_param_name)
                # dated routines for finding monomer pairs, triplets in Topology module
                energy_param=prop.attributes['value'].value.encode('ascii','ignore')
                my_energy=energy_param_name
                i_en=1
                while my_energy in cluster.params.iterkeys():
                    i_en+=1
                    my_energy='_'.join([energy_param_name,str(i_en)])
                cluster.params[my_energy] = float(energy_param) * HARTREE
                if prop_method == energy_from:
                    cluster.params['Energy']=float(energy_param) * HARTREE
                    energy_found=True
            elif extract_dipole and prop_name=='Dipole moment':
                dipole_param_name="_".join([prop_method,prop_name])
                dipole_param_name=dipole_param_name.replace(" ","_")
                log.info("found dipole moment: "+dipole_param_name)
                dipole_param=prop.attributes['value'].value.encode('ascii','ignore')
                cluster.params[dipole_param_name]=dipole_param

    if not energy_found:
        log.critical("couldn't find energy from "+energy_from+" prop method : "+prop_method)
                      
        
                
    # read gradients if requested
    if extract_forces:
        if not cluster.has_property('force'):
            cluster.add_property('force', 0.0, n_cols=3)

        grads = dom.documentElement.getElementsByTagName('gradient')
        force_matrix = grads[0].childNodes[0].data.split('\n')
        force_matrix = [str(i).split() for i in force_matrix]
        for i in force_matrix:
            try:
                force_matrix.remove([])
            except ValueError:
                break
        force_matrix = [[(-1.0 * HARTREE / BOHR) * float(j) for j in i]
                        for i in force_matrix]
       
        cluster.force[:] =farray(force_matrix).T

        if len(grads) != 1:
            for k in range(1,len(grads)):
                my_force='force%s'%str(k+1)
                force_matrix = grads[k].childNodes[0].data.split('\n')
                force_matrix = [str(i).split() for i in force_matrix]
                for i in force_matrix:
                    try:
                        force_matrix.remove([])
                    except ValueError:
                        break
                force_matrix = [[(-1.0 * HARTREE / BOHR) * float(j) for j in i]
                                for i in force_matrix]
                cluster.add_property(my_force,farray(force_matrix).T)

    return cluster
Esempio n. 7
0
def CubeReader(f, property_name='charge', discard_repeat=True, format=None):
    def convert_line(line, *fmts):
        return (f(s) for f, s in zip(fmts, line.split()))

    if type(f) == type(''):
        f = open(f)
        opened = True

    # First two lines are comments
    comment1 = f.readline()
    comment2 = f.readline()

    # Now number of atoms and origin
    n_atom, origin_x, origin_y, origin_z = convert_line(
        f.readline(), int, float, float, float)
    origin = farray([origin_x, origin_y, origin_z]) * BOHR

    # Next three lines define number of voxels and shape of each element
    shape = [0, 0, 0]
    voxel = fzeros((3, 3))
    for i in (1, 2, 3):
        shape[i - 1], voxel[1, i], voxel[2, i], voxel[3, i] = convert_line(
            f.readline(), int, float, float, float)

    at = Atoms(n=n_atom, lattice=voxel * BOHR * shape)
    at.add_property(property_name, 0.0)
    prop_array = getattr(at, property_name)

    # Now there's one line per atom
    for i in frange(at.n):
        at.z[i], prop_array[i], at.pos[1,
                                       i], at.pos[2,
                                                  i], at.pos[3,
                                                             i] = convert_line(
                                                                 f.readline(),
                                                                 int, float,
                                                                 float, float,
                                                                 float)
        at.pos[:, i] *= BOHR
    at.set_atoms(at.z)

    # Rest of file is volumetric data
    data = np.fromiter((float(x) for x in f.read().split()), float, count=-1)
    if data.size != shape[0] * shape[1] * shape[2]:
        raise IOError("Bad array length - expected shape %r, but got size %d" %
                      (shape, data.size))

    # Save volumetric data in at.data
    data = farray(data.reshape(shape))

    # Discard periodic repeats?
    if discard_repeat:
        at.data = data[:-1, :-1, :-1]
        shape = [s - 1 for s in shape]
        at.set_lattice(voxel * BOHR * shape, False)

    at.params['comment1'] = comment1
    at.params['comment2'] = comment2
    at.params['origin'] = origin
    at.params['shape'] = shape

    # save grids in at.grid_x, at.grid_y, at.grid_z
    if at.is_orthorhombic:
        at.grid_x, at.grid_y, at.grid_z = np.mgrid[origin[1]:origin[1] +
                                                   at.lattice[1, 1]:shape[0] *
                                                   1j, origin[2]:origin[2] +
                                                   at.lattice[2, 2]:shape[1] *
                                                   1j, origin[3]:origin[3] +
                                                   at.lattice[3, 3]:shape[2] *
                                                   1j]

    if opened:
        f.close()

    yield at
Esempio n. 8
0
                                                lotf_spring_hops=2,
                                                test_mode=params.test_mode,
                                                save_clusters=params.save_clusters,
                                                force_restart=params.force_restart)
                                                
     atoms.set_calculator(qmmm_pot)
 system_timer('init_fm_pot')
 
 if not params.continuation and not params.classical:
     print 'Finding initial Quantum Core positions...'
     if geom =='disloc':
       atoms  = set_quantum_disloc(atoms)
     elif geom=='crack':
       mom   = [3.0 for at in range(len(atoms))]
       atoms.set_initial_magnetic_moments(mom)
       atoms.add_property('hybrid', 0, overwrite=True)
       atoms.add_property('hybrid_vec', 0, overwrite=True)
       atoms.add_property('hybrid_1', 0)
       atoms.add_property('hybrid_mark_1', 0)
       crackpos = atoms.info['CrackPos']
       qm_list_old = []
       qm_list   = update_hysteretic_qm_region(atoms, [], crackpos, params.qm_inner_radius,
                                               params.qm_outer_radius,
                                               update_marks=False)
       atoms.hybrid[qm_list] = HYBRID_ACTIVE_MARK
       atoms.hybrid_vec[qm_list] = HYBRID_ACTIVE_MARK
       atoms.hybrid_1[qm_list] = HYBRID_ACTIVE_MARK
       atoms.hybrid_mark_1[qm_list] = HYBRID_ACTIVE_MARK
       atoms.params['core'] = crackpos
       print HYBRID_ACTIVE_MARK
       print 'Core Found. No. Quantum Atoms:', sum(atoms.hybrid[:])
Esempio n. 9
0
File: asap.py Progetto: xielm12/QUIP
def PosCelReader(basename=None,
                 pos='pos.in',
                 cel='cel.in',
                 force='force.in',
                 energy='energy.in',
                 stress='stress.in',
                 species_map={
                     'O': 1,
                     'Si': 2
                 },
                 cel_angstrom=False,
                 pos_angstrom=False,
                 rydberg=True,
                 format=None):

    if basename is not None:
        basename = os.path.splitext(basename)[0]
        pos = '%s.pos' % basename
        cel = '%s.cel' % basename
        energy = '%s.ene' % basename
        stress = '%s.str' % basename
        force = '%s.for' % basename

    doenergy = os.path.exists(energy)
    doforce = os.path.exists(force)
    dostress = os.path.exists(stress)

    if isinstance(pos, str): pos = open(pos)
    if isinstance(cel, str): cel = open(cel)
    if doenergy and isinstance(energy, str): energy = open(energy)
    if doforce and isinstance(force, str): force = open(force)
    if dostress and isinstance(stress, str): stress = open(stress)

    pos = iter(pos)
    cel = iter(cel)
    if doenergy: energy = iter(energy)
    if doforce: force = iter(force)
    if dostress: stress = iter(stress)

    pos.next()  # throw away blank line at start
    if doforce: force.next()

    rev_species_map = dict(zip(species_map.values(), species_map.keys()))

    while True:

        poslines = list(
            itertools.takewhile(
                lambda L: L.strip() != '' and not L.strip().startswith('STEP'),
                pos))
        if poslines == []:
            break

        cellines = list(itertools.islice(cel, 4))
        #lattice = farray([ [float(x) for x in L.split()] for L in cellines[1:4] ]).T
        lattice = fzeros((3, 3))
        for i in (1, 2, 3):
            lattice[:, i] = [float(x) for x in cellines[i].split()]
        if not cel_angstrom: lattice *= BOHR

        at = Atoms(n=len(poslines), lattice=lattice)
        at.pos[:] = farray([[float(x) for x in L.split()[0:3]]
                            for L in poslines]).T
        if not pos_angstrom: at.pos[:] *= BOHR
        species = [rev_species_map[int(L.split()[3])] for L in poslines]
        elements = [
            not el.isdigit() and atomic_number(el) or el for el in species
        ]
        at.set_atoms(elements)

        if doenergy:
            at.params['energy'] = float(energy.next().split()[0])
            if rydberg:
                at.params['energy'] *= RYDBERG

        if dostress:
            stress_lines = list(itertools.islice(stress, 4))
            virial = farray([[float(x) for x in L.split()]
                             for L in stress_lines[1:4]])
            virial *= at.cell_volume() / (10.0 * GPA)
            at.params['virial'] = virial

        if doforce:
            at.add_property('force', 0.0, n_cols=3)
            force_lines = list(
                itertools.takewhile(lambda L: L.strip() != '', force))
            if len(force_lines) != at.n:
                raise ValueError("len(force_lines) (%d) != at.n (%d)" %
                                 (len(force_lines), at.n))
            at.force[:] = farray([[float(x) for x in L.split()[0:3]]
                                  for L in force_lines]).T
            if rydberg:
                at.force[:] *= RYDBERG / BOHR

        yield at
Esempio n. 10
0
class Potential(_potential.Potential):
    __doc__ = update_doc_string(
        _potential.Potential.__doc__,
        r"""
The :class:`Potential` class also implements the ASE
:class:`ase.calculators.interface.Calculator` interface via the
the :meth:`get_forces`, :meth:`get_stress`, :meth:`get_stresses`,
:meth:`get_potential_energy`, :meth:`get_potential_energies`
methods. This simplifies calculation since there is no need
to set the cutoff or to call :meth:`~quippy.atoms.Atoms.calc_connect`,
as this is done internally. The example above reduces to::

    atoms = diamond(5.44, 14)
    atoms.rattle(0.01)
    atoms.set_calculator(pot)
    forces = atoms.get_forces()
    print forces

Note that the ASE force array is the transpose of the QUIP force
array, so has shape (len(atoms), 3) rather than (3, len(atoms)).

The optional arguments `pot1`, `pot2` and `bulk_scale` are
used by ``Sum`` and ``ForceMixing`` potentials (see also
wrapper class :class:`ForceMixingPotential`)

An :class:`quippy.mpi_context.MPI_context` object can be
passed as the `mpi_obj` argument to restrict the
parallelisation of this potential to a subset of the

The `callback` argument is used to implement the calculation of
the :class:`Potential` in a Python function: see :meth:`set_callback` for
an example.

In addition to the builtin QUIP potentials, it is possible to
use any ASE calculator as a QUIP potential by passing it as
the `calculator` argument to the :class:`Potential` constructor, e.g.::

   from ase.calculators.morse import MorsePotential
   pot = Potential(calculator=MorsePotential)

`cutoff_skin` is used to set the :attr:`cutoff_skin` attribute.

`atoms` if given, is used to set the calculator associated
with `atoms` to the new :class:`Potential` instance, by calling
:meth:'.Atoms.set_calculator`.

.. note::

    QUIP potentials do not compute stress and per-atom stresses
    directly, but rather the virial tensor which has units of stress
    :math:`\times` volume, i.e. energy. If the total stress is
    requested, it is computed by dividing the virial by the atomic
    volume, obtained by calling :meth:`.Atoms.get_volume`. If per-atom
    stresses are requested, a per-atom volume is needed. By default
    this is taken to be the total volume divided by the number of
    atoms. In some cases, e.g. for systems containing large amounts of
    vacuum, this is not reasonable. The ``vol_per_atom`` calc_arg can
    be used either to give a single per-atom volume, or the name of an
    array in :attr:`.Atoms.arrays` containing volumes for each atom.

""",
        signature=
        'Potential(init_args[, pot1, pot2, param_str, param_filename, bulk_scale, mpi_obj, callback, calculator, cutoff_skin, atoms])'
    )

    callback_map = {}

    def __init__(self,
                 init_args=None,
                 pot1=None,
                 pot2=None,
                 param_str=None,
                 param_filename=None,
                 bulk_scale=None,
                 mpi_obj=None,
                 callback=None,
                 calculator=None,
                 cutoff_skin=1.0,
                 atoms=None,
                 fpointer=None,
                 finalise=True,
                 error=None,
                 **kwargs):

        self.atoms = None
        self._prev_atoms = None
        self.energy = None
        self.energies = None
        self.forces = None
        self.stress = None
        self.stresses = None
        self.elastic_constants = None
        self.unrelaxed_elastic_constants = None
        self.numeric_forces = None
        self._calc_args = {}
        self._default_quantities = []
        self.cutoff_skin = cutoff_skin

        if callback is not None or calculator is not None:
            if init_args is None:
                init_args = 'callbackpot'

        if param_filename is not None:
            param_str = open(param_filename).read()

        if init_args is None and param_str is None:
            raise ValueError('Need one of init_args,param_str,param_filename')

        if init_args is not None:
            if init_args.lower().startswith('callbackpot'):
                if not 'label' in init_args:
                    init_args = init_args + ' label=%d' % id(self)
            else:
                # if param_str missing, try to find default set of QUIP params
                if param_str is None and pot1 is None and pot2 is None:
                    param_str = quip_xml_parameters(init_args)

        if kwargs != {}:
            if init_args is not None:
                init_args = init_args + ' ' + dict_to_args_str(kwargs)
            else:
                init_args = dict_to_args_str(kwargs)

        _potential.Potential.__init__(self,
                                      init_args,
                                      pot1=pot1,
                                      pot2=pot2,
                                      param_str=param_str,
                                      bulk_scale=bulk_scale,
                                      mpi_obj=mpi_obj,
                                      fpointer=fpointer,
                                      finalise=finalise,
                                      error=error)

        if init_args is not None and init_args.lower().startswith(
                'callbackpot'):
            _potential.Potential.set_callback(self, Potential.callback)

            if callback is not None:
                self.set_callback(callback)

            if calculator is not None:
                self.set_callback(calculator_callback_factory(calculator))

        if atoms is not None:
            atoms.set_calculator(self)

    __init__.__doc__ = _potential.Potential.__init__.__doc__

    def calc(self,
             at,
             energy=None,
             force=None,
             virial=None,
             local_energy=None,
             local_virial=None,
             args_str=None,
             error=None,
             **kwargs):

        if not isinstance(args_str, basestring):
            args_str = dict_to_args_str(args_str)

        kw_args_str = dict_to_args_str(kwargs)

        args_str = ' '.join((self.get_calc_args_str(), kw_args_str, args_str))

        if isinstance(energy, basestring):
            args_str = args_str + ' energy=%s' % energy
            energy = None
        if isinstance(energy, bool) and energy:
            args_str = args_str + ' energy'
            energy = None

        if isinstance(force, basestring):
            args_str = args_str + ' force=%s' % force
            force = None
        if isinstance(force, bool) and force:
            args_str = args_str + ' force'
            force = None

        if isinstance(virial, basestring):
            args_str = args_str + ' virial=%s' % virial
            virial = None
        if isinstance(virial, bool) and virial:
            args_str = args_str + ' virial'
            virial = None

        if isinstance(local_energy, basestring):
            args_str = args_str + ' local_energy=%s' % local_energy
            local_energy = None
        if isinstance(local_energy, bool) and local_energy:
            args_str = args_str + ' local_energy'
            local_energy = None

        if isinstance(local_virial, basestring):
            args_str = args_str + ' local_virial=%s' % local_virial
            local_virial = None
        if isinstance(local_virial, bool) and local_virial:
            args_str = args_str + ' local_virial'
            local_virial = None

        potlog.debug(
            'Potential invoking calc() on n=%d atoms with args_str "%s"' %
            (len(at), args_str))
        _potential.Potential.calc(self, at, energy, force, virial,
                                  local_energy, local_virial, args_str, error)

    calc.__doc__ = update_doc_string(
        _potential.Potential.calc.__doc__,
        """In Python, this method is overloaded to set the final args_str to
          :meth:`get_calc_args_str`, followed by any keyword arguments,
          followed by an explicit `args_str` argument if present. This ordering
          ensures arguments explicitly passed to :meth:`calc` will override any
          default arguments.""")

    @staticmethod
    def callback(at_ptr):
        from quippy import Atoms
        at = Atoms(fpointer=at_ptr, finalise=False)
        if 'label' not in at.params or at.params[
                'label'] not in Potential.callback_map:
            raise ValueError('Unknown Callback label %s' % at.params['label'])
        Potential.callback_map[at.params['label']](at)

    def set_callback(self, callback):
        """
        For a :class:`Potential` of type `CallbackPot`, this method is
        used to set the callback function. `callback` should be a Python
        function (or other callable, such as a bound method or class
        instance) which takes a single argument, of type
        :class:`~quippy.atoms.Atoms`. Information about which quantities should be
        computed can be obtained from the `calc_energy`, `calc_local_e`,
        `calc_force`, and `calc_virial` keys in `at.params`. Results
        should be returned either as `at.params` entries (for energy and
        virial) or by adding new atomic properties (for forces and local
        energy).

        Here's an example implementation of a simple callback::

          def example_callback(at):
              if at.calc_energy:
                 at.params['energy'] = ...

              if at.calc_force:
                 at.add_property('force', 0.0, n_cols=3)
                 at.force[:,:] = ...

          p = Potential('CallbackPot')
          p.set_callback(example_callback)
          p.calc(at, energy=True)
          print at.energy
          ...
        """
        Potential.callback_map[str(id(self))] = callback

    def wipe(self):
        """
        Mark all quantities as needing to be recalculated
        """
        self.energy = None
        self.energies = None
        self.forces = None
        self.stress = None
        self.stresses = None
        self.numeric_forces = None
        self.elastic_constants = None
        self.unrelaxed_elastic_constants = None

    def update(self, atoms):
        """
        Set the :class:`~quippy.atoms.Atoms` object associated with this :class:`Potential` to `atoms`.

        Called internally by :meth:`get_potential_energy`,
        :meth:`get_forces`, etc.  Only a weak reference to `atoms` is
        kept, to prevent circular references.  If `atoms` is not a
        :class:`quippy.atoms.Atoms` instance, then a copy is made and a
        warning will be printed.
        """
        # we will do the calculation in place, to minimise number of copies,
        # unless atoms is not a quippy Atoms
        if isinstance(atoms, Atoms):
            self.atoms = weakref.proxy(atoms)
        else:
            potlog.debug(
                'Potential atoms is not quippy.Atoms instance, copy forced!')
            self.atoms = Atoms(atoms)

        # check if atoms has changed since last call
        if self._prev_atoms is not None and self._prev_atoms.equivalent(
                self.atoms):
            return

        # Mark all quantities as needing to be recalculated
        self.wipe()

        # do we need to reinitialise _prev_atoms?
        if self._prev_atoms is None or len(self._prev_atoms) != len(
                self.atoms) or not self.atoms.connect.initialised:
            self._prev_atoms = Atoms()
            self._prev_atoms.copy_without_connect(self.atoms)
            self._prev_atoms.add_property('orig_pos', self.atoms.pos)
        else:
            # _prev_atoms is OK, update it in place
            self._prev_atoms.z[...] = self.atoms.z
            self._prev_atoms.pos[...] = self.atoms.pos
            self._prev_atoms.lattice[...] = self.atoms.lattice

        # do a calc_connect(), setting cutoff_skin so full reconnect will only be done when necessary
        self.atoms.set_cutoff(self.cutoff(), cutoff_skin=self.cutoff_skin)
        potlog.debug(
            'Potential doing calc_connect() with cutoff %f cutoff_skin %r' %
            (self.atoms.cutoff, self.cutoff_skin))
        self.atoms.calc_connect()

    # Synonyms for `update` for compatibility with ASE calculator interface
    def initialize(self, atoms):
        self.update(atoms)

    def set_atoms(self, atoms):
        self.update(atoms)

    def calculation_required(self, atoms, quantities):
        self.update(atoms)
        for quantity in quantities:
            if getattr(self, quantity) is None:
                return True
        return False

    def calculate(self, atoms, quantities=None):
        """
        Perform a calculation of `quantities` for `atoms` using this Potential.

        Automatically determines if a new calculation is required or if previous
        results are still appliciable (i.e. if the atoms haven't moved since last call)
        Called internally by :meth:`get_potential_energy`, :meth:`get_forces`, etc.
        """
        if quantities is None:
            quantities = ['energy', 'forces', 'stress']

        # Add any default quantities
        quantities = set(self.get_default_quantities() + quantities)

        if len(quantities) == 0:
            raise RuntimeError('Nothing to calculate')

        if not self.calculation_required(atoms, quantities):
            return

        args_map = {
            'energy': {
                'energy': None
            },
            'energies': {
                'local_energy': None
            },
            'forces': {
                'force': None
            },
            'stress': {
                'virial': None
            },
            'numeric_forces': {
                'force': 'numeric_force',
                'force_using_fd': True,
                'force_fd_delta': 1.0e-5
            },
            'stresses': {
                'local_virial': None
            },
            'elastic_constants': {},
            'unrelaxed_elastic_constants': {}
        }

        # list of quantities that require a call to Potential.calc()
        calc_quantities = [
            'energy', 'energies', 'forces', 'numeric_forces', 'stress',
            'stresses'
        ]

        # list of other quantities we know how to calculate
        other_quantities = ['elastic_constants', 'unrelaxed_elastic_constants']

        calc_args = {}
        calc_required = False
        for quantity in quantities:
            if quantity in calc_quantities:
                calc_required = True
                calc_args.update(args_map[quantity])
            elif quantity not in other_quantities:
                raise RuntimeError(
                    "Don't know how to calculate quantity '%s'" % quantity)

        if calc_required:
            self.calc(self.atoms, args_str=dict_to_args_str(calc_args))

        if 'energy' in quantities:
            self.energy = float(self.atoms.energy)
        if 'energies' in quantities:
            self.energies = self.atoms.local_energy.view(np.ndarray)
        if 'forces' in quantities:
            self.forces = self.atoms.force.view(np.ndarray).T
        if 'numeric_forces' in quantities:
            self.numeric_forces = self.atoms.numeric_force.view(np.ndarray).T
        if 'stress' in quantities:
            stress = -self.atoms.virial.view(
                np.ndarray) / self.atoms.get_volume()
            # convert to 6-element array in Voigt order
            self.stress = np.array([
                stress[0, 0], stress[1, 1], stress[2, 2], stress[1, 2],
                stress[0, 2], stress[0, 1]
            ])
        if 'stresses' in quantities:
            lv = np.array(self.atoms.local_virial)  # make a copy
            vol_per_atom = self.get('vol_per_atom',
                                    self.atoms.get_volume() / len(atoms))
            if isinstance(vol_per_atom, basestring):
                vol_per_atom = self.atoms.arrays[vol_per_atom]
            self.stresses = -lv.T.reshape(
                (len(atoms), 3, 3), order='F') / vol_per_atom

        if 'elastic_constants' in quantities:
            cij_dx = self.get('cij_dx', 1e-2)
            cij = fzeros((6, 6))
            self.calc_elastic_constants(self.atoms,
                                        fd=cij_dx,
                                        args_str=self.get_calc_args_str(),
                                        c=cij,
                                        relax_initial=False,
                                        return_relaxed=False)
            if not get_fortran_indexing():
                cij = cij.view(np.ndarray)
            self.elastic_constants = cij

        if 'unrelaxed_elastic_constants' in quantities:
            cij_dx = self.get('cij_dx', 1e-2)
            c0ij = fzeros((6, 6))
            self.calc_elastic_constants(self.atoms,
                                        fd=cij_dx,
                                        args_str=self.get_calc_args_str(),
                                        c0=c0ij,
                                        relax_initial=False,
                                        return_relaxed=False)
            if not get_fortran_indexing():
                c0ij = c0ij.view(np.ndarray)
            self.unrelaxed_elastic_constants = c0ij

    def get_potential_energy(self, atoms):
        """
        Return potential energy of `atoms` calculated with this Potential
        """
        self.calculate(atoms, ['energy'])
        return self.energy

    def get_potential_energies(self, atoms):
        """
        Return array of atomic energies calculated with this Potential
        """
        self.calculate(atoms, ['energies'])
        return self.energies.copy()

    def get_forces(self, atoms):
        """
        Return forces on `atoms` calculated with this Potential
        """
        self.calculate(atoms, ['forces'])
        return self.forces.copy()

    def get_numeric_forces(self, atoms):
        """
        Return forces on `atoms` computed with finite differences of the energy
        """
        self.calculate(atoms, ['numeric_forces'])
        return self.numeric_forces.copy()

    def get_stress(self, atoms):
        """
        Return stress tensor for `atoms` computed with this Potential

        Result is a 6-element array in Voigt notation:
           [sigma_xx, sigma_yy, sigma_zz, sigma_yz, sigma_xz, sigma_xy]
        """
        self.calculate(atoms, ['stress'])
        return self.stress.copy()

    def get_stresses(self, atoms):
        """
        Return the per-atoms virial stress tensors for `atoms` computed with this Potential
        """
        self.calculate(atoms, ['stresses'])
        return self.stresses.copy()

    def get_elastic_constants(self, atoms):
        """
        Calculate elastic constants of `atoms` using this Potential.

        Returns  6x6 matrix :math:`C_{ij}` of elastic constants.

        The elastic contants are calculated as finite difference
        derivatives of the virial stress tensor using positive and
        negative strains of magnitude the `cij_dx` entry in
        ``calc_args``.
        """
        self.calculate(atoms, ['elastic_constants'])
        return self.elastic_constants.copy()

    def get_unrelaxed_elastic_constants(self, atoms):
        """
        Calculate unrelaxed elastic constants of `atoms` using this Potential

        Returns 6x6 matrix :math:`C^0_{ij}` of unrelaxed elastic constants.

        The elastic contants are calculated as finite difference
        derivatives of the virial stress tensor using positive and
        negative strains of magnitude the `cij_dx` entry in
        :attr:`calc_args`.
        """
        self.calculate(atoms, ['unrelaxed_elastic_constants'])
        return self.unrelaxed_elastic_constants.copy()

    def get_default_quantities(self):
        "Get the list of quantities to be calculated by default"
        return self._default_quantities[:]

    def set_default_quantities(self, quantities):
        "Set the list of quantities to be calculated by default"
        self._default_quantities = quantities[:]

    def get(self, param, default=None):
        """
        Get the value of a ``calc_args`` parameter for this :class:`Potential`

        Returns ``None`` if `param` is not in the current ``calc_args`` dictionary.

        All calc_args are passed to :meth:`calc` whenever energies,
        forces or stresses need to be re-computed.
        """
        return self._calc_args.get(param, default)

    def set(self, **kwargs):
        """
        Set one or more calc_args parameters for this Potential

        All calc_args are passed to :meth:`calc` whenever energies,
        forces or stresses need to be computed.

        After updating the calc_args, :meth:`set` calls :meth:`wipe`
        to mark all quantities as needing to be recaculated.
        """
        self._calc_args.update(kwargs)
        self.wipe()

    def get_calc_args(self):
        """
        Get the current ``calc_args``
        """
        return self._calc_args.copy()

    def set_calc_args(self, calc_args):
        """
        Set the ``calc_args`` to be used subsequent :meth:`calc` calls
        """
        self._calc_args = calc_args.copy()

    def get_calc_args_str(self):
        """
        Get the ``calc_args`` to be passed to :meth:`calc` as a string
        """
        return dict_to_args_str(self._calc_args)

    def get_cutoff_skin(self):
        return self._cutoff_skin

    def set_cutoff_skin(self, cutoff_skin):
        self._cutoff_skin = cutoff_skin
        self._prev_atoms = None  # force a recalculation

    cutoff_skin = property(get_cutoff_skin,
                           set_cutoff_skin,
                           doc="""
                           The `cutoff_skin` attribute is only relevant when the ASE-style
                           interface to the Potential is used, via the :meth:`get_forces`,
                           :meth:`get_potential_energy` etc. methods. In this case the
                           connectivity of the :class:`~quippy.atoms.Atoms` object for which
                           the calculation is requested is automatically kept up to date by
                           using a neighbour cutoff of :meth:`cutoff` + `cutoff_skin`, and
                           recalculating the neighbour lists whenever the maximum displacement
                           since the last :meth:`Atoms.calc_connect` exceeds `cutoff_skin`.
                           """)