coulomb_kmc.kmc_fmm module

The KMCFMM class is the user facing interface to this implementation of the FMM-KMC algorithm. The underlying algorithm of this implementation is described in Fast electrostatic solvers for kinetic Monte Carlo simulations.

class coulomb_kmc.kmc_fmm.KMCFMM(positions, charges, domain, N=None, boundary_condition='pbc', r=None, shell_width=0.0, energy_unit=1.0, _debug=False, l=None, max_move=None, cuda_direct=False, mirror_direction=None)

Bases: coulomb_kmc.kmc_fmm._PY_KMCFMM, coulomb_kmc.kmc_inject_extract.InjectorExtractor

This class provides the methods required for the electrostatic energy component in a Kinetic Monte Carlo simulation. It is applicable to systems that are cubic with free space, fully periodic and (plate like) Dirichlet boundary conditions (work in progress).

A KMCFMM instance should be initialised, which performs the initial FMM solve, by calling the initialise method. Furthermore a KMCFMM instance should be freed by calling the free method.

After initialise is called the current system energy is given by the energy property. Moves can be proposed with the propose or propose_with_dats methods. A move can be accepted by calling the accept method. When a move is accepted with the accept method the move is also enacted in the position PositionDat that the KMCFMM instance is initialised with.

Parameters
  • positions (PositionDat) – PositionDat(ncomp=3, dtype=ctypes.c_double) containing charge positions.

  • charges (ParticleDat) – ParticleDat(ncomp=1, dtype=ctypes.c_int64) containing charge values.

  • N (int) – Number of charges.

  • boundary_condition (str) – Specific boundary condition as one of: ‘pbc’ for periodic boundary conditions, ‘free_space’ for the simulation cell in a vacuum or ‘27’ for the simulation cell and the nearest 26 periodic images in a vacuum.

  • r (int) – Number of FMM levels (overrides heuristic choice using N).

  • shell_width (float) – Additional argument for FMM solver, can be left at the default of 0.0.

  • energy_unit (float) – Energy unit to use, defaults to 1.0.

  • _debug (bool) – Debugging flag for FMM solver, can be ignored.

  • l (int) – Number of expansion terms to use for FMM solve and propose/accept operations.

  • max_move (float) – Maximum distance of a proposed move.

  • cuda_direct (bool) – Enable CUDA offloading of direct interactions (work in progress).

  • mirror_direction (tuple) – Tuple of bools, e.g. (True, False, False) to enable mirror charge handling for Dirichlet boundary conditions (work in progress).

accept(move=None)

Accept a move passed as a tuple (int: px, np.array: new_pos), where px is the current particle local id, and new_pos is the new position e.g. (46, np.array(1.0, 2.0, 3.0))

Note this will move the particle in the position dat used in the intialiser.

If the KMC instance is defined with a mirror direction move should be a tuple of moves. (42, np.array(1.0, 2.0, 3.0)), (84, np.array(-1.0, 2.0, 3.0))

Parameters

move (tuple) – move to accept

energy

Total system energy, valid after a call to initialise.

eval_field(points)

For a given (N, 3) NumPy, dtype=ctypes.c_double array of positions computes the potential field at the given positions. N.B. The value of energy_unit does not effect the output of this function.

free()

Free the KMCFMM instance. Should be called before program termination. Failure to call this method may result in subsequent MPI race conditions.

get_old_energy_with_dats(masks, energy)

Get the energy of a set of particles, useful for computing extraction and recombination rates. Currently implemented as a proposed hop where the new charge has a value of zero.

initialise()

Perform the initial FMM solve and initialise the KMCFMM instance. Must be called before propose or accept can be called.

propose(moves)

Propose moves by providing the local index of the particle and proposed new sites. NB. Returns system energy of proposed moves.

Parameters

moves (tuple) – Tuple of moves where each element is a pair (id, moves_id) where id is a particle id and moves_id is a NumPy array of new positions.

propose_with_dats(site_max_counts, current_sites, prop_positions, prop_masks, prop_energy_diffs, diff=True, prop_charges=None)

Compute the energy difference that would occur if a proposed move was performed. For this interface each charge exists at a site type. Each site type has a maximum number of proposed moves associated with it. The proposed moves are passed in a ParticleDat along with a mask. The change in system energy is returned in the ParticleDat prop_energy_diffs.

A detailed explanation of this interface, including example bookkeeping algorithms, is given in the paper Fast electrostatic solvers for kinetic Monte Carlo simulations.

Parameters
  • site_max_counts (ScalarArray) – ScalarArray(dtype=c_int64) (Input). Max moves associated with site type.

  • current_sites (ParticleDat) – ParticleDat(dtype=c_int64) (Input). Store which site each charge is currently on.

  • prop_positions (ParticleDat) – ParticleDat(dtype=c_double) (Input). Set of proposed positions for each charge.

  • prop_masks (ParticleDat) – ParticleDat(dtype=c_int64) (Input). Energy is computed for moves with mask > 0.

  • prop_energy_diffs (ParticleDat) – ParticleDat(dtype=c_double) (Output. Energy (difference) of each proposed move with mask > 0.

  • diff (bool) – default True. If diff is True then the prop_energy_diffs ParticleDat is populated with energy differences.

  • prop_charges (ParticleDat) – ParticleDat(dtype=c_double) (Input). (optional) Set of proposed new charges.