fullrmc.Core package

Collection

It contains a collection of methods and classes that are useful for the package.

fullrmc.Core.Collection.get_caller_frames(engine, frame, subframeToAll, caller)

Get list of frames for a function caller.

Parameters:
  1. engine (None, Engine): The stochastic engine in consideration.
  2. frame (None, string): The frame name. If engine is given as None, only None will be accepted as frame value.
  3. subframeToAll (boolean): If frame is a subframe then all multiframe subframes must be considered.
  4. caller (string): Caller name for logging and debugging purposes.
Returns:
  1. usedIncluded (boolean): Whether engine used frame is included in the built frames list. If frame is given as None, True will always be returned.
  2. frame (string): The given frame in parameters. If subframe is given and subframeToAll is True, then multiframe is returned.
  3. allFrames (list): List of all frames.
import inspect
from fullrmc.Core.Collection import get_caller_frames

# Assuming self is a constraint and get_caller_frames is called from within a method ...
usedIncluded, frame, allFrames = get_caller_frames(engine=self.engine,
                                                   frame='frame_name',
                                                   subframeToAll=True,
                                                   caller="%s.%s"%(self.__class__.__name__,inspect.stack()[0][3]) )
fullrmc.Core.Collection.get_real_elements_weight(elements, weightsDict, weighting)

Get elements weights given a dictionary of weights and a weighting scheme. If element weight is not defined in weightsDict then weight is fetched from pdbparser elements database using weighting scheme.

Parameters:
  1. elements (list): List of elements.
  2. weightsDict (dict): Dictionary of fixed weights.
  3. weighting (str): Weighting scheme.
Returns:
  1. elementsWeight (dict): Elements weights got from weightsDict and completed using weighting scheme,
fullrmc.Core.Collection.raise_if_collected(func)

Constraints method decorator that raises an error whenever the method is called and the system has atoms that were removed.

fullrmc.Core.Collection.reset_if_collected_out_of_date(func)

Constraints method decorator that resets the constraint whenever the method is called and the system has atoms that were removed.

fullrmc.Core.Collection.is_number(number)

Check if number is convertible to float.

Parameters:
  1. number (str, number): Input number.
Returns:
  1. result (bool): True if convertible, False otherwise
fullrmc.Core.Collection.is_integer(number, precision=1e-09)

Check if number is convertible to integer.

Parameters:
  1. number (str, number): Input number.
  2. precision (number): To avoid floating errors, a precision should be given.
Returns:
  1. result (bool): True if convertible, False otherwise.
fullrmc.Core.Collection.get_elapsed_time(start, format='%d days, %d hours, %d minutes, %d seconds')

Get formatted time elapsed.

Parameters:
  1. start (time.time): A time instance.
  2. format (string): The format string. must contain exactly four ‘%d’.
Returns:
  1. time (string): The formatted elapsed time.
fullrmc.Core.Collection.get_memory_usage()

Get current process memory usage. This is method requires psutils to be installed.

Returns:
  1. memory (float, None): The memory usage in Megabytes. When psutils is not installed, None is returned.
fullrmc.Core.Collection.get_path(key=None)

Get all paths information needed about the running script and python executable path.

Parameters:

#. key (None, string): the path to return. If not None is given, it can take any of the following:

  1. cwd: current working directory
  2. script: the script’s total path
  3. exe: python executable path
  4. script_name: the script name
  5. relative_script_dir: the script’s relative directory path
  6. script_dir: the script’s absolute directory path
  7. fullrmc: fullrmc package path
Returns:
  1. path (dictionary, value): If key is not None it returns the value of paths dictionary key. Otherwise all the dictionary is returned.
fullrmc.Core.Collection.rebin(data, bin=0.05, check=False)

Re-bin 2D data of shape (N,2). In general, fullrmc requires equivalently spaced experimental data bins. This function can be used to recompute any type of experimental data according to a set bin size.

Parameters:
  1. data (numpy.ndarray): The (N,2) shape data where first column is considered experimental data space values (e.g. r, q) and second column experimental data values.
  2. bin (number): New desired bin size.
  3. check (boolean): whether to check arguments before rebining.
Returns:
  1. X (numpy.ndarray): First column re-binned.
  2. Y (numpy.ndarray): Second column re-binned.
fullrmc.Core.Collection.smooth(data, winLen=11, window='hanning', check=False)

Smooth 1D data using window function and length.

Parameters:
  1. data (numpy.ndarray): the 1D numpy data.
  2. winLen (integer): the smoothing window length.
  3. window (str): The smoothing window type. Can be anything among ‘flat’, ‘hanning’, ‘hamming’, ‘bartlett’ and ‘blackman’.
  4. check (boolean): whether to check arguments before smoothing data.
Returns:
  1. smoothed (numpy.ndarray): the smoothed 1D data array.
fullrmc.Core.Collection.get_random_perpendicular_vector(vector)

Get random normalized perpendicular vector to a given vector.

Parameters:
  1. vector (numpy.ndarray, list, set, tuple): Given vector to compute a random perpendicular vector to it.
Returns:
  1. perpVector (numpy.ndarray): Perpendicular vector of type fullrmc.Globals.FLOAT_TYPE
fullrmc.Core.Collection.get_principal_axis(coordinates, weights=None)

Calculate principal axis of a set of atoms coordinates.

Parameters:
  1. coordinates (np.ndarray): Atoms (N,3) coordinates array.
  2. weights (numpy.ndarray, None): List of weights to compute the weighted Center Of Mass (COM) calculation. Must be a numpy.ndarray of numbers of the same length as indexes. None is accepted for equivalent weighting.
Returns:
  1. center (numpy.ndarray): the weighted COM of the atoms.
  2. eval1 (fullrmc.Globals.FLOAT_TYPE): Biggest eigen value.
  3. eval2 (fullrmc.Globals.FLOAT_TYPE): Second biggest eigen value.
  4. eval3 (fullrmc.Globals.FLOAT_TYPE): Smallest eigen value.
  5. axis1 (numpy.ndarray): Principal axis corresponding to the biggest eigen value.
  6. axis2 (numpy.ndarray): Principal axis corresponding to the second biggest eigen value.
  7. axis3 (numpy.ndarray): Principal axis corresponding to the smallest eigen value.
fullrmc.Core.Collection.get_rotation_matrix(rotationVector, angle)

Calculate the rotation (3X3) matrix about an axis (rotationVector) by a rotation angle.

Parameters:
  1. rotationVector (list, tuple, numpy.ndarray): Rotation axis coordinates.
  2. angle (float): Rotation angle in rad.
Returns:
  1. rotationMatrix (numpy.ndarray): Computed (3X3) rotation matrix
fullrmc.Core.Collection.rotate(xyzArray, rotationMatrix)

Rotate (N,3) numpy.array using a rotation matrix. The array itself will be rotated and not a copy of it.

Parameters:
  1. indexes (numpy.ndarray): the xyz (N,3) array to rotate.
  2. rotationMatrix (numpy.ndarray): the (3X3) rotation matrix.
fullrmc.Core.Collection.get_orientation_matrix(arrayAxis, alignToAxis)

Get the rotation matrix that aligns arrayAxis to alignToAxis

Parameters:
  1. arrayAxis (list, tuple, numpy.ndarray): xyzArray axis.
  2. alignToAxis (list, tuple, numpy.ndarray): The axis to align to.
fullrmc.Core.Collection.orient(xyzArray, arrayAxis, alignToAxis)

Rotates xyzArray using the rotation matrix that rotates and aligns arrayAxis to alignToAXis.

Parameters:
  1. xyzArray (numpy.ndarray): The xyz (N,3) array to rotate.
  2. arrayAxis (list, tuple, numpy.ndarray): xyzArray axis.
  3. alignToAxis (list, tuple, numpy.ndarray): The axis to align to.
fullrmc.Core.Collection.get_superposition_transformation(refArray, array, check=False)

Calculate the rotation tensor and the translations that minimizes the root mean square deviation between an array of vectors and a reference array.

Parameters:
  1. refArray (numpy.ndarray): The NX3 reference array to superpose to.
  2. array (numpy.ndarray): The NX3 array to calculate the transformation of.
  3. check (boolean): Whether to check arguments before generating points.
Returns:
  1. rotationMatrix (numpy.ndarray): The 3X3 rotation tensor.
  2. refArrayCOM (numpy.ndarray): The 1X3 vector center of mass of refArray.
  3. arrayCOM (numpy.ndarray): The 1X3 vector center of mass of array.
  4. rms (number)
fullrmc.Core.Collection.superpose_array(refArray, array, check=False)

Superpose arrays by calculating the rotation matrix and the translations that minimize the root mean square deviation between and array of vectors and a reference array.

Parameters:
  1. refArray (numpy.ndarray): the NX3 reference array to superpose to.
  2. array (numpy.ndarray): the NX3 array to calculate the transformation of.
  3. check (boolean): whether to check arguments before generating points.
Returns:
  1. superposedArray (numpy.ndarray): the NX3 array to superposed array.
fullrmc.Core.Collection.generate_random_vector(minAmp, maxAmp)

Generate random vector in 3D.

Parameters:
  1. minAmp (number): Vector minimum amplitude.
  2. maxAmp (number): Vector maximum amplitude.
Returns:
  1. vector (numpy.ndarray): the vector [X,Y,Z] array
fullrmc.Core.Collection.generate_points_on_sphere(thetaFrom, thetaTo, phiFrom, phiTo, npoints=1, check=False)

Generate random points on a sphere of radius 1. Points are generated using spherical coordinates arguments as in figure below. Theta [0,Pi] is the angle between the generated point and Z axis. Phi [0,2Pi] is the angle between the generated point and x axis.

_images/sphericalCoordsSystem.png
Parameters:
  1. thetaFrom (number): The minimum theta value.
  2. thetaTo (number): The maximum theta value.
  3. phiFrom (number): The minimum phi value.
  4. phiTo (number): The maximum phi value.
  5. npoints (integer): The number of points to generate
  6. check (boolean): whether to check arguments before generating points.
Returns:
  1. x (numpy.ndarray): The (npoints,1) numpy array of all generated points x coordinates.
  2. y (numpy.ndarray): The (npoints,1) numpy array of all generated points y coordinates.
  3. z (numpy.ndarray): The (npoints,1) numpy array of all generated points z coordinates.
fullrmc.Core.Collection.find_extrema(x, max=True, min=True, strict=False, withend=False)

Get a vector extrema indexes and values.

Parameters:
  1. max (boolean): Whether to index the maxima.
  2. min (boolean): Whether to index the minima.
  3. strict (boolean): Whether not to index changes to zero gradient.
  4. withend (boolean): Whether to always include x[0] and x[-1].
Returns:
  1. indexes (numpy.ndarray): Extrema indexes.
  2. values (numpy.ndarray): Extrema values.
fullrmc.Core.Collection.convert_Gr_to_gr(Gr, minIndex)

Converts G(r) to g(r) by computing the following \(g(r)=1+(\frac{G(r)}{4 \pi \rho_{0} r})\)

Parameters:
  1. Gr (numpy.ndarray): The G(r) numpy array of shape (number of points, 2)
  2. minIndex (int, tuple): The minima indexes to compute the number density rho0. It can be a single peak or a list of peaks to compute the mean slope instead.
Returns:
  1. minimas (numpy.ndarray): The minimas array found using minIndex and used to compute the slope and therefore \(\rho_{0}\).
  2. slope (float): The computed slope from the minimas.
  3. rho0 (float): The number density of the material.
  4. g(r) (numpy.ndarray): the computed g(r).

To visualize convertion

# peak indexes can be different, adjust according to your data
minPeaksIndex = [1,3,4]
minimas, slope, rho0, gr =  convert_Gr_to_gr(Gr, minIndex=minPeaksIndex)
print('slope: %s --> rho0: %s'%(slope,rho0))
import matplotlib.pyplot as plt
line = np.transpose( [[0, Gr[-1,0]], [0, slope*Gr[-1,0]]] )
plt.plot(Gr[:,0],Gr[:,1], label='G(r)')
plt.plot(minimas[:,0], minimas[:,1], 'o', label='minimas')
plt.plot(line[:,0], line[:,1], label='density')
plt.plot(gr[:,0],gr[:,1], label='g(r)')
plt.legend()
plt.show()
fullrmc.Core.Collection.generate_vectors_in_solid_angle(direction, maxAngle, numberOfVectors=1, check=False)

Generate random vectors that satisfy angle condition with a direction vector. Angle between any generated vector and direction must be smaller than given maxAngle.

_images/100randomVectors30deg.png

a) 100 vectors generated around OX axis within a maximum angle separation of 30 degrees.

_images/200randomVectors45deg.png

b) 200 vectors generated around [1,-1,1] axis within a maximum angle separation of 45 degrees.

_images/500randomVectors100deg.png

b) 500 vectors generated around [2,5,1] axis within a maximum angle separation of 100 degrees.

Parameters:
  1. direction (number): The direction around which to create the vectors.
  2. maxAngle (number): The maximum angle allowed.
  3. numberOfVectors (integer): The number of vectors to generate.
  4. check (boolean): whether to check arguments before generating vectors.
Returns:
  1. vectors (numpy.ndarray): The (numberOfVectors,3) numpy array of generated vectors.
fullrmc.Core.Collection.gaussian(x, center=0, FWHM=1, normalize=True, check=True)

Compute the normal distribution or gaussian distribution of a given vector. The probability density of the gaussian distribution is: \(f(x,\mu,\sigma) = \frac{1}{\sigma\sqrt{2\pi}} e^{\frac{-(x-\mu)^{2}}{2\sigma^2}}\)

Where:

  • \(\mu\) is the center of the gaussian, it is the mean or expectation of the distribution it is called the distribution’s median or mode.
  • \(\sigma\) is its standard deviation.
  • \(FWHM=2\sqrt{2 ln 2} \sigma\) is the Full Width at Half Maximum of the gaussian.
Parameters:
  1. x (numpy.ndarray): The vector to compute the gaussian
  2. center (number): The center of the gaussian.
  3. FWHM (number): The Full Width at Half Maximum of the gaussian.
  4. normalize(boolean): Whether to normalize the generated gaussian by \(\frac{1}{\sigma\sqrt{2\pi}}\) so the integral is equal to 1.
  5. check (boolean): whether to check arguments before generating vectors.
fullrmc.Core.Collection.step_function(x, center=0, FWHM=0.1, height=1, check=True)

Compute a step function as the cumulative summation of a gaussian distribution of a given vector.

Parameters:
  1. x (numpy.ndarray): The vector to compute the gaussian. gaussian is computed as a function of x.
  2. center (number): The center of the step function which is the the center of the gaussian.
  3. FWHM (number): The Full Width at Half Maximum of the gaussian.
  4. height (number): The height of the step function.
  5. check (boolean): whether to check arguments before generating vectors.
class fullrmc.Core.Collection.ListenerBase

Bases: object

All listeners base class.

listenerId

Listener unique id set at initialization

listen(message, argument=None)

Listens to any message sent from the Broadcaster.

Parameters:
  1. message (object): Any python object to send to constraint’s listen method.
  2. arguments (object): Any python object.
class fullrmc.Core.Collection.Broadcaster

Bases: object

A broadcaster broadcasts a message to all registered listener.

listeners

Listeners list copy.

add_listener(listener)

Add listener to the list of listeners.

Parameters:
  1. listener (object): Any python object having a listen method.
remove_listener(listener)

Remove listener to the list of listeners.

Parameters:
  1. listener (object): The listener object to remove.
broadcast(message, arguments=None)

Broadcast a message to all the listeners

Parameters:
  1. message (object): Any type of message object to pass to the listeners.
  2. arguments (object): Any type of argument to pass to the listeners.
class fullrmc.Core.Collection.RandomFloatGenerator(lowerLimit, upperLimit)

Bases: object

Generate random float number between a lower and an upper limit.

Parameters:
  1. lowerLimit (number): The lower limit allowed.
  2. upperLimit (number): The upper limit allowed.
lowerLimit

Lower limit of the number generation.

upperLimit

Upper limit of the number generation.

rang

Range defined as upperLimit-lowerLimit.

set_lower_limit(lowerLimit)

Set lower limit.

Parameters:
  1. lowerLimit (number): Lower limit allowed.
set_upper_limit(upperLimit)

Set upper limit.

Parameters:
  1. upperLimit (number): Upper limit allowed.
generate()

Generate a random float number between lowerLimit and upperLimit.

class fullrmc.Core.Collection.BiasedRandomFloatGenerator(lowerLimit, upperLimit, weights=None, biasRange=None, biasFWHM=None, biasHeight=1, unbiasRange=None, unbiasFWHM=None, unbiasHeight=None, unbiasThreshold=1)

Bases: fullrmc.Core.Collection.RandomFloatGenerator

Generate biased random float number between a lower and an upper limit. To bias the generator at a certain number, a bias gaussian is added to the weights scheme at the position of this particular number.

_images/biasedFloatGenerator.png
Parameters:
  1. lowerLimit (number): The lower limit allowed.
  2. upperLimit (number): The upper limit allowed.
  3. weights (None, list, numpy.ndarray): The weights scheme. The length defines the number of bins and the edges. The length of weights array defines the resolution of the biased numbers generation. If None is given, ones array of length 10000 is automatically generated.
  4. biasRange(None, number): The bias gaussian range. It must be smaller than half of limits range which is equal to (upperLimit-lowerLimit)/2. If None is given, it will be automatically set to (upperLimit-lowerLimit)/5
  5. biasFWHM(None, number): The bias gaussian Full Width at Half Maximum. It must be smaller than half of biasRange. If None, it will be automatically set to biasRange/10
  6. biasHeight(number): The bias gaussian maximum intensity.
  7. unbiasRange(None, number): The bias gaussian range. It must be smaller than half of limits range which is equal to (upperLimit-lowerLimit)/2. If None is given, it will be automatically set to biasRange.
  8. unbiasFWHM(None, number): The bias gaussian Full Width at Half Maximum. It must be smaller than half of biasRange. If None is given, it will be automatically set to biasFWHM.
  9. unbiasHeight(number): The unbias gaussian maximum intensity. If None is given, it will be automatically set to biasHeight.
  10. unbiasThreshold(number): unbias is only applied at a certain position only when the position weight is above unbiasThreshold. It must be a positive number.
originalWeights

Original weights as initialized.

weights

Current value weights vector.

scheme

Numbers generation scheme.

bins

Number of bins that is equal to the length of weights vector.

binWidth

Bin width defining the resolution of the biased random number generation.

bias

Bias step-function.

biasGuassian

Bias gaussian function.

biasRange

Bias gaussian extent range.

biasBins

Bias gaussian number of bins.

biasFWHM

Bias gaussian Full Width at Half Maximum.

biasFWHMBins

Bias gaussian Full Width at Half Maximum number of bins.

unbias

Unbias step-function.

unbiasGuassian

Unbias gaussian function.

unbiasRange

Unbias gaussian extent range.

unbiasBins

Unbias gaussian number of bins.

unbiasFWHM

Unbias gaussian Full Width at Half Maximum.

unbiasFWHMBins

Unbias gaussian Full Width at Half Maximum number of bins.

set_weights(weights=None)

Set generator’s weights.

Parameters:
  1. weights (None, list, numpy.ndarray): The weights scheme. The length defines the number of bins and the edges. The length of weights array defines the resolution of the biased numbers generation. If None is given, ones array of length 10000 is automatically generated.
set_bias(biasRange, biasFWHM, biasHeight)

Set generator’s bias gaussian function

Parameters:
  1. biasRange(None, number): Bias gaussian range. It must be smaller than half of limits range which is equal to (upperLimit-lowerLimit)/2. If None is given, it will be automatically set to (upperLimit-lowerLimit)/5.
  2. biasFWHM(None, number): Bias gaussian Full Width at Half Maximum. It must be smaller than half of biasRange. If None is given, it will be automatically set to biasRange/10.
  3. biasHeight(number): Bias gaussian maximum intensity.
set_unbias(unbiasRange, unbiasFWHM, unbiasHeight, unbiasThreshold)

Set generator’s unbias gaussian function

Parameters:
  1. unbiasRange(None, number): The bias gaussian range. It must be smaller than half of limits range which is equal to (upperLimit-lowerLimit)/2. If None, it will be automatically set to biasRange.
  2. unbiasFWHM(None, number): The bias gaussian Full Width at Half Maximum. It must be smaller than half of biasRange. If None is given, it will be automatically set to biasFWHM.
  3. unbiasHeight(number): The unbias gaussian maximum intensity. If None is given, it will be automatically set to biasHeight.
  4. unbiasThreshold(number): unbias is only applied at a certain position only when the position weight is above unbiasThreshold. It must be a positive number.
bias_scheme_by_index(index, scaleFactor=None, check=True)

Bias the generator’s scheme using the defined bias gaussian function at the given index.

Parameters:
  1. index(integer): The index of the position to bias
  2. scaleFactor(None, number): Whether to scale the bias gaussian before biasing the scheme. If None is given, bias gaussian is used as defined.
  3. check(boolean): Whether to check arguments.
bias_scheme_at_position(position, scaleFactor=None, check=True)

Bias the generator’s scheme using the defined bias gaussian function at the given number.

Parameters:
  1. position(number): The number to bias.
  2. scaleFactor(None, number): Whether to scale the bias gaussian before biasing the scheme. If None is given, bias gaussian is used as defined.
  3. check(boolean): Whether to check arguments.
unbias_scheme_by_index(index, scaleFactor=None, check=True)

Unbias the generator’s scheme using the defined bias gaussian function at the given index.

Parameters:
  1. index(integer): The index of the position to unbias.
  2. scaleFactor(None, number): Whether to scale the unbias gaussian before unbiasing the scheme. If None is given, unbias gaussian is used as defined.
  3. check(boolean): Whether to check arguments.
unbias_scheme_at_position(position, scaleFactor=None, check=True)

Unbias the generator’s scheme using the defined bias gaussian function at the given number.

Parameters:
  1. position(number): The number to unbias.
  2. scaleFactor(None, number): Whether to scale the unbias gaussian before unbiasing the scheme. If None is given, unbias gaussian is used as defined.
  3. check(boolean): Whether to check arguments.
generate()

Generate a random float number between the biased range lowerLimit and upperLimit.

class fullrmc.Core.Collection.RandomIntegerGenerator(lowerLimit, upperLimit)

Bases: object

Generate random integer number between a lower and an upper limit.

Parameters:
  1. lowerLimit (number): Lower limit allowed.
  2. upperLimit (number): Upper limit allowed.
lowerLimit

Lower limit of the number generation.

upperLimit

Upper limit of the number generation.

rang

The range defined as upperLimit-lowerLimit

set_lower_limit(lowerLimit)

Set lower limit.

Parameters:
  1. lowerLimit (number): The lower limit allowed.
set_upper_limit(upperLimit)

Set upper limit.

Parameters:
  1. upperLimit (number): The upper limit allowed.
generate()

Generate a random integer number between lowerLimit and upperLimit.

class fullrmc.Core.Collection.BiasedRandomIntegerGenerator(lowerLimit, upperLimit, weights=None, biasHeight=1, unbiasHeight=None, unbiasThreshold=1)

Bases: fullrmc.Core.Collection.RandomIntegerGenerator

Generate biased random integer number between a lower and an upper limit. To bias the generator at a certain number, a bias height is added to the weights scheme at the position of this particular number.

_images/biasedIntegerGenerator.png
Parameters:
  1. lowerLimit (integer): The lower limit allowed.
  2. upperLimit (integer): The upper limit allowed.
  3. weights (None, list, numpy.ndarray): The weights scheme. The length must be equal to the range between lowerLimit and upperLimit. If None is given, ones array of length upperLimit-lowerLimit+1 is automatically generated.
  4. biasHeight(number): The weight bias intensity.
  5. unbiasHeight(None, number): The weight unbias intensity. If None, it will be automatically set to biasHeight.
  6. unbiasThreshold(number): unbias is only applied at a certain position only when the position weight is above unbiasThreshold. It must be a positive number.
originalWeights

Original weights as initialized.

weights

Current value weights vector.

scheme

Numbers generation scheme.

bins

Number of bins that is equal to the length of weights vector.

set_weights(weights)

Set the generator integer numbers weights.

Parameters:
  1. weights (None, list, numpy.ndarray): The weights scheme. The length must be equal to the range between lowerLimit and upperLimit. If None is given, ones array of length upperLimit-lowerLimit+1 is automatically generated.
set_bias_height(biasHeight)

Set weight bias intensity.

Parameters:
  1. biasHeight(number): Weight bias intensity.
set_unbias_height(unbiasHeight)

Set weight unbias intensity.

Parameters:
  1. unbiasHeight(None, number): The weight unbias intensity. If None, it will be automatically set to biasHeight.
set_unbias_threshold(unbiasThreshold)

Set weight unbias threshold.

Parameters:
  1. unbiasThreshold(number): unbias is only applied at a certain position only when the position weight is above unbiasThreshold. It must be a positive number.
bias_scheme_by_index(index, scaleFactor=None, check=True)

Bias the generator’s scheme at the given index.

Parameters:
  1. index(integer): The index of the position to bias
  2. scaleFactor(None, number): Whether to scale the bias gaussian before biasing the scheme. If None, bias gaussian is used as defined.
  3. check(boolean): Whether to check arguments.
bias_scheme_at_position(position, scaleFactor=None, check=True)

Bias the generator’s scheme at the given number.

Parameters:
  1. position(number): The number to bias.
  2. scaleFactor(None, number): Whether to scale the bias gaussian before biasing the scheme. If None, bias gaussian is used as defined.
  3. check(boolean): Whether to check arguments.
unbias_scheme_by_index(index, scaleFactor=None, check=True)

Unbias the generator’s scheme at the given index.

Parameters:
  1. index(integer): The index of the position to unbias
  2. scaleFactor(None, number): Whether to scale the unbias gaussian before unbiasing the scheme. If None is given, unbias gaussian is used as defined.
  3. check(boolean): Whether to check arguments.
unbias_scheme_at_position(position, scaleFactor=None, check=True)

Unbias the generator’s scheme using the defined bias gaussian function at the given number.

Parameters:
  1. position(number): The number to unbias.
  2. scaleFactor(None, number): Whether to scale the unbias gaussian before unbiasing the scheme. If None is given, unbias gaussian is used as defined.
  3. check(boolean): Whether to check arguments.
generate()

Generate a random intger number between the biased range lowerLimit and upperLimit.

fullrmc.Core.Collection.generate_random_float()

random() -> x in the interval [0, 1).

Constraint

Constraint contains parent classes for all constraints. A Constraint is used to set certain rules for the stochastic engine to evolve the atomic system. Therefore it has become possible to fully customize and set any possibly imaginable rule.

Inheritance diagram of fullrmc.Core.Constraint
class fullrmc.Core.Constraint.Constraint

Bases: fullrmc.Core.Collection.ListenerBase

A constraint is used to direct the evolution of the atomic configuration towards the desired and most meaningful one.

constraintId

Constraint unique ID create at instanciation time.

constraintName

Constraints unique name in engine given when added to engine.

engine

Stochastic fullrmc’s engine instance.

usedFrame

Get used frame in engine. If None then engine is not defined yet

computationCost

Computation cost number.

state

Constraint’s state.

tried

Constraint’s number of tried moves.

accepted

Constraint’s number of accepted moves.

used

Constraint’s used flag. Defines whether constraint is used in the stochastic engine at runtime or set inactive.

varianceSquared

Constraint’s varianceSquared used in the stochastic engine at runtime to calculate the total constraint’s standard error.

standardError

Constraint’s standard error value.

originalData

Constraint’s original data calculated upon initialization.

data

Constraint’s current calculated data.

activeAtomsDataBeforeMove

Constraint’s current calculated data before last move.

activeAtomsDataAfterMove

Constraint’s current calculated data after last move.

afterMoveStandardError

Constraint’s current calculated StandardError after last move.

amputationData

Constraint’s current calculated data after amputation.

amputationStandardError

Constraint’s current calculated StandardError after amputation.

multiframeWeight

Get constraint weight towards total in a multiframe system.

is_in_engine(engine)

Get whether constraint is already in defined and added to engine. It can be the same exact instance or a repository pulled instance of the same constraintId

Parameters:
  1. engine (stochastic fullrmc engine): Engine instance.
Returns:
  1. result (boolean): Whether constraint exists in engine.
set_variance_squared(value, frame=None)

Set constraint’s variance squared that is used in the computation of the total stochastic engine standard error.

Parameters:
  1. value (number): Any positive non zero number.
  2. frame (None, string): Target frame name. If None, engine used frame is used. If multiframe is given, all subframes will be targeted. If subframe is given, all other multiframe subframes will be targeted.
set_computation_cost(value, frame=None)

Set constraint’s computation cost value. This is used at stochastic engine runtime to minimize computations and enhance performance by computing less costly constraints first. At every step, constraints will be computed in order starting from the less to the most computationally costly. Therefore upon rejection of a step because of an unsatisfactory rigid constraint, the left un-computed constraints at this step are guaranteed to be the most time coslty ones.

Parameters:
  1. value (number): computation cost.
  2. frame (None, string): Target frame name. If None, engine used frame is used. If multiframe is given, all subframes will be targeted. If subframe is given, all other multiframe subframes will be targeted.
set_used(*args, **kwargs)

Set used flag.

Parameters:
  1. value (boolean): True to use this constraint in stochastic engine runtime.
  2. frame (None, string): Target frame name. If None, engine used frame is used. If multiframe is given, all subframes will be targeted. If subframe is given, all other multiframe subframes will be targeted.
set_state(value)

Set constraint’s state. When constraint’s state and stochastic engine’s state don’t match, constraint’s data must be re-calculated.

Parameters:
  1. value (object): Constraint state value.
set_tried(value)

Set constraint’s number of tried moves.

Parameters:
  1. value (integer): Constraint tried moves value.
increment_tried()

Increment number of tried moves.

set_accepted(value)

Set constraint’s number of accepted moves.

Parameters:
  1. value (integer): Constraint’s number of accepted moves.
increment_accepted()

Increment constraint’s number of accepted moves.

set_standard_error(value)

Set constraint’s standardError value.

Parameters:
  1. value (number): standard error value.
set_data(value)

Set constraint’s data value.

Parameters:
  1. value (number): constraint’s data.
set_active_atoms_data_before_move(value)

Set constraint’s before move happens active atoms data value.

Parameters:
  1. value (number): Data value.
set_active_atoms_data_after_move(value)

Set constraint’s after move happens active atoms data value.

Parameters:
  1. value (number): data value.
set_after_move_standard_error(value)

Set constraint’s standard error value after move happens.

Parameters:
  1. value (number): standard error value.
set_amputation_data(value)

Set constraint’s after amputation data.

Parameters:
  1. value (number): data value.
set_amputation_standard_error(value)

Set constraint’s standardError after amputation.

Parameters:
  1. value (number): standard error value.
reset_constraint(reinitialize=True, flags=False, data=False, frame=None)

Reset constraint.

Parameters:
  1. reinitialize (boolean): If set to True, it will override the rest of the flags and will completely reinitialize the constraint.
  2. flags (boolean): Reset the state, tried and accepted flags of the constraint.
  3. data (boolean): Reset the constraints computed data.
  4. frame (None, string): Target frame name. If None, engine used frame is used. If multiframe is given, all subframes will be targeted. If subframe is given, rest of multiframe subframes will not be targeted.
update_standard_error()

Compute and set constraint’s standard error by calling compute_standard_error method and passing constraint’s data.

get_constraints_properties(frame)

Get a dictionary look up table of constraint’s properties

Parameters:
  1. frame (string): frame to pull and build contraint data. It can be a traditional frame, a multiframe or a subframe
Returns:
  1. propertiesLUT (dictionary): properties value look up table. Keys are described herein. All keys start with ‘frames-‘ and values are list of properties for every and each frame.

    • frames-name: list of all frames name
    • frames-weight: list of all frames weight
    • frames-number_of_removed_atoms: list of number of removed atoms from each frame
    • frames-constraint: list of constraint copy
    • frames-data: list of constraint data
    • frames-standard_error: list of all frames standard error
get_constraint_value()

Method must be overloaded in children classes.

get_constraint_original_value()

Method must be overloaded in children classes.

compute_standard_error()

Method must be overloaded in children classes.

compute_data(update=True)

Method must be overloaded in children classes.

compute_before_move(realIndexes, relativeIndexes)

Method must be overloaded in children classes.

compute_after_move(realIndexes, relativeIndexes, movedBoxCoordinates)

Method must be overloaded in children classes.

accept_move(realIndexes, relativeIndexes)

Method must be overloaded in children classes.

reject_move(realIndexes, relativeIndexes)

Method must be overloaded in children classes.

compute_as_if_amputated(realIndex, relativeIndex)

Method must be overloaded in children classes.

compute_as_if_inserted(realIndex, relativeIndex)

Method must be overloaded in children classes.

accept_amputation(realIndex, relativeIndex)

Method must be overloaded in children classes.

reject_amputation(realIndex, relativeIndex)

Method must be overloaded in children classes.

accept_insertion(realIndex, relativeIndex)

Method must be overloaded in children classes.

reject_insertion(realIndex, relativeIndex)

Method must be overloaded in children classes.

export(fileName, frame=None, format='%s', delimiter='\t', comments='#')

Export constraint data to text file or to an archive of files.

Parameters:
  1. fileName (path): full file name and path.
  2. frame (None, string): frame name to export data from. If multiframe is given, multiple files will be created with subframe name appended to the end.
  3. format (string): string format to export the data. format is as follows (%[flag]width[.precision]specifier)
  4. delimiter (string): String or character separating columns.
  5. comments (string): String that will be prepended to the header.
plot(frame=None, axes=None, subAdParams={'bottom': None, 'hspace': 0.4, 'left': None, 'right': None, 'top': None, 'wspace': None}, dataParams={'label': 'Y', 'linewidth': 2}, xlabelParams={'size': 10, 'xlabel': 'Core-Shell atoms'}, ylabelParams={'size': 10, 'ylabel': 'Coordination number'}, xticksParams={'fontsize': 8, 'rotation': 90}, yticksParams={'fontsize': 8, 'rotation': 0}, legendParams={'fontsize': 8, 'frameon': False, 'loc': 'upper right', 'ncol': 1}, titleParams={'fontsize': 10, 'label': '@{frame} (${numberOfRemovedAtoms:.1f}$ $rem.$ $at.$) $Std.Err.={standardError:.3f}$ - $used$ $({used})$'}, gridParams=None, show=True, **paramsKwargs)

Plot constraint data. This can be overloaded in children classes.

Parameters:
  1. frame (None, string): The frame name to plot. If None, used frame will be plotted.
  2. ax (None, matplotlib Axes): matplotlib Axes instance to plot in. If None is given a new plot figure will be created.
  3. subAdParams (None, dict): matplotlib.artist.Artist.subplots_adjust parameters subplots adjust parameters.
  4. dataParams (None, dict): constraint data plotting parameters
  5. xlabelParams (None, dict): matplotlib.axes.Axes.set_xlabel parameters.
  6. ylabelParams (None, dict): matplotlib.axes.Axes.set_ylabel parameters.
  7. legendParams (None, dict):matplotlib.axes.Axes.legend parameters.
  8. xticksParams (None, dict):matplotlib.axes.Axes.set_xticklabels parameters.
  9. yticksParams (None, dict):matplotlib.axes.Axes.set_yticklabels parameters.
  10. titleParams (None, dict): matplotlib.axes.Axes.set_title parameters
  11. gridParams (None, dict): matplotlib.axes.Axes.grid parameters
  12. show (boolean): Whether to render and show figure before returning.
Returns:
  1. figure (matplotlib Figure): matplotlib used figure.
  2. axes (matplotlib Axes): matplotlib used axes.
class fullrmc.Core.Constraint.ExperimentalConstraint(experimentalData, dataWeights=None, scaleFactor=1.0, adjustScaleFactor=(0, 0.8, 1.2))

Bases: fullrmc.Core.Constraint.Constraint

Experimental constraint is any constraint related to experimental data.

Parameters:
  1. engine (None, fullrmc.Engine): Constraint’s stochastic engine.

  2. experimentalData (numpy.ndarray, string): Experimental data goiven as numpy.ndarray or string path to load data using numpy.loadtxt method.

  3. dataWeights (None, numpy.ndarray): Weights array of the same number of points of experimentalData used in the constraint’s standard error computation. Therefore particular fitting emphasis can be put on different data points that might be considered as more or less important in order to get a reasonable and plausible modal.

    If None is given, all data points are considered of the same importance in the computation of the constraint’s standard error.

    If numpy.ndarray is given, all weights must be positive and all zeros weighted data points won’t contribute to the total constraint’s standard error. At least a single weight point is required to be non-zeros and the weights array will be automatically scaled upon setting such as the the sum of all the weights is equal to the number of data points.

  4. scaleFactor (number): A normalization scale factor used to normalize the computed data to the experimental ones.

  5. adjustScaleFactor (list, tuple): Used to adjust fit or guess the best scale factor during stochastic engine runtime. It must be a list of exactly three entries.

    1. The frequency in number of generated moves of finding the best scale factor. If 0 frequency is given, it means that the scale factor is fixed.
    2. The minimum allowed scale factor value.
    3. The maximum allowed scale factor value.

NB: If adjustScaleFactor first item (frequency) is 0, the scale factor will remain untouched and the limits minimum and maximum won’t be checked.

experimentalData

Experimental data of the constraint.

dataWeights

Experimental data points weight

multiframeWeight

Get constraint weight towards total in a multiframe system.

multiframePrior

Get constraint multiframe prior array.

scaleFactor

Constraint’s scaleFactor.

adjustScaleFactor

Adjust scale factor tuple.

adjustScaleFactorFrequency

Scale factor adjustment frequency.

adjustScaleFactorMinimum

Scale factor adjustment minimum number allowed.

adjustScaleFactorMaximum

Scale factor adjustment maximum number allowed.

limits

Used data X limits.

limitsIndexStart

Used data start index as calculated from limits.

limitsIndexEnd

Used data end index as calculated from limits.

set_scale_factor(scaleFactor)

Set the scale factor. This method doesn’t allow specifying frames. It will target used frame only.

Parameters:
  1. scaleFactor (number): A normalization scale factor used to normalize the computed data to the experimental ones.
set_adjust_scale_factor(adjustScaleFactor)

Set adjust scale factor. This method doesn’t allow specifying frames. It will target used frame only.

Parameters:
  1. adjustScaleFactor (list, tuple): Used to adjust fit or guess the best scale factor during stochastic engine runtime. It must be a list of exactly three entries.
    1. The frequency in number of generated moves of finding the best scale factor. If 0 frequency is given, it means that the scale factor is fixed.
    2. The minimum allowed scale factor value.
    3. The maximum allowed scale factor value.
set_experimental_data(experimentalData)

Set the constraint’s experimental data. This method will raise an error if called after adding constraint to stochastic engine.

Parameters:
  1. experimentalData (numpy.ndarray, string, list, tuple): Experimental data as numpy.ndarray or string path to load data using numpy.loadtxt method. If list or tuple are given, they will be automatically converted to a numpy array by calling numpy.array(experimentalData). Finally experimental data type will be converted to fullrmc.Globals.FLOAT_TYPE
set_data_weights(dataWeights, frame=None)

Set experimental data points weight. Data weights will be automatically normalized.

Parameters:
  1. dataWeights (None, numpy.ndarray): Weights array of the same number of points of experimentalData used in the constraint’s standard error computation. Therefore particular fitting emphasis can be put on different data points that might be considered as more or less important in order to get a reasonable and plausible model.

    If None is given, all data points are considered of the same importance in the computation of the constraint’s standard error.

    If numpy.ndarray is given, all weights must be positive and all zeros weighted data points won’t contribute to the total constraint’s standard error. At least a single weight point is required to be non-zeros and the weights array will be automatically scaled upon setting such as the the sum of all the weights is equal to the number of data points.

  2. frame (None, string): Target frame name. If None, engine used frame is used. If multiframe is given, all subframes will be targeted. If subframe is given, rest of multiframe subframes will not be targeted.

check_experimental_data(experimentalData)

Checks the constraint’s experimental data This method must be overloaded in all experimental constraint sub-classes.

Parameters:
  1. experimentalData (numpy.ndarray): Experimental data numpy.ndarray.
fit_scale_factor(experimentalData, modelData, dataWeights)

The best scale factor value is computed by minimizing \(E=sM\).

Where:
  1. \(E\) is the experimental data.
  2. \(s\) is the scale factor.
  3. \(M\) is the model constraint data.

This method doesn’t allow specifying frames. It will target used frame only.

Parameters:
  1. experimentalData (numpy.ndarray): Experimental data.
  2. modelData (numpy.ndarray): Constraint modal data.
  3. dataWeights (None, numpy.ndarray): Data points weights to compute the scale factor. If None is given, all data points will be considered as having the same weight.
Returns:
  1. scaleFactor (number): The new scale factor fit value.

NB: This method won’t update the internal scale factor value of the constraint. It always computes the best scale factor given experimental and atomic model data.

get_adjusted_scale_factor(experimentalData, modelData, dataWeights)

Checks if scale factor should be updated according to the given scale factor frequency and engine’s accepted steps. If adjustment is due, a new scale factor will be computed using fit_scale_factor method, otherwise the the constraint’s scale factor will be returned.

Parameters:
  1. experimentalData (numpy.ndarray): the experimental data.
  2. modelData (numpy.ndarray): the constraint modal data.
  3. dataWeights (None, numpy.ndarray): the data points weights to compute the scale factor. If None is given, all data points will be considered as having the same weight.
Returns:

#. scaleFactor (number): Constraint’s scale factor or the new scale factor fit value.

NB: This method WILL NOT UPDATE the internal scale factor value of the constraint.

compute_standard_error(experimentalData, modelData)

Compute the squared deviation between modal computed data and the experimental ones.

\[SD = \sum \limits_{i}^{N} W_{i}(Y(X_{i})-F(X_{i}))^{2}\]

Where:

\(N\) is the total number of experimental data points.

\(W_{i}\) is the data point weight. It becomes equivalent to 1 when dataWeights is set to None.

\(Y(X_{i})\) is the experimental data point \(X_{i}\).

\(F(X_{i})\) is the computed from the model data \(X_{i}\).

Parameters:
  1. experimentalData (numpy.ndarray): Experimental data.
  2. modelData (numpy.ndarray): The data to compare with the experimental one and compute the squared deviation.
Returns:
  1. standardError (number): The calculated standard error of the constraint.
get_constraints_properties(frame)

Get a dictionary look up table of constraint’s properties that are needed to plot or export

Parameters:
  1. frame (string): frame to pull and build contraint data. It can be a traditional frame, a multiframe or a subframe
Returns:
  1. propertiesLUT (dictionary): properties value look up table. Keys are described herein. Values of keys that start with ‘frames-‘ are a list for all frames. Values of keys that start with ‘weighted-‘ are weighted values for all frames

    • frames-name: list of all frames name.
    • frames-weight: list of all frames weight.
    • frames-number_of_removed_atoms: list of number of removed atoms from each frame.
    • frames-experimental_x: list of numpy array of experimental x data.
    • frames-experimental_y: list of numpy array of experimental y data.
    • frames-output: list of frames dictionary constraint output data
    • frames-model_x: list of numpy array of model x data.
    • frames-shape_array: list of system shape function (numpy array) of all frames.
    • frames-window_function: list of window function (numpy array) of all frames.
    • frames-scale_factor: list of all frames scale factor.
    • frames-standard_error: list of all frames standard error.
    • weighted-output: dictionary of all frames weighted constraint data using ‘frames-weight’
    • weighted-number_of_removed_atoms: All frames averaged number of removed atoms using ‘frames-weight’
    • weighted-scale_factor: All frames averaged scale factor using ‘frames-weight’
    • weighted-standard_error: All frames weighted standard error using ‘frames-weight’
plot(frame=None, axes=None, asMesoscopic=False, intra=True, inter=True, shapeFunc=True, subAdParams={'bottom': None, 'hspace': 0.4, 'left': None, 'right': None, 'top': None, 'wspace': None}, totParams={'color': 'black', 'label': 'total', 'linewidth': 3.0, 'zorder': 1}, expParams={'color': 'red', 'label': 'experimental', 'marker': 'o', 'markersize': 7.5, 'markevery': 1, 'zorder': 0}, noWParams={'color': 'black', 'label': 'total - no window', 'linewidth': 1.0, 'zorder': 1}, shaParams={'color': 'black', 'label': 'shape function', 'linestyle': 'dashed', 'linewidth': 1.0, 'zorder': 2}, parParams={'linewidth': 1.0, 'markersize': 5, 'markevery': 5, 'zorder': 3}, xlabelParams={'size': 10, 'xlabel': 'X'}, ylabelParams={'size': 10, 'ylabel': 'Y'}, xticksParams={'fontsize': 8, 'rotation': 0}, yticksParams={'fontsize': 8, 'rotation': 0}, legendParams={'fontsize': 8, 'frameon': False, 'loc': 'upper right', 'ncol': 2}, titleParams={'fontsize': 10, 'label': '@{frame} (${numberOfRemovedAtoms:.1f}$ $rem.$ $at.$) $Std.Err.={standardError:.3f}$\n$scale$ $factor$=${scaleFactor:.2f}$ - $multiframe$ $weight$=${multiframeWeight:.3f}$ - $used$ $({used})$'}, gridParams=None, show=True, **paramsKwargs)

Plot constraint data. This can be overloaded in children classes.

Parameters:
  1. frame (None, string): The frame name to plot. If None, used frame will be plotted.
  2. axes (None, matplotlib Axes): matplotlib Axes instance to plot in. If None is given a new plot figure will be created.
  3. asMesoscopic (boolean): If given frame is a multiframe is, when true, asMesoscopic considers all frames as a statistical average of all frames in a mesoscopic system. All subframes will be then plotted as a single weighted mesoscopic structure. If asMesocopic is False and given frame is a multiframe given ax will be disregarded
  4. intra (boolean): Whether to add intra-molecular pair distribution function features to the plot.
  5. inter (boolean): Whether to add inter-molecular pair distribution function features to the plot.
  6. shapeFunc (boolean): Whether to add shape function to the plot only when exists.
  7. subAdParams (None, dict): matplotlib.artist.Artist.subplots_adjust parameters subplots adjust parameters.
  8. totParams (None, dict): constraint total plotting parameters
  9. expParams (None, dict): constraint experimental data parameters
  10. noWParams (None, dict): constraint total without window parameters
  11. shaParams (None, dict): constraint shape function parameters
  12. parParams (None, dict): constraint partials parameters
  13. xlabelParams (None, dict): matplotlib.axes.Axes.set_xlabel parameters.
  14. ylabelParams (None, dict): matplotlib.axes.Axes.set_ylabel parameters.
  15. legendParams (None, dict):matplotlib.axes.Axes.legend parameters.
  16. xticksParams (None, dict):matplotlib.axes.Axes.set_xticklabels parameters.
  17. yticksParams (None, dict):matplotlib.axes.Axes.set_yticklabels parameters.
  18. titleParams (None, dict): matplotlib.axes.Axes.set_title parameters
  19. gridParams (None, dict): matplotlib.axes.Axes.grid parameters
  20. show (boolean): Whether to render and show figure before returning.
Returns:
  1. figure (matplotlib Figure): matplotlib used figure.
  2. axes (matplotlib Axes): matplotlib used axes.
plot_multiframe_weights(frame, ax=None, titleFormat='@{frame} [subframes probability distribution]', show=True)

plot multiframe subframes weight distribution histogram

Parameters:
  1. frame (string): multiframe name
  2. ax (None, matplotlib Axes): matplotlib Axes instance to plot in. If None is given a new plot figure will be created.
  3. titleFormat (string): title format. If empty string is given no title will be added to figure axes
  4. show (boolean): Whether to render and show figure before
    returning.
Returns:
  1. figure (matplotlib Figure): matplotlib used figure.
  2. axes (matplotlib Axes): matplotlib used axes.
export(fileName, frame=None, format='%12.5f', delimiter='\t', comments='#')

Export constraint data to text file or to an archive of files.

Parameters:
  1. fileName (path): full file name and path.
  2. frame (None, string): frame name to export data from. If multiframe is given, multiple files will be created with subframe name appended to the end.
  3. format (string): string format to export the data. format is as follows (%[flag]width[.precision]specifier)
  4. delimiter (string): String or character separating columns.
  5. comments (string): String that will be prepended to the header.
class fullrmc.Core.Constraint.SingularConstraint

Bases: fullrmc.Core.Constraint.Constraint

A singular constraint is a constraint that doesn’t allow multiple instances in the same engine.

is_singular(engine)

Get whether only one instance of this constraint type is present in the stochastic engine. True for only itself found, False for other instance of the same __class__.__name__ or constraintId.

Parameters:
  1. engine (stochastic fullrmc engine): Engine instance.
Returns:
  1. result (boolean): Whether constraint is singular in engine.
assert_singular(engine)

Checks whether only one instance of this constraint type is present in the stochastic engine. Raises Exception if multiple instances are present.

class fullrmc.Core.Constraint.RigidConstraint(rejectProbability)

Bases: fullrmc.Core.Constraint.Constraint

A rigid constraint is a constraint that doesn’t count into the total standard error of the stochastic Engine. But it’s internal standard error must monotonously decrease or remain the same from one engine step to another. If standard error of an rigid constraint increases the step will be rejected even before engine’s new standardError get computed.

Parameters:
  1. rejectProbability (Number): Rejecting probability of all steps where standard error increases. It must be between 0 and 1 where 1 means rejecting all steps where standardError increases and 0 means accepting all steps regardless whether standard error increases or not.
rejectProbability

Rejection probability.

set_reject_probability(rejectProbability)

Set the rejection probability. This method doesn’t allow specifying frames. It will target used frame only.

Parameters:
  1. rejectProbability (Number): rejecting probability of all steps where standard error increases. It must be between 0 and 1 where 1 means rejecting all steps where standardError increases and 0 means accepting all steps regardless whether standard error increases or not.
should_step_get_rejected(standardError)

Given a standard error, return whether to keep or reject new standard error according to the constraint reject probability.

Parameters:

#. standardError (number): The standard error to compare with the Constraint standard error

Return:
  1. result (boolean): True to reject step, False to accept
should_step_get_accepted(standardError)

Given a standard error, return whether to keep or reject new standard error according to the constraint reject probability.

Parameters:
  1. standardError (number): The standard error to compare with the Constraint standard error
Return:
  1. result (boolean): True to accept step, False to reject
fullrmc.Core.Constraint.randfloat()

random() -> x in the interval [0, 1).

Group

Group contains parent classes for all groups. A Group is a set of atoms indexes used to gather atoms and apply actions such as moves upon them. Therefore it has become possible to fully customize and separate atoms to groups and perform stochastic actions on groups rather than on single atoms.

Inheritance diagram of fullrmc.Core.Group
class fullrmc.Core.Group.Group(indexes, moveGenerator=None, refine=False, name='')

Bases: object

A Group is a set of atoms indexes container.

Parameters:
  1. indexes (np.ndarray, list, set, tuple): list of atoms indexes.
  2. moveGenerator (None, MoveGenerator): Move generator instance. If None is given AtomsRemoveGenerator is considered by default.
  3. refine (bool): The refinement flag used by the Engine.
  4. name (str): The group user defined name.
# import fullrmc modules
from fullrmc.Engine import Engine
from fullrmc.Core.Group import Group

# create engine
ENGINE = Engine(path='my_engine.rmc')

# set pdb file
ENGINE.set_pdb('system.pdb')

# Add constraints ...

# re-define groups as atoms.
groups = [Group([idx]) for idx in ENGINE.pdb.indexes]
ENGINE.set_groups( groups )

# Re-define groups generators as needed ... By default AtomsRemoveGenerator is used.
indexes

Atoms index array.

moveGenerator

Group’s move generator instance.

refine

Refine flag.

name

groud user defined name.

set_refine(refine)

Set the selector refine flag.

Parameters:
  1. refine (bool): The selector refinement flag.
set_name(name)

Set the group’s name.

Parameters:
  1. name (str): The group user defined name.
set_indexes(indexes)

Set group atoms index. Indexes redundancy is not checked and indexes order is preserved.

Parameters:
  1. indexes (list,set,tuple,np.ndarray): The group atoms indexes.
set_move_generator(generator)

Set group move generator.

Parameters:
  1. generator (None, MoveGenerator): Move generator instance. If None is given TranslationGenerator is considered by default.
class fullrmc.Core.Group.EmptyGroup(*args, **kwargs)

Bases: fullrmc.Core.Group.Group

Empty group is a special group that takes no atoms indexes. It’s mainly used to remove atoms from system upon fitting.

# import fullrmc modules
from fullrmc.Engine import Engine
from fullrmc.Core.Group import EmptyGroup

# create engine
ENGINE = Engine(path='my_engine.rmc')

# set pdb file
ENGINE.set_pdb('system.pdb')

# Add constraints ...

# re-define groups and set a single empty group
ENGINE.set_groups( EmptyGroup() )

# Re-define groups generators as needed ... By default RemoveGenerator is used.
moveGenerator

Group’s move generator instance.

indexes

Always returns None for EmptyGroup

set_move_generator(generator)

Set group move generator.

Parameters:
  1. generator (None, MoveGenerator): Move generator instance. If None is given TranslationGenerator is considered by default.
set_indexes(indexes)

Sets the group indexes. For an EmptyGroup, this method will disregard given indexes argument and will always set indexes property to None.

Parameters:
  1. indexes (object): The group atoms indexes. This argument will always be disregarded in this particular case.

GroupSelector

GroupSelector contains parent classes for all group selectors. A GroupSelector is used at the stochastic engine’s runtime to select groups upon which a move will be applied. Therefore it has become possible to fully customize the selection of groups of atoms and to choose when and how frequently a group can be chosen to perform a move upon.

Inheritance diagram of fullrmc.Core.GroupSelector
class fullrmc.Core.GroupSelector.GroupSelector(engine=None)

Bases: object

Group selector is the parent class that selects groups to perform moves at stochastic engine’s runtime.

Parameters:
  1. engine (None, fullrmc.Engine): Selector’s stochastic engine instance.
engine

Stochastic engine’s instance.

refine

Get refine flag value. It will always return False because refine is a property of RecursiveGroupSelector instances only.

explore

Get explore flag value. It will always return False because explore is a property of RecursiveGroupSelector instances only.

willSelect

Get whether next step a new selection is occur or still the same group is going to be selected again. It will always return True because recurrence is a property of RecursiveGroupSelector instances only.

willRecur

Get whether next step the same group will be returned. It will always return False because this is a property of RecursiveGroupSelector instances only.

willRefine

Get whether selection is recurring and refine flag is True. It will always return False because recurrence is a property of RecursiveGroupSelector instances only.

willExplore

Get whether selection is recurring and explore flag is True. It will always return False because recurrence is a property of RecursiveGroupSelector instances only.

isNewSelection

Get whether the last step a new selection was made. It will always return True because recurrence is a property of RecursiveGroupSelector instances only.

isRecurring

Get whether the last step the same group was returned. It will always return False because this is a property of RecursiveGroupSelector instances only.

isRefining

Get whether selection is recurring and refine flag is True. It will always return False because recurrence is a property of RecursiveGroupSelector instances only.

isExploring

Get whether selection is recurring and explore flag is True. It will always return False because recurrence is a property of RecursiveGroupSelector instances only.

set_engine(engine)

Set selector’s stochastic engine instance.

Parameters:
  1. engine (None, fullrmc.Engine): Selector’s stochastic engine.
select_index()

This method must be overloaded in every GroupSelector sub-class

Returns:
  1. index (integer): the selected group index in engine groups list.
move_accepted(index)

This method is called by the stochastic engine when a move generated on a group is accepted. This method is empty must be overloaded when needed.

Parameters:
  1. index (integer): the selected group index in engine groups list.
move_rejected(index)

This method is called by the stochastic engine when a move generated on a group is rejected. This method is empty must be overloaded when needed.

Parameters:
  1. index (integer): the selected group index in engine groups list.
class fullrmc.Core.GroupSelector.RecursiveGroupSelector(selector, recur=10, override=True, refine=False, explore=False)

Bases: fullrmc.Core.GroupSelector.GroupSelector

Recursive selector is the only selector that can use the recursive property on a selection. It is used as a wrapper around a GroupSelector instance.

Parameters:
  1. selector (fullrmc.Core.GroupSelector.GroupSelector): The selector instance to wrap.
  2. recur (integer): Set number of times to recur. It must be a positive integer.
  3. override (boolean): Override temporary recur value. recur value will be overridden only when selected group move generator is a PathGenerator instance. In this particular case, recur value will be temporary changed to the number of moves stored in the PathGenerator. If selected group move generator is not a PathGenerator instance, recur value will take back its original value.
  4. refine (boolean): Its an engine flag that is used to refine the position of a group until recurrence expires and a new group is selected. Refinement is done by applying moves upon the selected group always from its initial position at the time it was selected until recurrence expires, then the best position is kept.
  5. explore (boolean): Its an engine flag that is used to make a group explore the space around it until recurrence expires and a new group is selected. Exploring is done by applying moves upon the selected group starting from its initial position and evolving in a trajectory like way until recurrence expires, then the best position is kept.

NB: refine and explore flags can’t both be set to True at the same time. When this happens refine flag gets automatically switched to False. The usage of those flags is very important because they allow groups of atoms to go out of local minima in the energy surface. The way traditional reverse mote carlo works is by minimizing the total energy of the system (error) using gradient descent method. Using of those flags allows the system to go up hill in the energy surface searching for other lower minimas, while always conserving the lowest energy state found and not changing the system structure until a better structure with smaller error is found.

The following video compares the Reverse Monte Carlo traditional fitting mode with fullrmc's recursive selection one with explore flag set to True. From a potential point of view, exploring allows to cross forbidden unlikely energy barriers and going out of local minimas.

The following video is an example of refining the position of a molecule using RecursiveGroupSelector and setting refine flag to True. The molecule is always refined from its original position towards a new one generated by the move generator.

The following video is an example of exploring the space of a molecule using RecursiveGroupSelector and setting explore flag to True. The molecule explores the allowed space by wandering via its move generator and only moves enhancing the structure are stored.

# import fullrmc modules
from fullrmc.Engine import Engine
from fullrmc.Core.GroupSelector import RecursiveGroupSelector

# create engine
ENGINE = Engine(path='my_engine.rmc')

# set pdb file
ENGINE.set_pdb('system.pdb')

# Add constraints ...
# Re-define groups if needed ...
# Re-define groups selector if needed ...

##### Wrap engine group selector with a recursive group selector. #####
# create recursive group selector. Recurrence is set to 20 with explore flag set to True.
RGS = RecursiveGroupSelector(ENGINE.groupSelector, recur=20, refine=False, explore=True)
ENGINE.set_group_selector(RGS)
selector

The wrapped selector instance.

lastSelectedIndex

The last selected group index.

willSelect

Get whether next step a new selection is occur or still the same group is going to be selected again.

willRecur

Get whether next step the same group will be returned.

willRefine

Get whether next step the same group will be returned and refine flag is True.

willExplore

Get whether next step the same group will be returned and explore flag is True.

isNewSelection

Get whether this last step a new selection was made.

isRecurring

Get whether this last step the same group was returned.

isRefining

Get whether this last step the same group was returned and refine flag is True.

isExploring

Get whether this last step the same group was returned and explore flag is True.

override

Override flag value.

refine

Refine flag value.

explore

Explore flag value.

currentRecur

The current recur value which is selected group dependant when override flag is True.

recur

The current recur value. The set recur value can change during engine runtime if override flag is True. To get the recur value as set by set_recur method recurAsSet must be used.

recurAsSet

Get recur value as set but set_recur method.

position

Get the position of the selector in the path.

engine

Get the wrapped selector engine instance.

set_engine(engine)

Sets the wrapped selector stochastic engine instance.

Parameters:
  1. engine (None, fullrmc.Engine): The selector stochastic engine.
set_recur(recur)

Sets the recur value.

Parameters:
  1. recur (integer): Set the recur value. It must be a positive integer.
set_override(override)

Select override value.

Parameters:
  1. override (boolean): Override selector recur value only when selected group move generator is a PathGenerator instance. Overridden recur value is temporary and totally selected group dependant. If selected group move generator is not a PathGenerator instance, recur value will take back selector’s recur value.
set_refine(refine)

Select override value.

Parameters:
  1. refine (boolean): Its an engine flag that is used to refine the position of a group until recurrence expires and a new group is selected. Refinement is done by applying moves upon the selected group always from its initial position at the time it was selected until recurrence expires, then the best position is kept.
set_explore(explore)

Select override value.

Parameters:
  1. explore (boolean): Its an engine flag that is used to make a group explore the space around it until recurrence expires and a new group is selected. Exploring is done by applying moves upon the selected group starting from its initial position and evolving in a trajectory like way until recurrence expires, then the best position is kept.
select_index()

Select new index.

Returns:
  1. index (integer): the selected group index in engine groups list.

MoveGenerator

MoveGenerator contains parent classes for all move generators. A MoveGenerator sub-class is used at fullrmc’s stochastic engine runtime to generate moves upon selected groups. Every group has its own MoveGenerator class and definitions, therefore it is possible to fully customize how a group of atoms should move.

Inheritance diagram of fullrmc.Core.MoveGenerator
class fullrmc.Core.MoveGenerator.MoveGenerator(group=None)

Bases: object

It is the parent class of all moves generators. This class can’t be instantiated but its sub-classes might be.

Parameters:
  1. group (None, Group): The group instance.
group

Group instance.

set_group(group)

Set the MoveGenerator group.

Parameters:
  1. group (None, Group): Group instance.
check_group(group)

Check the generator’s group. This method must be overloaded in all MoveGenerator sub-classes.

Parameters:
  1. group (Group): the Group instance
transform_coordinates(coordinates, argument=None)

Transform coordinates. This method is called to move atoms. This method must be overloaded in all MoveGenerator sub-classes.

Parameters:
  1. coordinates (np.ndarray): The coordinates on which to apply the move.
  2. argument (object): Any other argument needed to perform the move. In General it’s not needed.
Returns:
  1. coordinates (np.ndarray): The new coordinates after applying the move.
move(coordinates)

Moves coordinates. This method must NOT be overloaded in MoveGenerator sub-classes.

Parameters:
  1. coordinates (np.ndarray): The coordinates on which to apply the transformation.
Returns:
  1. coordinates (np.ndarray): The new coordinates after applying the transformation.
class fullrmc.Core.MoveGenerator.RemoveGenerator(group=None, maximumCollected=None, allowFittingScaleFactor=False, atomsList=None)

Bases: fullrmc.Core.MoveGenerator.MoveGenerator

This is a very particular move generator that will not generate moves on atoms but removes them from the atomic configuration using a general collector mechanism. Remove generators must be used to create defects in the simulated system. When the standard error is high, removing atoms might reduce the total fit standard error but this can be illusional and very limiting because artificial non physical voids can get created in the system which will lead to an impossibility to finding a solution at the end. It’s strongly recommended to exhaust all ideas and possibilities in finding a good solution prior to start removing atoms unless structural defects is the goal of the simulation.

All removed or amputated atoms are collected by the engine and will become available to be re-inserted in the system if needed. But keep in mind, it might be physically easy to remove and atom but an impossibility to add it back especially if the created voids are smeared out.

Removers are called generators but they behave like selectors. Instead of applying a certain move on a group of atoms, they normally pick atoms from defined atoms list and apply no moves on those. ‘move’ and ‘transform_coordinates’ methods are not implemented in this class of generators and a usage error will be raised if called. ‘pick_from_list’ method is used instead and must be overloaded by all RemoveGenerator subclasses.

N.B. This class can’t be instantiated but its sub-classes might be.

Parameters:
  1. group (None, Group): The group instance which is this case must be fullrmc EmptyGroup.
  2. maximumCollected (None, Integer): The maximum number allowed of atoms to be removed and collected from atomic configuration by the stochastic engine. This property is general to the system and checks engine’s collected atoms not the number of removed atoms via this generator. If None is given, the remover will not check for the number of already removed atoms before attempting a remove.
  3. allowFittingScaleFactor (bool): Constraints and especially experimental ones have a scale factor constant that can be fit. Fitting a scale factor happens at stochastic engine’s runtime at a certain fitting frequency. If this flag set to True, then fitting the scale factor will be allowed upon removing atoms. When set to False, fitting the constraint scale factor will be forbidden upon removing atoms. By default, allowFittingScaleFactor is set to False because it’s more logical to allow removing only atoms that enhances the total standard error without rescaling the model’s data.
  4. atomsList (None,list,set,tuple,np.ndarray): The list of atomss index to chose and remove from.
atomsList

Atoms list from which atoms will be picked to attempt removal.

allowFittingScaleFactor

Whether to allow constraints to fit their scale factor upon removing atoms.

maximumCollected

Maximum collected atoms allowed.

check_group(group)

Check the generator’s group.

Parameters:
  1. group (Group): The group instance.
set_maximum_collected(maximumCollected)

Set maximum collected number of atoms allowed.

Parameters:
  1. maximumCollected (None, Integer): The maximum number allowed of atoms to be removed and collected from atomic configuration by the stochastic engine. This property is general to the system and checks engine’s collected atoms not the number of removed atoms via this generator. If None is given, the remover will not check for the number of already removed atoms before attempting a remove.
set_allow_fitting_scale_factor(allowFittingScaleFactor)

Set allow fitting scale factor flag.

Parameters:
  1. allowFittingScaleFactor (bool): Constraints and especially experimental ones have a scale factor constant that can be fit. Fitting a scale factor happens at stochastic engine’s runtime at a certain fitting frequency. If this flag set to True, then fitting the scale factor will be allowed upon removing atoms. When set to False, fitting the constraint scale factor will be forbidden upon removing atoms. By default, allowFittingScaleFactor is set to False because it’s more logical to allow removing only atoms that enhances the total standard error without rescaling the model’s data.
set_atoms_list(atomsList)

Set atoms index list from which atoms will be picked to attempt removal. This method must be overloaded and not be called from this class but from its children. Otherwise a usage error will be raised.

Parameters:
  1. atomsList (None, list,set,tuple,np.ndarray): The list of atoms index to chose and remove from.
move(coordinates)

Moves coordinates. This method must NOT be overloaded in MoveGenerator sub-classes.

Parameters:
  1. coordinates (np.ndarray): Not used here.
transform_coordinates(coordinates, argument)

This method must NOT be overloaded in MoveGenerator sub-classes.

Parameters:
  1. coordinates (np.ndarray): Not used here. the translation.
  2. argument (object): Not used here.
pick_from_list(engine)

This method must be overloaded in all RemoveGenerator sub-classes.

Parameters:
  1. engine (Engine): stochastic engine calling the method.
class fullrmc.Core.MoveGenerator.SwapGenerator(group=None, swapLength=1, swapList=None)

Bases: fullrmc.Core.MoveGenerator.MoveGenerator

It is a particular move generator that instead of generating a move upon a group of atoms, it will exchange the group atom positions with other atoms from a defined swapList. Because the swapList can be big, swapGenerator can be assigned to multiple groups at the same time under the condition of all groups having the same length.

During stochastic engine runtime, whenever a swap generator is encountered, all sophisticated selection recurrence modes such as (refining, exploring) will be reduced to simple recurrence.

This class can’t be instantiated but its sub-classes might be.

Parameters:
  1. group (None, Group): The group instance.
  2. swapLength (Integer): The swap length that defines the length of the group and the length of the every swap sub-list in swapList.
  3. swapList (None, List): List of atoms index. If None is given, no swapping or exchanging will be performed. If List is given, it must contain lists of atom indexes where every sub-list must have the same number of atoms as the group.
swapLength

Swap length.

swapList

Swap list.

groupAtomsIndexes

Last selected group atoms index.

swapAtomsIndexes

Last swap atoms index.

set_swap_length(swapLength)

Set swap length. it will empty and reset swaplist automatically.

Parameters:
  1. swapLength (Integer): The swap length that defines the length of the group and the length of the every swap sub-list in swapList.
set_group(group)

Set the MoveGenerator group.

Parameters:
  1. group (None, Group): group instance.
set_swap_list(swapList)

Set the swap-list to exchange atoms position from.

Parameters:
  1. swapList (None, List): The list of atoms.

    If None is given, no swapping or exchanging will be performed.

    If List is given, it must contain lists of atom indexes where every sub-list length must be equal to swapLength.

append_to_swap_list(subList)

Append a sub list to swap list.

Parameters:
  1. subList (List): The sub-list of atoms index to append to swapList.
get_ready_for_move(engine, groupAtomsIndexes)

Set the swap generator ready to perform a move. Unlike a normal move generator, swap generators will affect not only the selected atoms but other atoms as well. Therefore at stochastic engine runtime, selected atoms will be extended to all affected atoms by the swap.

This method is called automatically upon stochastic engine runtime to ensure that all affect atoms with the swap are updated.

Parameters:
  1. engine (fullrmc.Engine): The stochastic engine calling for the move.
  2. groupAtomsIndexes (numpy.ndarray): The atoms index to swap.
Returns:
  1. indexes (numpy.ndarray): All the atoms involved in the swap move including the given groupAtomsIndexes.
class fullrmc.Core.MoveGenerator.PathGenerator(group=None, path=None, randomize=False)

Bases: fullrmc.Core.MoveGenerator.MoveGenerator

PathGenerator is a MoveGenerator sub-class where moves definitions are pre-stored in a path and get pulled out at every move step.

This class can’t be instantiated but its sub-classes might be.

Parameters:
  1. group (None, Group): The group instance.
  2. path (None, list): The list of moves.
  3. randomize (boolean): Whether to pull moves randomly from path or pull moves in order at every step.
step

Current step number.

path

Path list of moves.

randomize

Randomize flag.

check_path(path)

Check the generator’s path.

This method must be overloaded in all PathGenerator sub-classes.

Parameters:
  1. path (list): The list of moves.
normalize_path(path)

Normalizes all path moves. It is called automatically upon set_path method is called.

This method can be overloaded in all MoveGenerator sub-classes.

Parameters:
  1. path (list): The list of moves.
Returns:
  1. path (list): The list of moves.
set_path(path)

Set the moves path.

Parameters:
  1. path (list): The list of moves.
set_randomize(randomize)

Set whether to randomize moves selection.

Parameters:
  1. randomize (boolean): Whether to pull moves randomly from path or pull moves in order at every step.
move(coordinates)

Move coordinates.

Parameters:
  1. coordinates (np.ndarray): The coordinates on which to apply the transformation.
Returns:
  1. coordinates (np.ndarray): The new coordinates after applying the transformation.
class fullrmc.Core.MoveGenerator.MoveGeneratorCombinator(group=None, combination=None, shuffle=False)

Bases: fullrmc.Core.MoveGenerator.MoveGenerator

MoveGeneratorCombinator combines all moves of a list of MoveGenerators and applies it at once.

Parameters:
  1. group (None, Group): The constraint stochastic engine.
  2. combination (list): The list of MoveGenerator instances.
  3. shuffle (boolean): Whether to shuffle generator instances at every move or to combine moves in the list order.
# import fullrmc modules
from fullrmc.Engine import Engine
from fullrmc.Core.MoveGenerator import MoveGeneratorCombinator
from fullrmc.Generators.Translations import TranslationGenerator
from fullrmc.Generators.Rotations import RotationGenerator

# create engine
ENGINE = Engine(path='my_engine.rmc')

# set pdb file
ENGINE.set_pdb('system.pdb')

# Add constraints ...
# Re-define groups if needed ...

##### Define each group move generator as a combination of a translation and a rotation. #####
# create recursive group selector. Recurrence is set to 20 with explore flag set to True.
# shuffle is set to True which means that at every selection the order of move generation
# is random. At one step a translation is performed prior to rotation and in another step
# the rotation is performed at first.
# selected from the collector.
for g in ENGINE.groups:
    # create translation generator
    TMG = TranslationGenerator(amplitude=0.2)
    # create rotation generator only when group length is bigger than 1.
    if len(g)>1:
        RMG = RotationGenerator(amplitude=2)
        MG  = MoveGeneratorCombinator(collection=[TMG,RMG],shuffle=True)
    else:
        MG  = MoveGeneratorCombinator(collection=[TMG],shuffle=True)
    g.set_move_generator( MG )
shuffle

Shuffle flag.

combination

Combination list of MoveGenerator instances.

check_group(group)

Checks the generator’s group. This methods always returns True because normally all combination MoveGenerator instances groups are checked.

This method must NOT be overloaded unless needed.

Parameters:
  1. group (Group): the Group instance
set_group(group)

Set the MoveGenerator group.

Parameters:
  1. group (None, Group): group instance.
set_combination(combination)

Set the generators combination list.

Parameters:
  1. combination (list): The list of MoveGenerator instances.
set_shuffle(shuffle)

Set whether to shuffle moves generator.

Parameters:
  1. shuffle (boolean): Whether to shuffle generator instances at every move or to combine moves in the list order.
move(coordinates)

Move coordinates.

Parameters:
  1. coordinates (np.ndarray): The coordinates on which to apply the transformation.
Returns:
  1. coordinates (np.ndarray): The new coordinates after applying the transformation.
class fullrmc.Core.MoveGenerator.MoveGeneratorCollector(group=None, collection=None, randomize=True, weights=None)

Bases: fullrmc.Core.MoveGenerator.MoveGenerator

MoveGeneratorCollector collects MoveGenerators instances and applies the move of one instance at every step.

Parameters:
  1. group (None, Group): The constraint stochastic engine.
  2. collection (list): The list of MoveGenerator instances.
  3. randomize (boolean): Whether to pull MoveGenerator instance randomly from collection list or in order.
  4. weights (None, list): Generators selections Weights list. It must be None for equivalent weighting or list of (generatorIndex, weight) tuples. If randomize is False, weights list is ignored upon generator selection from collection.
# import fullrmc modules
from fullrmc.Engine import Engine
from fullrmc.Core.MoveGenerator import MoveGeneratorCollector
from fullrmc.Generators.Translations import TranslationGenerator
from fullrmc.Generators.Rotations import RotationGenerator

# create engine
ENGINE = Engine(path='my_engine.rmc')

# set pdb file
ENGINE.set_pdb('system.pdb')

# Add constraints ...
# Re-define groups if needed ...

##### Define each group move generator as a combination of a translation and a rotation. #####
# create recursive group selector. Recurrence is set to 20 with explore flag set to True.
# randomize is set to True which means that at every selection a generator is randomly
# selected from the collector.
for g in ENGINE.groups:
    # create translation generator
    TMG = TranslationGenerator(amplitude=0.2)
    # create rotation generator only when group length is bigger than 1.
    if len(g)>1:
        RMG = RotationGenerator(amplitude=2)
        MG  = MoveGeneratorCollector(collection=[TMG,RMG],randomize=True)
    else:
        MG  = MoveGeneratorCollector(collection=[TMG],randomize=True)
    g.set_move_generator( MG )
randomize

Randomize flag.

collection

List of MoveGenerator instances.

generatorsWeight

Generators selection weights list.

selectionScheme

Selection scheme.

set_group(group)

Set the MoveGenerator group.

Parameters:
  1. group (None, Group): group instance.
check_group(group)

Check the generator’s group. This methods always returns True because normally all collection MoveGenerator instances groups are checked.

This method must NOT be overloaded unless needed.

Parameters:
  1. group (Group): the Group instance.
set_collection(collection)

Set the generators instances collection list.

Parameters:
  1. collection (list): The list of move generator instance.
set_randomize(randomize)

Set whether to randomize MoveGenerator instance selection from collection list.

Parameters:
  1. randomize (boolean): Whether to pull MoveGenerator instance randomly from collection list or in order.
set_weights(weights)

Set groups selection weighting scheme.

Parameters:
  1. weights (None, list): Generators selections Weights list. It must be None for equivalent weighting or list of (generatorIndex, weight) tuples. If randomize is False, weights list is ignored upon generator selection from collection.
set_selection_scheme()

Set selection scheme.

move(coordinates)

Move coordinates.

Parameters:
  1. coordinates (np.ndarray): The coordinates on which to apply the transformation.
Returns:
  1. coordinates (np.ndarray): The new coordinates after applying the transformation.
fullrmc.Core.MoveGenerator.generate_random_float()

random() -> x in the interval [0, 1).

boundary conditions collection

This is a C compiled module to compute boundary conditions related calculations

fullrmc.Core.boundary_conditions_collection.get_reciprocal_basis()

Computes reciprocal box matrix.

Arguments:
  1. basis (float32 array): The (3,3) box matrix
Returns:
  1. rbasis (float32 array): The (3,3) reciprocal box matrix.
fullrmc.Core.boundary_conditions_collection.transform_coordinates()

Transforms coordinates array using a transformation matrix.

Arguments:
  1. transMatrix (float32 array): The (3,3) transformation matrix
  2. coords (float32 array): The (N,3) coordinates array.
Returns:
  1. transCoords (float32 array): The (N,3) transformed coordinates array.
fullrmc.Core.boundary_conditions_collection.box_coordinates_real_distances()

Computes atomic real distances given box coordinates.

Arguments:
  1. atomIndex (int32): The index of atom to compute the distance from.
  2. indexes (int32 array): The list of atom indexes to compute the distance to
  3. boxCoords (float32 array): The (N,3) box coordinates array.
  4. basis (float32 array): The (3,3) box matrix
Returns:
  1. distances (float32 array): The (N,) distances array.

reciprocal space

This is a C compiled module to compute transformations from real to reciprocal space and vice versa.

fullrmc.Core.reciprocal_space.gr_to_sq()

Transform pair correlation function g(r) to static structure factor S(q).

Arguments:
  1. distances (float32 (n,) numpy.ndarray): The g(r) bins positions in real space.
  2. gr (float32 (n,) numpy.ndarray): The pair correlation function g(r) data.
  3. qrange (float32 (m,) numpy.ndarray): The S(q) bins positions in reciprocal space.
  4. rho (float32) [default=1]: The number density of the system.
Returns:
  1. sq (float32 (m,) numpy.ndarray): The static structure factor S(q) data.
fullrmc.Core.reciprocal_space.Gr_to_sq()

Transform pair distribution function G(r) to static structure factor S(q).

Arguments:
  1. distances (float32 (n,) numpy.ndarray): The G(r) bins positions in real space.
  2. Gr (float32 (n,) numpy.ndarray): The pair correlation function Gr) data.
  3. qrange (float32 (m,) numpy.ndarray): The S(q) bins positions in reciprocal space.
Returns:
  1. sq (float32 (m,) numpy.ndarray): The static structure factor S(q) data.
fullrmc.Core.reciprocal_space.sq_to_Gr()

Transform static structure factor S(q) to pair distribution function G(r).

Arguments:
  1. qValues (float32 (m,) numpy.ndarray): The S(q) bins positions in reciprocal space.
  2. rValues (float32 (n,) numpy.ndarray): The G(r) bins positions in real space.
  3. sq (float32 (m,) numpy.ndarray): The static structure factor S(q) data.
Returns:
  1. Gr (float32 (n,) numpy.ndarray): The pair correlation function Gr) data.

pairs distances

This is a C compiled module to compute atomic pair distances.

fullrmc.Core.pairs_distances.from_to_points_differences()

Compute point to point vector difference between two atomic coordinates arrays taking into account periodic or infinite boundary conditions. Difference is calculated as the following:

\[differences[i,:] = boundaryConditions( pointsTo[i,:] - pointsFrom[i,:] )\]
Arguments:
  1. pointsFrom (float32 (n,3) numpy.ndarray): The first atomic coordinates array of the same shape as pointsTo.
  2. pointsTo (float32 (n,3) numpy.ndarray): The second atomic coordinates array of the same shape as pointsFrom.
  3. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. differences (float32 (n,3) numpy.ndarray): The computed differences array.
fullrmc.Core.pairs_distances.pairs_differences_to_point()

Compute differences between one atomic coordinates arrays to a point coordinates taking into account periodic or infinite boundary conditions. Difference is calculated as the following:

\[differences[i,:] = boundaryConditions( point[0,:] - coords[i,:] )\]
Arguments:
  1. point (float32 (1,3) numpy.ndarray): The atomic coordinates point.
  2. coords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  3. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. differences (float32 (n,3) numpy.ndarray): The computed differences array.
fullrmc.Core.pairs_distances.pairs_differences_to_indexcoords()

Compute differences between one atomic coordinates arrays to a point coordinates given its index in the coordinates array and taking into account periodic or infinite boundary conditions. Difference is calculated as the following:

\[differences[i,:] = boundaryConditions( coords[atomIndex,:] - coords[i,:] )\]
Arguments:
  1. atomIndex (int32): The index of the atomic coordinates point.
  2. coords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  3. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. differences (float32 (n,3) numpy.ndarray): The computed differences array.
fullrmc.Core.pairs_distances.pairs_differences_to_multi_points()

Compute differences between one atomic coordinates arrays to a multiple points coordinates taking into account periodic or infinite boundary conditions. Difference is calculated as the following:

\[differences[i,:,k] = boundaryConditions( points[k,:] - coords[i,:] )\]
Arguments:
  1. points (float32 (k,3) numpy.ndarray): The multiple atomic coordinates points.
  2. coords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  3. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. differences (float32 (n,3,k) numpy.ndarray): The computed differences array.
fullrmc.Core.pairs_distances.pairs_differences_to_multi_indexcoords()

Compute differences between one atomic coordinates arrays to a points coordinates given their indexes in the coordinates array and taking into account periodic or infinite boundary conditions. Difference is calculated as the following:

\[differences[i,:,k] = boundaryConditions( coords[indexes[k],:] - coords[i,:] )\]
Arguments:
  1. indexes (int32 (k,3) numpy.ndarray): The atomic coordinates indexes array.
  2. coords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  3. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. differences (float32 (n,3,k) numpy.ndarray): The computed differences array.
fullrmc.Core.pairs_distances.pairs_distances_to_point()

Compute distances between one atomic coordinates arrays to a point coordinates taking into account periodic or infinite boundary conditions. Distances is calculated as the following:

\[distances[i] = \sqrt{ \sum_{d}^{3}{ boundaryConditions( point[0,d] - coords[i,d] )^{2}} }\]
Arguments:
  1. point (float32 (1,3) numpy.ndarray): The atomic coordinates point.
  2. coords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  3. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. distances (float32 (n,) numpy.ndarray): The computed distances array.
fullrmc.Core.pairs_distances.pairs_distances_to_indexcoords()

Compute distances between one atomic coordinates arrays to a points coordinates given their indexes in the coordinates array and taking into account periodic or infinite boundary conditions. Distances is calculated as the following:

\[distances[i] = \sqrt{ \sum_{d}^{3}{ boundaryConditions( coords[atomIndex[i],d] - coords[i,d] )^{2}} }\]
Arguments:
  1. point (float32 (1,3) numpy.ndarray): The atomic coordinates point.
  2. coords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  3. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. distances (float32 (n,) numpy.ndarray): The computed distances array.
fullrmc.Core.pairs_distances.pairs_distances_to_multi_points()

Compute distances between one atomic coordinates arrays to a multiple points coordinates taking into account periodic or infinite boundary conditions. Distances is calculated as the following:

\[distances[i,k] = \sqrt{ \sum_{d}^{3}{ boundaryConditions( points[k,d] - coords[i,d] )^{2}} }\]
Arguments:
  1. points (float32 (k,3) numpy.ndarray): The multiple atomic coordinates points.
  2. coords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  3. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. distances (float32 (n,) numpy.ndarray): The computed distances array.
fullrmc.Core.pairs_distances.pairs_distances_to_multi_indexcoords()

Compute distances between one atomic coordinates arrays to a points coordinates given their indexes in the coordinates array and taking into account periodic or infinite boundary conditions. Distances is calculated as the following:

\[distances[i,k] = \sqrt{ \sum_{d}^{3}{ boundaryConditions( coords[indexes[k],:] - coords[i,d] )^{2}} }\]
Arguments:
  1. indexes (int32 (k,3) numpy.ndarray): The atomic coordinates indexes array.
  2. coords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  3. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. distances (float32 (n,) numpy.ndarray): The computed distances array.

atomic coordination number

This is a C compiled module to compute atomic bonds.

fullrmc.Core.atomic_coordination.single_atom_single_shell_subdists()
fullrmc.Core.atomic_coordination.single_atom_single_shell_totdists()
fullrmc.Core.atomic_coordination.single_atom_single_shell_coords()
fullrmc.Core.atomic_coordination.single_atom_multi_shells_totdists()
fullrmc.Core.atomic_coordination.single_atom_multi_shells_coords()
fullrmc.Core.atomic_coordination.single_atom_coord_number_totdists()
fullrmc.Core.atomic_coordination.single_atom_coord_number_coords()
fullrmc.Core.atomic_coordination.multi_atoms_coord_number_totdists()
fullrmc.Core.atomic_coordination.multi_atoms_coord_number_coords()
fullrmc.Core.atomic_coordination.all_atoms_coord_number_totdists()
fullrmc.Core.atomic_coordination.all_atoms_coord_number_coords()

atomic distances

This is a C compiled module to compute atomic inter-molecular distances.

fullrmc.Core.atomic_distances.single_atomic_distances_dists()

Computes the inter-molecular distances constraint of a single atom given a distances array.

Arguments:
  1. atomIndex (int32): The index of the atom.
  2. distances (float32 array): The distances array of the atom with the rest of atoms.
  3. moleculeIndex (int32 array): The molecule’s index array, assigning a molecule index for every atom.
  4. elementIndex (int32 array): The element’s index array, assigning an element index for every atom.
  5. dintra (float32 array): The (numberOfElements,numberOfElements,1) array for intra-molecular counted distances.
  6. dinter (float32 array): The (numberOfElements,numberOfElements,1) array for inter-molecular counted distances.
  7. nintra (float32 array): The (numberOfElements,numberOfElements,1) array for intra-molecular counted elements.
  8. ninter (float32 array): The (numberOfElements,numberOfElements,1) array for inter-molecular counted elements.
  9. lowerLimit (float32 array): The (numberOfElements,numberOfElements,1) array of lower distance limits.
  10. upperLimit (float32 array): The (numberOfElements,numberOfElements,1) array of upper distance limits.
  11. interMolecular (bool): Whether to consider inter-molecular distances. DEFAULT: True
  12. intraMolecular (bool): Whether to consider intra-molecular distances. DEFAULT: True
  13. countWithinLimits (bool): Whether to count distances and atoms found within the lower and upper limits or outside. DEFAULT: True
  14. reduceDistanceToUpper (bool): Whether to reduce counted distances to the difference between the found distance and the upper limit. When True, this flag has the higher priority. DEFAULT: False
  15. reduceDistanceToLower (bool): Whether to reduce counted distances to the difference between the found distance and the lower limit. When True, this flag may lose its priority for reduceDistanceToUpper if the later is True. DEFAULT: False
  16. reduceDistance (bool): Whether to reduce counted distances to the difference between the found distance and the closest limit. When True, this flag may lose its priority if any of reduceDistanceToLower or reduceDistanceToUpper is True. DEFAULT: False
  17. allAtoms (bool): Perform the calculation over all the atoms. If False calculation starts from the given atomIndex. DEFAULT: True
  18. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. dintra (float32 array): The updated (numberOfElements,numberOfElements,1) array for intra-molecular counted distances.
  2. dinter (float32 array): The updated (numberOfElements,numberOfElements,1) array for inter-molecular counted distances.
  3. nintra (float32 array): The updated (numberOfElements,numberOfElements,1) array for intra-molecular counted elements.
  4. ninter (float32 array): The updated (numberOfElements,numberOfElements,1) array for inter-molecular counted elements.
fullrmc.Core.atomic_distances.multiple_atomic_distances_coords()

Computes multiple atoms inter-molecular distances constraint given coordinates.

Arguments:
  1. indexes (int32 array): The atoms indexes array.
  2. boxCoords (float32 array): The whole system box coordinates.
  3. basis (float32 array): The box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. moleculeIndex (int32 array): The molecule’s index array, assigning a molecule index for every atom.
  6. elementIndex (int32 array): The element’s index array, assigning an element index for every atom.
  7. numberOfElements (int32): The number of elements in the system.
  8. lowerLimit (float32 array): The (numberOfElements,numberOfElements,1) array of lower distance limits.
  9. upperLimit (float32 array): The (numberOfElements,numberOfElements,1) array of upper distance limits.
  10. interMolecular (bool): Whether to consider inter-molecular distances. DEFAULT: True
  11. intraMolecular (bool): Whether to consider intra-molecular distances. DEFAULT: True
  12. countWithinLimits (bool): Whether to count distances and atoms found within the lower and upper limits or outside.
  13. reduceDistanceToUpper (bool): Whether to reduce counted distances to the difference between the found distance and the upper limit. When True, this flag has the higher priority. DEFAULT: False
  14. reduceDistanceToLower (bool): Whether to reduce counted distances to the difference between the found distance and the lower limit. When True, this flag may lose its priority for reduceDistanceToUpper if the later is True. DEFAULT: False
  15. reduceDistance (bool): Whether to reduce counted distances to the difference between the found distance and the closest limit. When True, this flag may lose its priority if any of reduceDistanceToLower or reduceDistanceToUpper is True. DEFAULT: False
  16. allAtoms (bool): Perform the calculation over all the atoms. If False calculation starts from the given atomIndex. DEFAULT: True
  17. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. dintra (float32 array): The created (numberOfElements,numberOfElements,1) array for intra-molecular counted distances.
  2. dinter (float32 array): The created (numberOfElements,numberOfElements,1) array for inter-molecular counted distances.
  3. nintra (float32 array): The created (numberOfElements,numberOfElements,1) array for intra-molecular counted elements.
  4. ninter (float32 array): The created (numberOfElements,numberOfElements,1) array for inter-molecular counted elements.
fullrmc.Core.atomic_distances.multiple_atomic_distances_dists()

Computes multiple atoms inter-molecular distances constraint given distances.

Arguments:
  1. indexes (int32 array): The atoms indexes array.
  2. distances (float32 array): The distances array of the atoms with the rest of atoms.
  3. moleculeIndex (int32 array): The molecule’s index array, assigning a molecule index for every atom.
  4. elementIndex (int32 array): The element’s index array, assigning an element index for every atom.
  5. numberOfElements (int32): The number of elements in the system.
  6. lowerLimit (float32 array): The (numberOfElements,numberOfElements,1) array of lower distance limits.
  7. upperLimit (float32 array): The (numberOfElements,numberOfElements,1) array of upper distance limits.
  8. interMolecular (bool): Whether to consider inter-molecular distances. DEFAULT: True
  9. intraMolecular (bool): Whether to consider intra-molecular distances. DEFAULT: True
  10. countWithinLimits (bool): Whether to count distances and atoms found within the lower and upper limits or outside.
  11. reduceDistanceToUpper (bool): Whether to reduce counted distances to the difference between the found distance and the upper limit. When True, this flag has the higher priority. DEFAULT: False
  12. reduceDistanceToLower (bool): Whether to reduce counted distances to the difference between the found distance and the lower limit. When True, this flag may lose its priority for reduceDistanceToUpper if the later is True. DEFAULT: False
  13. reduceDistance (bool): Whether to reduce counted distances to the difference between the found distance and the closest limit. When True, this flag may lose its priority if any of reduceDistanceToLower or reduceDistanceToUpper is True. DEFAULT: False
  14. allAtoms (bool): Perform the calculation over all the atoms. If False calculation starts from the given atomIndex. DEFAULT: True
  15. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. dintra (float32 array): The created (numberOfElements,numberOfElements,1) array for intra-molecular counted distances.
  2. dinter (float32 array): The created (numberOfElements,numberOfElements,1) array for inter-molecular counted distances.
  3. nintra (float32 array): The created (numberOfElements,numberOfElements,1) array for intra-molecular counted elements.
  4. ninter (float32 array): The created (numberOfElements,numberOfElements,1) array for inter-molecular counted elements.
fullrmc.Core.atomic_distances.full_atomic_distances_coords()

Computes all atoms inter-molecular distances constraint given coordinates.

Arguments:
  1. boxCoords (float32 array): The whole system box coordinates.
  2. basis (float32 array): The box vectors.
  3. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  4. moleculeIndex (int32 array): The molecule’s index array, assigning a molecule index for every atom.
  5. elementIndex (int32 array): The element’s index array, assigning an element index for every atom.
  6. numberOfElements (int32): The number of elements in the system.
  7. lowerLimit (float32 array): The (numberOfElements,numberOfElements,1) array of lower distance limits.
  8. upperLimit (float32 array): The (numberOfElements,numberOfElements,1) array of upper distance limits.
  9. interMolecular (bool): Whether to consider inter-molecular distances. DEFAULT: True
  10. intraMolecular (bool): Whether to consider intra-molecular distances. DEFAULT: True
  11. countWithinLimits (bool): Whether to count distances and atoms found within the lower and upper limits or outside.
  12. reduceDistanceToUpper (bool): Whether to reduce counted distances to the difference between the found distance and the upper limit. When True, this flag has the higher priority. DEFAULT: False
  13. reduceDistanceToLower (bool): Whether to reduce counted distances to the difference between the found distance and the lower limit. When True, this flag may lose its priority for reduceDistanceToUpper if the later is True. DEFAULT: False
  14. reduceDistance (bool): Whether to reduce counted distances to the difference between the found distance and the closest limit. When True, this flag may lose its priority if any of reduceDistanceToLower or reduceDistanceToUpper is True. DEFAULT: False
  15. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. dintra (float32 array): The created (numberOfElements,numberOfElements,1) array for intra-molecular counted distances.
  2. dinter (float32 array): The created (numberOfElements,numberOfElements,1) array for inter-molecular counted distances.
  3. nintra (float32 array): The created (numberOfElements,numberOfElements,1) array for intra-molecular counted elements.
  4. ninter (float32 array): The created (numberOfElements,numberOfElements,1) array for inter-molecular counted elements.
fullrmc.Core.atomic_distances.full_atomic_distances_dists()

Computes all atoms inter-molecular distances constraint given distances.

Arguments:
  1. distances (float32 array): The distances array of the atoms with the rest of atoms.
  2. moleculeIndex (int32 array): The molecule’s index array, assigning a molecule index for every atom.
  3. elementIndex (int32 array): The element’s index array, assigning an element index for every atom.
  4. numberOfElements (int32): The number of elements in the system.
  5. lowerLimit (float32 array): The (numberOfElements,numberOfElements,1) array of lower distance limits.
  6. upperLimit (float32 array): The (numberOfElements,numberOfElements,1) array of upper distance limits.
  7. interMolecular (bool): Whether to consider inter-molecular distances. DEFAULT: True
  8. intraMolecular (bool): Whether to consider intra-molecular distances. DEFAULT: True
  9. countWithinLimits (bool): Whether to count distances and atoms found within the lower and upper limits or outside.
  10. reduceDistanceToUpper (bool): Whether to reduce counted distances to the difference between the found distance and the upper limit. When True, this flag has the higher priority. DEFAULT: False
  11. reduceDistanceToLower (bool): Whether to reduce counted distances to the difference between the found distance and the lower limit. When True, this flag may lose its priority for reduceDistanceToUpper if the later is True. DEFAULT: False
  12. reduceDistance (bool): Whether to reduce counted distances to the difference between the found distance and the closest limit. When True, this flag may lose its priority if any of reduceDistanceToLower or reduceDistanceToUpper is True. DEFAULT: False
  13. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. dintra (float32 array): The created (numberOfElements,numberOfElements,1) array for intra-molecular counted distances.
  2. dinter (float32 array): The created (numberOfElements,numberOfElements,1) array for inter-molecular counted distances.
  3. nintra (float32 array): The created (numberOfElements,numberOfElements,1) array for intra-molecular counted elements.
  4. ninter (float32 array): The created (numberOfElements,numberOfElements,1) array for inter-molecular counted elements.

bonds

This is a C compiled module to compute atomic bonds.

fullrmc.Core.bonds.full_bonds_coords()

It calculates the bonds constraint of box coordinates.

Arguments:
  1. idx1 (int32 (n,) numpy.ndarray): First atoms index array
  2. idx2 (int32 (n,) numpy.ndarray): Second atoms index array
  3. lowerLimit (float32 (n,) numpy.ndarray): Lower limit or minimum bond length allowed.
  4. upperLimit (float32 (n,) numpy.ndarray): Upper limit or minimum bond length allowed.
  5. boxCoords (float32 (n,3) numpy.ndarray): The atomic coordinates array.
  6. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  7. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  8. reduceDistanceToUpper (bool): Whether to reduce bonds length found out of limits to the difference between the bond length and the upper limit. When True, this flag has the higher priority. DEFAULT: False
  9. reduceDistanceToLower (bool): Whether to reduce bonds length found out of limits to the difference between the bond length and the lower limit. When True, this flag may lose its priority for reduceDistanceToUpper if the later is True. DEFAULT: False
  10. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. bondsLength: The calculated bonds length
  2. reducedLengths: The calculated reduced bonds length

angles

This is a C compiled module to compute bonded atoms angle.

fullrmc.Core.angles.full_angles_coords()

Computes the angles constraint given bonded atoms vectors.

Arguments:
  1. central (int32 (n,) numpy.ndarray): The central atom indexes.
  2. left (int32 (n,) numpy.ndarray): The left atom indexes.
  3. right (int32 (n,) numpy.ndarray): The right atom indexes.
  4. lowerLimit (float32 (n,) numpy.ndarray): The angles lower limits.
  5. upperLimit (float32 (n,) numpy.ndarray): The angles upper limits.
  6. boxCoords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  7. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  8. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  9. reduceAngleToUpper (bool): Whether to reduce angle found out of limits to the difference between the angle and the upper limit. When True, this flag has the higher priority. DEFAULT: False
  10. reduceAngleToLower (bool): Whether to reduce angle found out of limits to the difference between the angle and the lower limit. When True, this flag may lose its priority for reduceAngleToUpper if the later is True. DEFAULT: False
  11. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. angles (float32 (n,) numpy.ndarray): The calculated angles (rad).
  2. reducedAngles (float32 (n,) numpy.ndarray): The reduced angles (rad).

dihedral angles

This is a C compiled module to compute improper angles.

fullrmc.Core.dihedral_angles.full_dihedral_angles_coords()

Computes the improper angles constraint between an improper atom and a plane atoms. The plane normal vector is calculated using the right-hand rule where (thumb=ox vector), (index=oy vector) hence (oz=normal=second finger)

Arguments:
  1. indexes1 (int32 (n,) numpy.ndarray): Diherdral first atom indexes.
  2. indexes2 (int32 (n,) numpy.ndarray): Diherdral second atom indexes.
  3. indexes3 (int32 (n,) numpy.ndarray): Diherdral third atom indexes.
  4. indexes4 (int32 (n,) numpy.ndarray): Diherdral fourth atom indexes.
  5. lowerLimit1 (float32 (n,) numpy.ndarray): First shells lower limit.
  6. upperLimit1 (float32 (n,) numpy.ndarray): First shells upper limits.
  7. lowerLimit2 (float32 (n,) numpy.ndarray): Second shells lower limit.
  8. upperLimit2 (float32 (n,) numpy.ndarray): Second shells upper limits.
  9. lowerLimit3 (float32 (n,) numpy.ndarray): Third shells lower limit.
  10. upperLimit3 (float32 (n,) numpy.ndarray): Third shells upper limits.
  11. boxCoords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  12. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  13. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  14. reduceAngleToUpper (bool): Whether to reduce angle found out of limits to the difference between the angle and the upper limit. When True, this flag has the higher priority. DEFAULT: False
  15. reduceAngleToLower (bool): Whether to reduce angle found out of limits to the difference between the angle and the lower limit. When True, this flag may lose its priority for reduceAngleToUpper if the later is True. DEFAULT: False
  16. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. angles: The calculated angles (rad).
  2. reducedAngles: The reduced angles (rad)

improper angles

This is a C compiled module to compute improper angles.

fullrmc.Core.improper_angles.full_improper_angles_coords()

Computes the improper angles constraint between an improper atom and a plane atoms. The plane normal vector is calculated using the right-hand rule where (thumb=ox vector), (index=oy vector) hence (oz=normal=second finger)

Arguments:
  1. improperIdxs (int32 (n,) numpy.ndarray): The improper atom indexes.
  2. oIdxs (int32 (n,) numpy.ndarray): The O atom indexes.
  3. xIdxs (int32 (n,) numpy.ndarray): The x atom indexes.
  4. yIdxs (int32 (n,) numpy.ndarray): The y atom indexes.
  5. lowerLimit (float32 (n,) numpy.ndarray): The angles lower limits.
  6. upperLimit (float32 (n,) numpy.ndarray): The angles upper limits.
  7. boxCoords (float32 (n,3) numpy.ndarray): The atomic coordinates array of the same shape as pointsFrom.
  8. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  9. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  10. reduceAngleToUpper (bool): Whether to reduce angle found out of limits to the difference between the angle and the upper limit. When True, this flag has the higher priority. DEFAULT: False
  11. reduceAngleToLower (bool): Whether to reduce angle found out of limits to the difference between the angle and the lower limit. When True, this flag may lose its priority for reduceAngleToUpper if the later is True. DEFAULT: False
  12. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. angles: The calculated angles (rad).
  2. reducedAngles: The reduced angles (rad)

pairs histogram

This is a C compiled module to compute pair distances histograms.

fullrmc.Core.pairs_histograms.single_pairs_histograms()

Computes the pair distribution histograms of a single atom given a distances array.

Arguments:
  1. atomIndex (int32): The index of the atom.
  2. distances (float32 array): The distances array.
  3. moleculeIndex (int32 array): The molecule’s index array, assigning a molecule index for every atom.
  4. elementIndex (int32 array): The element’s index array, assigning an element index for every atom.
  5. hintra (float32 array): The (numberOfElements,numberOfElements,1) array for intra-molecular distances histograms.
  6. hinter (float32 array): The (numberOfElements,numberOfElements,1) array for inter-molecular distances histograms.
  7. minDistance (float32): The minimum distance to be counted in the histogram.
  8. maxDistance (float32): The maximum distance to be counted in the histogram.
  9. bin (float32): The histogram bin size.
  10. allAtoms (bool): Perform the calculation over all the atoms. If False calculation starts from the given atomIndex. DEFAULT: True
  11. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. hintra (float32 array): The updated (numberOfElements,numberOfElements,1) array for intra-molecular distances histograms.
  2. hinter (float32 array): The updated (numberOfElements,numberOfElements,1) array for inter-molecular distances histograms.
fullrmc.Core.pairs_histograms.multiple_pairs_histograms_coords()

Computes the pair distribution histograms of multiple atoms given atomic coordinates.

Arguments:
  1. indexes (int32 (k,3) numpy.ndarray): The atomic coordinates indexes array.
  2. coords (float32 (n,3) numpy.ndarray): The atomic coordinates array.
  3. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  4. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  5. moleculeIndex (int32 array): The molecule’s index array, assigning a molecule index for every atom.
  6. elementIndex (int32 array): The element’s index array, assigning an element index for every atom.
  7. numberOfElements (int32): The number of elements in the system.
  8. minDistance (float32): The minimum distance to be counted in the histogram.
  9. maxDistance (float32): The maximum distance to be counted in the histogram.
  10. bin (float32): The histogram bin size.
  11. histSize(int32): The histograms size.
  12. allAtoms (bool): Perform the calculation over all the atoms. If False calculation starts from the given atomIndex. DEFAULT: True
  13. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. hintra (float32 array): The updated (numberOfElements,numberOfElements,1) array for intra-molecular distances histograms.
  2. hinter (float32 array): The updated (numberOfElements,numberOfElements,1) array for inter-molecular distances histograms.
fullrmc.Core.pairs_histograms.multiple_pairs_histograms_dists()

Computes the pair distribution histograms of multiple atoms given atomic distances.

Arguments:
  1. indexes (int32 (k,3) numpy.ndarray): The atomic coordinates indexes array.
  2. distances (float32 array): The distances array.
  3. moleculeIndex (int32 array): The molecule’s index array, assigning a molecule index for every atom.
  4. elementIndex (int32 array): The element’s index array, assigning an element index for every atom.
  5. numberOfElements (int32): The number of elements in the system.
  6. minDistance (float32): The minimum distance to be counted in the histogram.
  7. maxDistance (float32): The maximum distance to be counted in the histogram.
  8. bin (float32): The histogram bin size.
  9. histSize(int32): The histograms size.
  10. allAtoms (bool): Perform the calculation over all the atoms. If False calculation starts from the given atomIndex. DEFAULT: True
  11. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. hintra (float32 array): The updated (numberOfElements,numberOfElements,1) array for intra-molecular distances histograms.
  2. hinter (float32 array): The updated (numberOfElements,numberOfElements,1) array for inter-molecular distances histograms.
fullrmc.Core.pairs_histograms.full_pairs_histograms_coords()

Computes the pair distribution histograms of multiple atoms given atomic coordinates.

Arguments:
  1. boxCoords (float32 (n,3) numpy.ndarray): The atomic coordinates array.
  2. basis (float32 (3,3) numpy.ndarray): The (3x3) boundary conditions box vectors.
  3. isPBC (bool): Whether it is a periodic boundary conditions or infinite.
  4. moleculeIndex (int32 array): The molecule’s index array, assigning a molecule index for every atom.
  5. elementIndex (int32 array): The element’s index array, assigning an element index for every atom.
  6. numberOfElements (int32): The number of elements in the system.
  7. minDistance (float32): The minimum distance to be counted in the histogram.
  8. maxDistance (float32): The maximum distance to be counted in the histogram.
  9. bin (float32): The histogram bin size.
  10. histSize(int32): The histograms size.
  11. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. hintra (float32 array): The updated (numberOfElements,numberOfElements,1) array for intra-molecular distances histograms.
  2. hinter (float32 array): The updated (numberOfElements,numberOfElements,1) array for inter-molecular distances histograms.
fullrmc.Core.pairs_histograms.full_pairs_histograms_dists()

Computes the pair distribution histograms of multiple atoms given atomic distances.

Arguments:
  1. distances (float32 array): The distances array.
  2. moleculeIndex (int32 array): The molecule’s index array, assigning a molecule index for every atom.
  3. elementIndex (int32 array): The element’s index array, assigning an element index for every atom.
  4. numberOfElements (int32): The number of elements in the system.
  5. minDistance (float32): The minimum distance to be counted in the histogram.
  6. maxDistance (float32): The maximum distance to be counted in the histogram.
  7. bin (float32): The histogram bin size.
  8. histSize(int32): The histograms size.
  9. ncores (int32) [default=1]: The number of cores to use.
Returns:
  1. hintra (float32 array): The updated (numberOfElements,numberOfElements,1) array for intra-molecular distances histograms.
  2. hinter (float32 array): The updated (numberOfElements,numberOfElements,1) array for inter-molecular distances histograms.