Table Of Contents

Previous topic

pysic_interface (PyInterface.f90)

Next topic

potentials (Potentials.f90)

This Page

pysic_core (Core.f90)

Core, true to its name, is the heart of the Fortran core of Pysic. It contains the data structures defining the simulation geometry and interactions and defines the central routines for calculating the total energy of and the forces acting on the system.

Many of the routines in pysic_interface (PyInterface.f90) which f2py interfaces to Python are simply calling routines here.

List of subroutines in pysic_core

Full documentation of global variables in pysic_core

atoms

type(atom) pointer size(:)

an array of atom objects representing the system

atoms_created

logical scalar

initial value = .false.

logical tag indicating if atom storing arrays have been created

bo_factors

double precision pointer size(:)

bo_gradients

double precision pointer size(:, :, :)

bo_scaling

logical pointer size(:)

bo_sums

double precision pointer size(:)

bo_temp

double precision pointer size(:)

bond_factors

type(bond_order_parameters) pointer size(:)

an array of bond_order_parameters objects representing bond order factors modifying the potentials

bond_factors_allocated

logical scalar

initial value = .false.

logical tag indicating if bond order parameter storing arrays have been allocated

bond_storage_allocated

logical scalar

initial value = .false.

logical tag indicating if bond order factor storing arrays have been allocated

cell

type(supercell) scalar

a supercell object representing the simulation cell

electronegativity_evaluation_index

integer scalar parameter

initial value = 3

energy_evaluation_index

integer scalar parameter

initial value = 1

evaluate_ewald

logical scalar

initial value = .false.

switch for enabling Ewald summation of coulomb interactions

ewald_allocated

logical scalar

initial value = .false.

ewald_cutoff

double precision scalar

ewald_epsilon

double precision scalar

ewald_k_cutoffs

integer size(3)

ewald_k_radius

double precision scalar

ewald_scaler

double precision pointer size(:)

ewald_sigma

double precision scalar

force_evaluation_index

integer scalar parameter

initial value = 2

group_index_save_slot

integer pointer size(:)

interactions

type(potential) pointer size(:)

an array of potential objects representing the interactions

multipliers

type(potential) allocatable size(:)

a temporary array for storing multiplying potentials before associating them with a master potential

n_bond_factors

integer scalar

initial value = 0

n_interactions

integer scalar

initial value = 0

number of potentials

n_multi

integer scalar

initial value = 0

number of temporary product potentials

n_nbs

integer pointer size(:)

n_saved_bond_order_factors

integer scalar

initial value = 0

number of saved bond order factors

number_of_atoms

integer scalar

potentials_allocated

logical scalar

initial value = .false.

logical tag indicating if potential storing arrays have been allocated

saved_bond_order_factors

double precision pointer size(:, :)

Array for storing calculated bond order factors. Indexing: (atom index, group_index_save_slot(group index))

saved_bond_order_gradients

double precision pointer size(:, :, :, :)

Array for storing calculated bond order gradients. Indexing: (xyz, atom index, group_index_save_slot(group index), target index)

saved_bond_order_sums

double precision pointer size(:, :)

Array for storing calculated bond order sums. Indexing: (atom index, group_index_save_slot(group index))

saved_bond_order_virials

double precision pointer size(:, :, :)

Array for storing calculated bond order virials. Indexing: (xyz, group_index_save_slot(group index), target index)

temp_enegs

double precision pointer size(:)

temp_forces

double precision pointer size(:, :)

temp_gradient

double precision pointer size(:, :, :)

total_n_nbs

integer pointer size(:)

use_saved_bond_order_factors

logical scalar

initial value = .false.

Logical tag which enables / disables bond order saving. If true, bond order calculation routines try to find the precalculated factors in the saved bond order arrays instead of calculating.

use_saved_bond_order_gradients

integer pointer size(:, :)

Array storing the atom index of the bond gradient stored for indices (group index, target index). Since gradients are needed for all factors (N) with respect to moving all atoms (N), storing them all would require an N x N matrix. Therefore only some are stored. This array is used for searching the stroage to see if the needed gradient is there or needs to be calculated.

Full documentation of subroutines in pysic_core

core_add_bond_order_factor(n_targets, n_params, n_split, bond_name, parameters, param_split, cutoff, smooth_cut, elements, orig_elements, group_index, success)

Creates one additional bond_order_factor in the core. The routine assumes that adequate memory has been allocated already using core_allocate_bond_order_factors.

When the bond order parameters in the Python interface are imported to the Fortran core, the target specifiers (elements) are permutated to create all equivalent bond order parameters. That is, if we have parameters for Si-O, both Si-O and O-Si parameters are created. This is because the energy and force calculation loops only deal with atom pairs A-B once (so only A-B or B-A is considered, not both) and if, say, the loop only finds an O-Si pair, it is important to apply the Si-O parameters also on that pair. In some cases, such as with the tersoff factor affecting triplets (A-B-C), the contribution is not symmetric for all the atoms. Therefore it is necessary to also store the original targets of the potential as specified in the Python interface. These are to be given in the ‘orig_elements’ lists.

called from PyInterface: add_bond_order_factor()

Parameters:

n_targets: integer intent(in) scalar
number of targets (interacting bodies)
n_params: integer intent(in) scalar
number of parameters
n_split: integer intent(in) scalar
number of subsets in the list of parameters, should equal n_targets
bond_name: character(len=*) intent(in) scalar
bond order factor names
parameters: double precision intent(in) size(n_params)
numeric parameters
param_split: integer intent(in) size(n_split)
the numbers of parameters for 1-body, 2-body etc.
cutoff: double precision intent(in) scalar
interaction hard cutoff
smooth_cut: double precision intent(in) scalar
interaction soft cutoff
elements: character(len=label_length) intent(in) size(n_targets)
atomic symbols specifying the elements the interaction acts on
orig_elements: character(len=label_length) intent(in) size(n_targets)
original atomic symbols specifying the elements the interaction acts on
group_index: integer intent(in) scalar
index denoting the potential to which the factor is connected
success: logical intent(out) scalar
logical tag specifying if creation of the factor succeeded
core_add_bond_order_forces(group_index, atom_index, prefactor, forces, stress)

Parameters:

group_index: integer intent(in) scalar

atom_index: integer intent(in) scalar

prefactor: double precision intent(in) scalar

forces: double precision intent(inout) size(:, :)

stress: double precision intent(inout) size(6)

core_add_pair_bond_order_forces(index1, index2, prefactor, separation, direction, distance, group_index, pair_bo_sums, pair_bo_factors, forces, stress)

Evaluates the local force affecting two atoms from bond order factors.

Parameters:

index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2

prefactor: double precision intent(in) scalar

separation: double precision intent(in) size(3)

direction: double precision intent(in) size(3)

distance: double precision intent(in) scalar

group_index: integer intent(in) scalar

pair_bo_sums: double precision intent(in) size(2)

pair_bo_factors: double precision intent(in) size(2)

forces: double precision intent(inout) size(:, :)
calculated forces
stress: double precision intent(inout) size(6)
calculated stress
core_add_potential(n_targets, n_params, pot_name, parameters, cutoff, smooth_cut, elements, tags, indices, orig_elements, orig_tags, orig_indices, pot_index, is_multiplier, success)

Creates one additional potential in the core. The routine assumes that adequate memory has been allocated already using core_allocate_potentials.

When the potentials in the Python interface are imported to the Fortran core, the target specifiers (elements, tags, indices) are permutated to create all equivalent potentials. That is, if we have a potential for Si-O, both Si-O and O-Si potentials are created. This is because the energy and force calculation loops only deal with atom pairs A-B once (so only A-B or B-A is considered, not both) and if, say, the loop only finds an O-Si pair, it is important to apply the Si-O interaction also on that pair. In some cases, such as with the bond-bending potential affecting triplets (A-B-C), the interaction is not symmetric for all the atoms. Therefore it is necessary to also store the original targets of the potential as specified in the Python interface. These are to be given in the ‘orig_*’ lists.

If product potentials are created, all but the first one of the potentials are created with is_multiplier == .true.. This leads to the potentials being stored in the global temporary array multipliers. The last potential of a group should be created with is_multiplier = .false. and the stored multipliers are attached to it. The list of multipliers is not cleared automatically, since usually one creates copies of the same potential with permutated targets and all of these need the same multipiers. Instead the multipliers are cleared with a call of clear_potential_multipliers().

called from PyInterface: add_potential()

Parameters:

n_targets: integer intent(in) scalar
number of targets (interacting bodies)
n_params: integer intent(in) scalar
number of parameters
pot_name: character(len=*) intent(in) scalar
potential names
parameters: double precision intent(in) size(n_params)
numeric parameters
cutoff: double precision intent(in) scalar
interaction hard cutoff
smooth_cut: double precision intent(in) scalar
interaction soft cutoff
elements: character(len=label_length) intent(in) size(n_targets)
atomic symbols specifying the elements the interaction acts on
tags: integer intent(in) size(n_targets)
tags specifying the atoms the interaction acts on
indices: integer intent(in) size(n_targets)
indices specifying the atoms the interaction acts on
orig_elements: character(len=label_length) intent(in) size(n_targets)
original atomic symbols specifying the elements the interaction acts on
orig_tags: integer intent(in) size(n_targets)
original tags specifying the atoms the interaction acts on
orig_indices: integer intent(in) size(n_targets)
original indices specifying the atoms the interaction acts on
pot_index: integer intent(in) scalar
index of the potential
is_multiplier: logical intent(in) scalar
logical tag specifying if this potential should be treated as a multiplier
success: logical intent(out) scalar
logical tag specifying if creation of the potential succeeded
core_allocate_bond_order_factors(n_bond_factors)

Allocates pointers for storing bond order factors.

called from PyInterface: allocate_bond_order_factors()

Parameters:

n_bond_factors: integer intent(in) scalar

core_allocate_bond_order_storage(n_atoms, n_groups, n_factors)

Allocates arrays for storing precalculated values of bond order factors and gradients.

called from PyInterface: allocate_bond_order_factors()

Parameters:

n_atoms: integer intent(in) scalar
number of atoms
n_groups: integer intent(in) scalar
number of bond order groups
n_factors: integer intent(in) scalar
number of bond order parameters
core_allocate_potentials(n_pots)

Allocates pointers for storing potentials.

called from PyInterface: allocate_potentials()

Parameters:

n_pots: integer intent(in) scalar
number of potentials
core_assign_bond_order_factor_indices()

This routine finds for each atom the potentials for which the atom is an accepted target at the first position. First position here means that for instance in an A-B-C triplet. A is in first position. Being an accepted target means that the atom has the correct element.

called from PyInterface: create_bond_order_factor_list()

core_assign_potential_indices()

This routine finds for each atom the potentials for which the atom is an accepted target at the first position. First position here means that for instance in an A-B-C triplet. A is in first position. Being an accepted target means that the atom has the correct element, index or tag (one that the potential targets).

called from PyInterface: create_potential_list()

core_build_neighbor_lists(cutoffs)

Builds the neighbor lists in the core. The simulation cell must be partitioned with core_create_space_partitioning() before this routine can be called.

Parameters:

cutoffs: double precision intent(in) size(:)
list of cutoffs, atom by atom
core_calculate_bond_order_factors(group_index, total_bond_orders)

Calculates the bond order sums of all atoms for the given group.

For a factor such as

\[b_i = f(\sum_j c_{ij})\]

The routine calculates

\[\sum_j c_{ij}.\]

The full bond order factor is then obtained by applying the scaling function \(f\). This is done with core_post_process_bond_order_factors().

Parameters:

group_index: integer intent(in) scalar
an index denoting the potential to which the factor is connected
total_bond_orders: double precision intent(inout) size(:)
the calculated bond order sums
core_calculate_bond_order_gradients(group_index, atom_index, raw_sums, total_gradient, total_virial, for_factor)

Returns the gradients of bond order factors.

For a factor such as

\[b_i = f(\sum_j c_{ij})\]

The routine calculates

\[\nabla_\alpha b_i = f'(\sum_j c_{ij}) \nabla_\alpha \sum_j c_{ij}.\]

By default, the gradients of all factors \(i\) are calculated with respect to moving the given atom \(\alpha\). If for_factor is .true., the gradients of the bond factor of the given atom are calculated with respect to moving all atoms.

Parameters:

group_index: integer intent(in) scalar
an index denoting the potential to which the factor is connected
atom_index: integer intent(in) scalar
index of the atom with respect to which the factors are differentiated (\(\alpha\)), or the atoms whose factor is differentiated (\(i\)) if for_factor is .true.
raw_sums: double precision intent(in) size(:)
precalculated bond order sums, \(\sum_j c_{ij}\), in the above example.
total_gradient: double precision intent(inout) size(:, :)
the calculated bond order gradients \(\nabla_\alpha b_i\)
total_virial: double precision intent(inout) size(6)
the components of the virial due to the bond order gradients
for_factor: logical intent(in) scalar optional
a switch for requesting the gradients for a given \(i\) instead of a given \(\alpha\)
core_calculate_bond_order_gradients_of_factor(group_index, atom_index, raw_sums, total_gradient, total_virial)

Returns the gradients of one bond order factor with respect to moving all atoms.

This calls core_calculate_bond_order_gradients() with for_factor = .true.

For a factor such as

\[b_i = f(\sum_j c_{ij})\]

The routine calculates

\[\nabla_\alpha b_i = f'(\sum_j c_{ij}) \nabla_\alpha \sum_j c_{ij}.\]

The gradients of the bond factor of the given atom \(i\) are calculated with respect to moving all atoms \(\alpha\).

Parameters:

group_index: integer intent(in) scalar
an index denoting the potential to which the factor is connected
atom_index: integer intent(in) scalar
index of the atom whose factor is differentiated (\(i\))
raw_sums: double precision intent(in) size(:)
precalculated bond order sums, \(\sum_j c_{ij}\), in the above example.
total_gradient: double precision intent(inout) size(:, :)
the calculated bond order gradients \(\nabla_\alpha b_i\)
total_virial: double precision intent(inout) size(6)
the components of the virial due to the bond order gradient
core_calculate_electronegativities(total_enegs)

Calculates electronegativity forces acting on all atomic charges of the system.

The routine calculates the electronegativities

\[\chi_{\alpha} = -\frac{\partial V}{\partial q_\alpha}\]

for all atoms \(\alpha\). This is done according to the the structure and potentials allocated in the core, so the routine does not accept arguments. Instead, the core modifying routines such as core_generate_atoms() must be called first to set up the calculation.

called from PyInterface: calculate_electronegativities()

Parameters:

total_enegs: double precision intent(inout) size(:)
an array containing the calculated charge forces for all atoms
core_calculate_energy(total_energy)

Calculates the total potential energy of the system.

This is done according to the the structure and potentials allocated in the core, so the routine does not accept arguments. Instead, the core modifying routines such as core_generate_atoms() must be called first to set up the calculation.

called from PyInterface: calculate_energy()

Parameters:

total_energy: double precision intent(out) scalar
calculated total potential energy
core_calculate_forces(total_forces, total_stress)

Calculates forces acting on all atoms of the system.

The routine calculates the potential gradient

\[\mathbf{F}_\alpha = - \nabla_\alpha V\]

for all atoms \(\alpha\). This is done according to the the structure and potentials allocated in the core, so the routine does not accept arguments. Instead, the core modifying routines such as core_generate_atoms() must be called first to set up the calculation.

called from PyInterface: calculate_forces()

Parameters:

total_forces: double precision intent(inout) size(:, :)
an array containing the calculated forces for all atoms
total_stress: double precision intent(inout) size(6)
as array containing the calculated stress tensor
core_calculate_pair_bond_order_factor(atom_pair, separation, distance, direction, group_index, bond_order_sum)

Calculates the bond order sum for a given pair of atoms for the given group.

For a factor such as

\[b_ij = f(\sum_k c_{ijk})\]

The routine calculates

\[\sum_k c_{ijk}.\]

The full bond order factor is then obtained by applying the scaling function \(f\). This is done with core_post_process_bond_order_factors().

Parameters:

atom_pair: integer intent(in) size(2)

separation: double precision intent(in) size(3)

distance: double precision intent(in) scalar

direction: double precision intent(in) size(3)

group_index: integer intent(in) scalar
an index denoting the potential to which the factor is connected
bond_order_sum: double precision intent(out) size(2)
the calculated bond order sums
core_calculate_pair_bond_order_gradients(atom_pair, separation, distance, direction, group_index, raw_sums, total_gradient, total_virial)

Returns the gradients of a pair bond order factor.

For a factor such as

\[b_{ij} = f(\sum_k c_{ijk})\]

The routine calculates

\[\nabla_\alpha b_{ij} = f'(\sum_k c_{ijk}) \nabla_\alpha \sum_k c_{ijk}.\]

By default, the gradients the factor \(ij\) is calculated with respect to moving all atoms \(\alpha\).

Parameters:

atom_pair: integer intent(in) size(2)

separation: double precision intent(in) size(3)

distance: double precision intent(in) scalar

direction: double precision intent(in) size(3)

group_index: integer intent(in) scalar
an index denoting the potential to which the factor is connected
raw_sums: double precision intent(in) size(2)
precalculated bond order sums, \(\sum_j c_{ij}\), in the above example.
total_gradient: double precision intent(inout) size(:, :, :)
the calculated bond order gradients \(\nabla_\alpha b_i\)
total_virial: double precision intent(inout) size(6, 2)
the components of the virial due to the bond order gradient
core_clear_atoms()

Deallocates the array of atoms in the core, if allocated.

core_clear_bond_order_factors()

Deallocates pointers for bond order factors (the parameters)

core_clear_bond_order_storage()

Deallocates pointers for bond order factors (the precalculated factor values).

core_clear_ewald_arrays()
core_clear_potential_multipliers()
core_clear_potentials()

Deallocates pointers for potentials

core_create_cell(vectors, inverse, periodicity)

Creates a supercell for containing the calculation geometry.

called from PyInterface: create_cell()

Parameters:

vectors: double precision intent(in) size(3, 3)
A 3x3 matrix containing the vectors spanning the supercell. The first index runs over xyz and the second index runs over the three vectors.
inverse: double precision intent(in) size(3, 3)
A 3x3 matrix containing the inverse matrix of the one given in vectors, i.e. \(A*B = I\) for the two matrices. Since the latter represents a cell of non-zero volume, this inverse must exist. It is not tested that the given matrix actually is the inverse, the user must make sure it is.
periodicity: logical intent(in) size(3)
A 3-element vector containing logical tags specifying if the system is periodic in the directions of the three vectors spanning the supercell.
core_create_neighbor_list(n_nbors, atom_index, neighbors, offsets)

Assigns a precalculated neighbor list to a single atom of the given index. The neighbor list must be precalculated, this method only stores them in the core. The list must contain an array storing the indices of the neighboring atoms as well as the supercell offsets. The offsets are integer triplets showing how many times must the supercell vectors be added to the position of the neighbor to find the neighboring image in a periodic system. For example, let the supercell be:

[[1.0, 0, 0], [0, 1.0, 0], [0, 0, 1.0]],

i.e., a unit cube, with periodic boundaries. Now, if we have particles with coordinates:

a = [1.5, 0.5, 0.5]
b = [0.4, 1.6, 3.3]

the closest separation vector \(\mathbf{r}_b-\mathbf{r}_a\) between the particles is:

[-.1, .1, -.2]

obtained if we add the vector of periodicity:

[1.0, -1.0, -3.0]

to the coordinates of particle b. The offset vector (for particle b, when listing neighbors of a) is then:

[1, -1, -3]

Note that if the system is small, one atom can in principle appear several times in the neighbor list with different offsets.

called from PyInterface: create_neighbor_list()

Parameters:

n_nbors: integer intent(in) scalar

atom_index: integer intent(in) scalar
index of the atom for which the neighbor list is created
neighbors: integer intent(in) size(n_nbors)
An array containing the indices of the neighboring atoms
offsets: integer intent(in) size(3, n_nbors)
An array containing vectors specifying the offsets of the neighbors in periodic systems.
core_create_space_partitioning(max_cutoff)

Partitions the simulation volume in subvolumes for fast neighbor searching

Parameters:

max_cutoff: double precision intent(in) scalar
the maximum cutoff radius for neighbor search
core_debug_dump(forces)

Write atomic coordinates and other info in a file. This is only for debugging.

Parameters:

forces: double precision intent(in) size(:, :)

core_empty_bond_order_gradient_storage(index)

Clears bond order factor gradients (the precalculated gradient values) but does not deallocate the arrays. If an index is given, then only that column is emptied.

Parameters:

index: integer intent(in) scalar optional
the column to be emptied
core_empty_bond_order_storage()

Clears bond order factors (the precalculated factor values) but does not deallocate the arrays.

core_evaluate_local_doublet(n_atoms, atom_doublet, index1, index2, test_index1, interaction_indices, separations, directions, distances, calculation_type, energy, forces, enegs, stress, many_bodies_found)

Evaluates the interactions affecting two atoms.

Parameters:

n_atoms: integer intent(in) scalar
total number of atoms in the system
atom_doublet: type(atom) intent(in) size(2)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
test_index1: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2
interaction_indices: integer intent() pointer size(:)
the interactions targeting the given atoms
separations: double precision intent(in) size(3, 1)
distance vector from 1 to 2, as an array
directions: double precision intent(in) size(3, 1)
unit vector from 1 to 2, as an array
distances: double precision intent(in) size(1)
distance from 1 to 2, as an array
calculation_type: integer intent(in) scalar
the type of information requested
energy: double precision intent(inout) scalar
calculated energy
forces: double precision intent(inout) size(3, n_atoms)
calculated forces
enegs: double precision intent(inout) size(n_atoms)
calculated electronegativities
stress: double precision intent(inout) size(6)
calculated stress
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets
core_evaluate_local_doublet_electronegativities(n_atoms, atom_doublet, index1, index2, test_index1, interaction_indices, separations, directions, distances, enegs, many_bodies_found)

Evaluates the local electronegativity affecting two atoms.

Parameters:

n_atoms: integer intent(in) scalar

atom_doublet: type(atom) intent(in) size(2)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
test_index1: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2
interaction_indices: integer intent() pointer size(:)
the interactions targeting the given atoms
separations: double precision intent(in) size(3, 1)
distance vector from 1 to 2, as an array
directions: double precision intent(in) size(3, 1)
unit vector from 1 to 2, as an array
distances: double precision intent(in) size(1)
distance from 1 to 2, as an array
enegs: double precision intent(inout) size(n_atoms)
calculated electronegativities
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets
core_evaluate_local_doublet_electronegativities_B(atom_doublet, index1, index2, test_index1, interaction_indices, separations, directions, distances, enegs, many_bodies_found, manybody_indices, n_manybody)

Evaluates the local electronegativity affecting two atoms. (Rearranged internally.)

Parameters:

atom_doublet: type(atom) intent(in) size(2)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
test_index1: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2
interaction_indices: integer intent() pointer size(:)
the interactions targeting the given atoms
separations: double precision intent(in) size(3, 1)
distance vector from 1 to 2, as an array
directions: double precision intent(in) size(3, 1)
unit vector from 1 to 2, as an array
distances: double precision intent(in) size(1)
distance from 1 to 2, as an array
enegs: double precision intent(inout) size(:)
calculated electronegativities
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets

manybody_indices: integer intent() pointer size(:)

n_manybody: integer intent(out) scalar

core_evaluate_local_doublet_energy(n_atoms, atom_doublet, index1, index2, test_index1, interaction_indices, separations, directions, distances, energy, many_bodies_found)

Evaluates the local potential affecting two atoms.

Parameters:

n_atoms: integer intent(in) scalar

atom_doublet: type(atom) intent(in) size(2)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
test_index1: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2
interaction_indices: integer intent() pointer size(:)
the interactions targeting the given atoms
separations: double precision intent(in) size(3, 1)
distance vector from 1 to 2, as an array
directions: double precision intent(in) size(3, 1)
unit vector from 1 to 2, as an array
distances: double precision intent(in) size(1)
distance from 1 to 2, as an array
energy: double precision intent(inout) scalar
calculated energy
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets
core_evaluate_local_doublet_energy_B(atom_doublet, index1, index2, test_index1, interaction_indices, separations, directions, distances, energy, many_bodies_found, manybody_indices, n_manybody)

Evaluates the local potential affecting two atoms. (Rearranged internally compared to ‘A’.)

Parameters:

atom_doublet: type(atom) intent(in) size(2)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
test_index1: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2
interaction_indices: integer intent() pointer size(:)
the interactions targeting the given atoms
separations: double precision intent(in) size(3, 1)
distance vector from 1 to 2, as an array
directions: double precision intent(in) size(3, 1)
unit vector from 1 to 2, as an array
distances: double precision intent(in) size(1)
distance from 1 to 2, as an array
energy: double precision intent(inout) scalar
calculated energy
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets

manybody_indices: integer intent() pointer size(:)

n_manybody: integer intent(out) scalar

core_evaluate_local_doublet_forces(n_atoms, atom_doublet, index1, index2, test_index1, interaction_indices, separations, directions, distances, forces, stress, many_bodies_found)

Evaluates the local force affecting two atoms.

Parameters:

n_atoms: integer intent(in) scalar
total number of atoms in the system
atom_doublet: type(atom) intent(in) size(2)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
test_index1: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2
interaction_indices: integer intent() pointer size(:)
the interactions targeting the given atoms
separations: double precision intent(in) size(3, 1)
distance vector from 1 to 2, as an array
directions: double precision intent(in) size(3, 1)
unit vector from 1 to 2, as an array
distances: double precision intent(in) size(1)
distance from 1 to 2, as an array
forces: double precision intent(inout) size(3, n_atoms)
calculated forces
stress: double precision intent(inout) size(6)
calculated stress
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets
core_evaluate_local_doublet_forces_B(atom_doublet, index1, index2, test_index1, interaction_indices, separations, directions, distances, forces, stress, many_bodies_found, manybody_indices, n_manybody)

Evaluates the local force affecting two atoms. (Rearranged internally.)

Parameters:

atom_doublet: type(atom) intent(in) size(2)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
test_index1: integer intent(in) scalar
if 1, test if the interaction targets atom1; similarly for 2
interaction_indices: integer intent() pointer size(:)
the interactions targeting the given atoms
separations: double precision intent(in) size(3, 1)
distance vector from 1 to 2, as an array
directions: double precision intent(in) size(3, 1)
unit vector from 1 to 2, as an array
distances: double precision intent(in) size(1)
distance from 1 to 2, as an array
forces: double precision intent(inout) size(:, :)
calculated forces
stress: double precision intent(inout) size(6)
calculated stress
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets

manybody_indices: integer intent() pointer size(:)

n_manybody: integer intent(out) scalar

core_evaluate_local_quadruplet(n_atoms, atom_quadruplet, index1, index2, index3, index4, test_index1, test_index2, test_index3, interaction_indices, separations, directions, distances, calculation_type, energy, forces, enegs, stress, many_bodies_found)

Evaluates the interactions affecting four atoms.

Parameters:

n_atoms: integer intent(in) scalar
total number of atoms in the system
atom_quadruplet: type(atom) intent(in) size(4)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
index3: integer intent(in) scalar
index of the atom 3
index4: integer intent(in) scalar
index of the atom 4
test_index1: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2, 3
test_index2: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2, 3
test_index3: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2, 3
interaction_indices: integer intent() pointer size(:)
the interactions targeting the given atoms
separations: double precision intent(in) size(3, 3)
distance vector from 1 to 2, 2 to 3 and 3 to 4 as an array
directions: double precision intent(in) size(3, 3)
unit vector from 1 to 2, 2 to 3 and 3 to 4 as an array
distances: double precision intent(in) size(3)
distance from 1 to 2, 2 to 3 and 3 to 4 as an array
calculation_type: integer intent(in) scalar
the type of information requested
energy: double precision intent(out) scalar
calculated energy
forces: double precision intent(out) size(3, n_atoms)
calculated forces
enegs: double precision intent(out) size(n_atoms)
calculated electronegativities
stress: double precision intent(out) size(6)
calculated stress
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets
core_evaluate_local_quadruplet_B(atom_quadruplet, index1, index2, index3, index4, test_index1, test_index2, test_index3, separations, directions, distances, calculation_type, energy, forces, enegs, stress, many_bodies_found, manybody_indices, n_manybody)

Evaluates the interactions affecting four atoms. (Rearranged internally.)

Parameters:

atom_quadruplet: type(atom) intent(in) size(4)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
index3: integer intent(in) scalar
index of the atom 3
index4: integer intent(in) scalar
index of the atom 4
test_index1: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2, 3
test_index2: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2, 3
test_index3: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2, 3
separations: double precision intent(in) size(3, 3)
distance vector from 1 to 2, 2 to 3 and 3 to 4 as an array
directions: double precision intent(in) size(3, 3)
unit vector from 1 to 2, 2 to 3 and 3 to 4 as an array
distances: double precision intent(in) size(3)
distance from 1 to 2, 2 to 3 and 3 to 4 as an array
calculation_type: integer intent(in) scalar
the type of information requested
energy: double precision intent(inout) scalar
calculated energy
forces: double precision intent(inout) size(:, :)
calculated forces
enegs: double precision intent(inout) size(:)
calculated electronegativities
stress: double precision intent(inout) size(6)
calculated stress
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets

manybody_indices: integer intent() pointer size(:)

n_manybody: integer intent(in) scalar

core_evaluate_local_singlet(index1, atom_singlet, interaction_indices, calculation_type, energy, forces, stress, enegs)

Evaluates the local potential affecting a single atom

Parameters:

index1: integer intent(in) scalar
index of the atom
atom_singlet: type(atom) intent(in) scalar
the atom that is targeted
interaction_indices: integer intent() pointer size(:)
the interactions targeting the given atom
calculation_type: integer intent(in) scalar
specifies if we are evaluating the energy, forces, or electronegativities
energy: double precision intent(inout) scalar
calculated energy
forces: double precision intent(inout) size(:, :)
calculated forces
stress: double precision intent(inout) size(6)
calculated stress
enegs: double precision intent(inout) size(:)
calculated electronegativities
core_evaluate_local_triplet(n_atoms, atom_triplet, index1, index2, index3, test_index1, test_index2, interaction_indices, separations, directions, distances, calculation_type, energy, forces, enegs, stress, many_bodies_found)

Evaluates the interactions affecting three atoms.

Parameters:

n_atoms: integer intent(in) scalar
total number of atoms in the system
atom_triplet: type(atom) intent(in) size(3)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
index3: integer intent(in) scalar
index of the atom 3
test_index1: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2, 3
test_index2: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2, 3
interaction_indices: integer intent() pointer size(:)
the interactions targeting the given atoms
separations: double precision intent(in) size(3, 2)
distance vector from 1 to 2 and 2 to 3 as an array
directions: double precision intent(in) size(3, 2)
unit vector from 1 to 2 and 2 to 3 as an array
distances: double precision intent(in) size(2)
distance from 1 to 2 and 2 to 3 as an array
calculation_type: integer intent(in) scalar
the type of information requested
energy: double precision intent(out) scalar
calculated energy
forces: double precision intent(out) size(3, n_atoms)
calculated forces
enegs: double precision intent(out) size(n_atoms)
calculated electronegativities
stress: double precision intent(out) size(6)
calculated stress
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets
core_evaluate_local_triplet_B(atom_triplet, index1, index2, index3, test_index1, test_index2, separations, directions, distances, calculation_type, energy, forces, enegs, stress, many_bodies_found, manybody_indices, n_manybody)

Evaluates the interactions affecting three atoms. (Rearranged internally.)

Parameters:

atom_triplet: type(atom) intent(in) size(3)
the atoms that are targeted
index1: integer intent(in) scalar
index of the atom 1
index2: integer intent(in) scalar
index of the atom 2
index3: integer intent(in) scalar
index of the atom 3
test_index1: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2, 3
test_index2: integer intent(in) scalar
if 1, test if the ineraction targets atom1; similarly for 2, 3
separations: double precision intent(in) size(3, 2)
distance vector from 1 to 2 and 2 to 3 as an array
directions: double precision intent(in) size(3, 2)
unit vector from 1 to 2 and 2 to 3 as an array
distances: double precision intent(in) size(2)
distance from 1 to 2 and 2 to 3 as an array
calculation_type: integer intent(in) scalar
the type of information requested
energy: double precision intent(inout) scalar
calculated energy
forces: double precision intent(inout) size(:, :)
calculated forces
enegs: double precision intent(inout) size(:)
calculated electronegativities
stress: double precision intent(inout) size(6)
calculated stress
many_bodies_found: logical intent(out) scalar
returns true if the loop finds an interaction with 3 or more targets

manybody_indices: integer intent() pointer size(:)

n_manybody: integer intent(in) scalar

core_fill_bond_order_storage()

Fills the storage for bond order factors and bond order sums. This is meant to be called in the beginning of force and energy evaluation. The routine calculates all bond order factors (in parallel, if run in MPI) and stores them. Then during the energy or force calculation, it is sufficient to just look up the needed values in the arrays. The routine does not calculate and store bond factor gradients.

core_generate_atoms(n_atoms, masses, charges, positions, momenta, tags, elements)

Creates the atomic particles by invoking a subroutine in the geometry module.

called from PyInterface: create_atoms()

Parameters:

n_atoms: integer intent(in) scalar
number of atoms
masses: double precision intent(in) size(n_atoms)
masses of atoms
charges: double precision intent(in) size(n_atoms)
electric charges of atoms
positions: double precision intent(in) size(3, n_atoms)
coordinates of atoms
momenta: double precision intent(in) size(3, n_atoms)
momenta of atoms
tags: integer intent(in) size(n_atoms)
numeric tags for the atoms
elements: character(len=label_length) intent(in) size(n_atoms)
atomic symbols of the atoms
core_get_bond_order_factor_of_atom(group_index, atom_index, bond_order_factor)

Returns the bond order factors of the given atom for the given group.

Parameters:

group_index: integer intent(in) scalar
index for the bond order factor group
atom_index: integer intent(in) scalar
index of the atom whose bond order factor is returned
bond_order_factor: double precision intent(inout) scalar
the calculated bond order factor
core_get_bond_order_factors(group_index, bond_order_factors)

Returns the bond order factors of all atoms for the given group. The routines tries to find the values in the stored precalculated values first if use_saved_bond_order_factors is true, and saves the calculated values if it does not find them.

Parameters:

group_index: integer intent(in) scalar
index for the bond order factor group
bond_order_factors: double precision intent(inout) size(:)
the calculated bond order factors
core_get_bond_order_gradients(group_index, atom_index, slot_index, bond_order_gradients, bond_order_virial)

Returns the gradients of the bond order factor of the given atom with respect to moving all atoms, for the given group. The routine tries to find the values in the stored precalculated values first if use_saved_bond_order_factors is true, and saves the calculated values if it does not find them.

The slot index is the index of the atom in the interaction being evaluated (so for a triplet A-B-C, A would have slot 1, B slot 2, and C slot 3). This is only used for storing the values.

Parameters:

group_index: integer intent(in) scalar
index for the bond order factor group
atom_index: integer intent(in) scalar
index of the atom whose bond order factor is differentiated
slot_index: integer intent(in) scalar
index denoting the position of the atom in an interacting group (such as A-B-C triplet)
bond_order_gradients: double precision intent(inout) size(:, :)
the calculated gradients of the bond order factor
bond_order_virial: double precision intent(inout) size(6)
the components of the virial due to the bond order factors
core_get_bond_order_sums(group_index, bond_order_sums)

Returns the bond order sums of all atoms for the given group. By ‘bond order sum’, we mean the summation of local terms without per atom scaling. E.g., for \(b_i = 1 + \sum c_{ij}\), \(\sum c_{ij}\) is the sum. The routines tries to find the values in the stored precalculated values first if use_saved_bond_order_factors is true, and saves the calculated values if it does not find them.

Parameters:

group_index: integer intent(in) scalar
index for the bond order factor group
bond_order_sums: double precision intent(inout) size(:)
the calculated bond order sums
core_get_cell_vectors(vectors)

Returns the vectors defining the supercell stored in the core.

called from PyInterface: get_cell_vectors()

Parameters:

vectors: double precision intent(out) size(3, 3)
A 3x3 matrix containing the vectors spanning the supercell. The first index runs over xyz and the second index runs over the three vectors.
core_get_ewald_energy(real_cut, k_cut, reciprocal_cut, sigma, epsilon, energy)

Debug routine for Ewald

Parameters:

real_cut: double precision intent(in) scalar

k_cut: double precision intent(in) scalar

reciprocal_cut: integer intent(in) size(3)

sigma: double precision intent(in) scalar

epsilon: double precision intent(in) scalar

energy: double precision intent(out) scalar

core_get_neighbor_list_of_atom(atom_index, n_neighbors, neighbors, offsets)

Returns the list of neighbros for an atom

Parameters:

atom_index: integer intent(in) scalar
the index of the atom whose neighbors are returned
n_neighbors: integer intent(in) scalar
the number of neighbors
neighbors: integer intent(out) size(n_neighbors)
the indices of the neighboring atoms
offsets: integer intent(out) size(3, n_neighbors)
the offsets for periodic boundaries
core_get_number_of_atoms(n_atoms)

Returns the number of atoms in the array allocated in the core.

called from PyInterface: get_number_of_atoms()

Parameters:

n_atoms: integer intent(out) scalar
number of atoms
core_get_number_of_neighbors(atom_index, n_neighbors)

Returns the number of neighbors for an atom

Parameters:

atom_index: integer intent(in) scalar
the index of the atoms
n_neighbors: integer intent(out) scalar
the number of neighbors
core_loop_over_local_interactions(calculation_type, total_energy, total_forces, total_enegs, total_stress)

Loops over atoms, atomic pairs, atomic triplets, and atomic quadruplets and calculates the contributions from local potentials to energy, forces, or electronegativities. This routine is called from the routines

Parameters:

calculation_type: integer intent(in) scalar
index to specify if the loop calculates energies, forces, or e-negativities
total_energy: double precision intent(inout) scalar
calculated energy
total_forces: double precision intent(inout) size(:, :)
calculated forces
total_enegs: double precision intent(inout) size(:)
calculated electronegativities
total_stress: double precision intent(inout) size(6)
calculated stress
core_post_process_bond_order_factors(group_index, raw_sums, total_bond_orders)

Bond-order post processing, i.e., application of per-atom scaling functions.

By post processing, we mean any operations done after calculating the sum of pair- and many-body terms. That is, if a factor is, say,

\[b_i = f(\sum_j c_{ij}) = 1 + \sum_j c_{ij},\]

the \(\sum_j c_{ij}\) would have been calculated already (with core_calculate_bond_order_factors()) and the operation \(f(x) = 1 + x\) remains to be carried out. The post processing is done per atom regardless of if the bond factor is of a pair or many body type.

Parameters:

group_index: integer intent(in) scalar
an index denoting the potential to which the factor is connected
raw_sums: double precision intent(in) size(:)
precalculated bond order sums, \(\sum_j c_{ij}\), in the above example.
total_bond_orders: double precision intent(inout) size(:)
the calculated bond order factors \(b_i\)
core_post_process_bond_order_gradients(group_index, raw_sums, raw_gradients, total_bond_gradients, mpi_split)

Bond-order post processing, i.e., application of per-atom scaling functions. This routine does the scaling for all bond factors with the given bond order sums and gradients of these sums.

By post processing, we mean any operations done after calculating the sum of pair- and many-body terms. That is, if a factor is, say,

\[b_i = f(\sum_j c_{ij}) = 1 + \sum_j c_{ij},\]

the \(\sum_j c_{ij}\) would have been calculated already and the operation \(f(x) = 1 + x\) remains to be carried out. The post processing is done per atom regardless of if the bond factor is of a pair or many body type.

For gradients, one needs to evaluate

\[\nabla_\alpha b_i = f'(\sum_j c_{ij}) \nabla_\alpha \sum_j c_{ij}\]

Parameters:

group_index: integer intent(in) scalar
an index denoting the potential to which the factor is connected
raw_sums: double precision intent(in) size(:)
precalculated bond order sums, \(\sum_j c_{ij}\), in the above example
raw_gradients: double precision intent(in) size(:, :)
precalculated gradients of bond order sums, \(\nabla_\alpha \sum_j c_{ij}\), in the above example
total_bond_gradients: double precision intent(inout) size(:, :)
the calculated bond order gradients \(\nabla_\alpha b_i\)
mpi_split: logical intent(in) scalar optional
A switch for enabling MPI parallelization. By default the routine is sequential since the calculation may be called from within an already parallelized routine.
core_post_process_bond_order_gradients_of_factor(group_index, atom_index, raw_sum, raw_gradients, total_bond_gradients, raw_virial, total_virial, mpi_split)

Bond-order post processing, i.e., application of per-atom scaling functions. This routine does the scaling for the bond order factor of the given atom with respect to moving all atoms with the given bond order sum for the factor and the gradients of the sum with respect to moving all atoms.

By post processing, we mean any operations done after calculating the sum of pair- and many-body terms. That is, if a factor is, say,

\[b_i = f(\sum_j c_{ij}) = 1 + \sum_j c_{ij},\]

the \(\sum_j c_{ij}\) would have been calculated already and the operation \(f(x) = 1 + x\) remains to be carried out. The post processing is done per atom regardless of if the bond factor is of a pair or many body type.

For gradients, one needs to evaluate

\[\nabla_\alpha b_i = f'(\sum_j c_{ij}) \nabla_\alpha \sum_j c_{ij}\]

Parameters:

group_index: integer intent(in) scalar
an index denoting the potential to which the factor is connected
atom_index: integer intent(in) scalar
the index of the atom whose factor is differentiated (\(i\))
raw_sum: double precision intent(in) scalar
precalculated bond order sum for the given atom, \(\sum_j c_{ij}\), in the above example
raw_gradients: double precision intent(in) size(:, :)
precalculated gradients of bond order sums, \(\nabla_\alpha \sum_j c_{ij}\), in the above example
total_bond_gradients: double precision intent(inout) size(:, :)
the calculated bond order gradients \(\nabla_\alpha b_i\)
raw_virial: double precision intent(in) size(6)
the precalculated virial due to the bond order gradient
total_virial: double precision intent(inout) size(6)
the scaled virial due to the bond order gradient
mpi_split: logical intent(in) scalar optional
A switch for enabling MPI parallelization. By default the routine is sequential since the calculation may be called from within an already parallelized routine.
core_post_process_pair_bond_order_factor(atom1, group_index, raw_sum, total_bond_order)

Bond-order post processing, i.e., application of per-pair scaling functions.

By post processing, we mean any operations done after calculating the sum of pair- and many-body terms. That is, if a factor is, say,

\[b_{ij} = f(\sum_k c_{ijk}) = 1 + \sum_k c_{ijk},\]

the \(\sum_k c_{ijk}\) would have been calculated already (with core_calculate_pair_bond_order_factor()) and the operation \(f(x) = 1 + x\) remains to be carried out.

Parameters:

atom1: type(atom) intent(in) scalar
the central atom of the pair bond order factor
group_index: integer intent(in) scalar
an index denoting the potential to which the factor is connected
raw_sum: double precision intent(in) scalar
precalculated bond order sum, \(\sum_k c_{ijk}\), in the above example.
total_bond_order: double precision intent(out) scalar
the calculated bond order factor \(b_{ij}\)
core_post_process_pair_bond_order_gradients(group_index, atom1, raw_sum, raw_gradients, total_bond_gradients, raw_virial, total_virial, mpi_split)

Bond-order post processing, i.e., application of per-pair scaling functions. This routine does the scaling for the bond order factor of the given pair with respect to moving all atoms with the given bond order sum for the factor and the gradients of the sum with respect to moving all atoms.

By post processing, we mean any operations done after calculating the sum of pair- and many-body terms. That is, if a factor is, say,

\[b_{ij} = f(\sum_k c_{ijk}) = 1 + \sum_k c_{ijk},\]

the \(\sum_k c_{ijk}\) would have been calculated already and the operation \(f(x) = 1 + x\) remains to be carried out. The post processing is done per pair.

For gradients, one needs to evaluate

\[\nabla_\alpha b_{ij} = f'(\sum_k c_{ijk}) \nabla_\alpha \sum_k c_{ijk}\]

Parameters:

group_index: integer intent(in) scalar
an index denoting the potential to which the factor is connected
atom1: type(atom) intent(in) scalar
the central atom of the pair bond order factor
raw_sum: double precision intent(in) scalar
precalculated bond order sum for the given atom, \(\sum_j c_{ij}\), in the above example
raw_gradients: double precision intent(in) size(:, :)
precalculated gradients of bond order sums, \(\nabla_\alpha \sum_j c_{ij}\), in the above example
total_bond_gradients: double precision intent(out) size(:, :)
the calculated bond order gradients \(\nabla_\alpha b_i\)
raw_virial: double precision intent(in) size(6)
the precalculated virial due to the bond order gradient
total_virial: double precision intent(out) size(6)
the scaled virial due to the bond order gradient
mpi_split: logical intent(in) scalar optional
A switch for enabling MPI parallelization. By default the routine is sequential since the calculation may be called from within an already parallelized routine.
core_release_all_memory()

Release all allocated pointer arrays in the core.

core_set_ewald_parameters(real_cut, k_radius, reciprocal_cut, sigma, epsilon, scaler)

Sets the parameters for Ewald summation in the core.

Parameters:

real_cut: double precision intent(in) scalar
the real-space cutoff
k_radius: double precision intent(in) scalar
the k-space cutoff (in inverse length)
reciprocal_cut: integer intent(in) size(3)
the k-space cutoffs (in numbers of k-space cells)
sigma: double precision intent(in) scalar
the split parameter
epsilon: double precision intent(in) scalar
electric constant
scaler: double precision intent(in) size(:)
scaling factors for the individual charges
core_update_atom_charges(n_atoms, charges)

Updates the charges of atomic particles.

called from PyInterface: update_atom_charges()

Parameters:

n_atoms: integer intent(in) scalar
number of atoms
charges: double precision intent(in) size(n_atoms)
new charges for the atoms
core_update_atom_coordinates(n_atoms, positions, momenta)

Updates the positions and momenta of atomic particles.

called from PyInterface: update_atom_coordinates()

Parameters:

n_atoms: integer intent(in) scalar
number of atoms
positions: double precision intent(in) size(3, n_atoms)
new coordinates for the atoms
momenta: double precision intent(in) size(3, n_atoms)
new momenta for the atoms
expand_neighbor_storage(nbors_and_offsets, length, new_length, n_atoms)

Expands the allocated memory for storing neighbor lists

Parameters:

nbors_and_offsets: integer intent() pointer size(:, :, :)

length: integer intent(in) scalar

new_length: integer intent(in) scalar

n_atoms: integer intent(in) scalar

list_atoms()

Prints some information on the atoms stored in the core in stdout.

list_bonds()

Prints some information on the bond order factors stored in the core in stdout.

list_cell()

Prints some information on the supercell stored in the core in stdout.

list_interactions()

Prints some information on the potentials stored in the core in stdout.