CDTK.Tools.Geometry module

class CDTK.Tools.Geometry.Coordinate(**opts)[source]

Bases: object

List of numbers representing a coordinate vector

It can represent a normal mode, or in general anything that can be represented by a set of coordinates. Can also contain the masses related to those coordinates. Uses numpy.array to store coordinates and other information.

coordinates – array, the coordinates masses – masses corresponding to each coordinate, in AMU labels – list of labels representing each coordinate

property coordinates
getCoordinates()[source]
getLabels()[source]
getMasses()[source]
property labels
property masses
setCoordinates(x)[source]
setLabels(s)[source]
setMasses(x)[source]
class CDTK.Tools.Geometry.XYZGeometry(**opts)[source]

Bases: object

Stores a cartesian geometry.

property atomPositions
cartesianProjector(subgroups=None)[source]

Returns a (3N, 3N) matrix that projects out the subspace spanned by the vectors corresponding to translation, overall rotation, and subgroup rotation around the axis with the smallest moment of inertia. Can be used to project Cartesian Hessian matrices prior to normal mode analysis.

subgroups – array with arrays of subgroup indices

cartesianRotationVector(ix, subgroup=None)[source]

Returns the vector for minimal rotation around axis ix in unweighted cartesian coordinates

ix - index of the vector subgroup - index array of atoms forming a subgroup. If None (default),

the whole molecule is rotated. For subgroup rotations, regardless of the provided index, the rotation with the smallest moment of inertia is returned

cartesianTranslationVector(ix)[source]

Returns the vector for displacement along axis ix in unweighted cartesian coordinates

cartesianVectors(subgroups=None)[source]

Returns a (3N,6+s) matrix where s is the number of subgroups. The matrix contains, columnwise, orthonormal (3N) vectors. These vectors form a basis of the space spanned by the vectors corresponding to translation, overall rotation, and subgroup rotation around the axis with the smallest moment of inertia of the respective group. They can be used to project out normal modes corresponding to these movements from a Hessian matrix. The vectors are given in unweighted, cartesian coordinates.

subgroups - array with arrays of subgroup indices

centerGeometry()[source]

Move the center of mass of the molecule to the origin

centerOfMass()[source]

Return a 3 components vector with the center of mass of the molecule

getAtomPositions()[source]
getLabels()[source]
getMasses()[source]
getPbcDims()[source]
property labels
property masses
momentsOfInertia(sort=True, subgroup=None)[source]

Return a 3x3 matrix where the vectors v[:,i] are the moments of inertia that diagonalize the inertia tensor

sort - if True (default), return sorted eigenvalues and eigenvectors subgroup - index array, if it is given, the moments of inertia for a

rotation of the atoms specified by the indices is returned. Default: None

property pbcDims
radiusGyration()[source]

Returns the radius of gyration defined as

rg**2 = sum_i m_i * (r_i - r_COM)**2 / M

setAtomPositions(x)[source]
setLabels(s)[source]
setMasses(x)[source]
setPbcDims(x)[source]
translate(translation)[source]

Translate the system according translation translation – 3D vector

CDTK.Tools.Geometry.angularMomentumEckart(euler, *args)[source]

Calculate the Eckart conditions given 3 euler angles

This is the functions to be zeroed (minimized) if x is to be roted to the Eckart frame defined by xref

euler – euler angles a,b,c args[0] – structure; (Nx3) array args[1] – reference structure; (Nx3) array args[2] – mass vector; N array

Returns the square of the norm of the vector:

amv = sum_i { mass[i] * (dx[i] cross xref[i]) }

CDTK.Tools.Geometry.centerOfMass(geom, m)[source]

Return 3D array with the center of mass of geom

xs – list or array of 3D arrays with the positions of the atoms m – list or array with the mass of each atom in xs

CDTK.Tools.Geometry.eckartFrame(x, xref, m)[source]

Calculate the euler angles that rotate x to the Eckart frame w.r.t xref

x,xref – Nx3 array with atom coordinates m – array with masses, len(m) = N

CDTK.Tools.Geometry.inertiaTensor(geom, m)[source]

Return a 3x3 matrix with the inertia tensor corresponding to the configuration x and masses m

CDTK.Tools.Geometry.pbc_translation_matrix(pbc_dims)[source]

Periodic boundary condition: calculate translation matrix from box dimensions Minimum image convention

Parameters: pbc_dims – array containing the box dimensions along x/y/z

Returns: [3,3,3,3] matrix containing the translation vectors and their linear combinations

e.g. T_xyz_3D[0,2,1,:] is the vector that generates a translation in -x,+y,0z direction

CDTK.Tools.Geometry.pbc_translation_vectors(pbc_dims)[source]

Periodic boundary condition: calculate translation vectors from box dimensions

Parameters: pbc_dims – array containing the box dimensions along x/y/z

Returns: [3,3] matrix containing the translation vectors along x/y/z

CDTK.Tools.Geometry.rotateGeometry(x, euler)[source]