Javelin

Contents:

About

Javelin is a moderisation of DISCUS written in the Python programming language

Installation

Install the latest release

  • Update once we have a release

Requirements

Optional:

  • ASE (to use the ase atoms structue)

Development

A development enviroment can easily be set up with either conda of PyPI

Using Conda

conda env create
source activate javelin
python setup.py install

Using PyPI

Tests

The unit tests can be run with pytest

Install with conda

conda install pytest

Install with PyPI

pip install pytest

Tutorials

Modules

List of javelin modules:

fourier

This module define the Fourier class and other functions related to the fourier transformation.

class javelin.fourier.Fourier[source]

The Fourier class contains everything required to calculate the diffuse scattering. The only required thing to be set is javelin.fourier.Fourier.structure. There are defaults for all other options including grid, radiation, average structure subtraction and lots options.

Examples:
>>> from javelin.structure import Structure
>>> fourier = Fourier()
>>> fourier.structure = Structure()
>>> print(fourier)
Structure         : Structure(False, a=1.0, b=1.0, c=1.0, alpha=90.0, beta=90.0, gamma=90.0)
Radiation         : neutron
Fourier volume    : complete crystal
Aver. subtraction : False

Reciprocal layer  :
lower left  corner :     [0 0 0]
lower right corner :     [1 0 0]
upper left  corner :     [0 1 0]
top   left  corner :     [0 0 1]

hor. increment     :     [ 0.01  0.    0.  ]
vert. increment    :     [ 0.    0.01  0.  ]
top   increment    :     [0 0 1]

# of points        :     101 x 101 x 1
>>> results = fourier.calc()
>>> print(results) 
<xarray.DataArray 'Intensity' ([1 0 0]: 101, [0 1 0]: 101, [0 0 1]: 1)>
array([[[ 0.],
        [ 0.],
        ...,
        [ 0.],
        [ 0.]],

       [[ 0.],
        [ 0.],
        ...,
        [ 0.],
        [ 0.]],

       ...,
       [[ 0.],
        [ 0.],
        ...,
        [ 0.],
        [ 0.]],

       [[ 0.],
        [ 0.],
        ...,
        [ 0.],
        [ 0.]]])
Coordinates:
  * [1 0 0]  ([1 0 0]) float64 0.0 0.01 0.02 0.03 0.04 0.05 0.06 0.07 0.08 ...
  * [0 1 0]  ([0 1 0]) float64 0.0 0.01 0.02 0.03 0.04 0.05 0.06 0.07 0.08 ...
  * [0 0 1]  ([0 0 1]) float64 0.0
Attributes:
    units:    r.l.u
approximate

This sets the options of calculating the approximate scattering instead of exact. This is much quicker and is likely good enough for most cases.

Getter:Returns bool of approximate scattering option
Setter:Sets whether approximate sacttering is calculated
Type:bool
average

This sets the options of calculating average structure and subtracted it from the simulated scattering

Getter:Returns bool of average structure subtraction option
Setter:Sets whether average stucture is subtracted
Type:bool
calc()[source]

Calculates the fourier transform

Returns:DataArray containing calculated diffuse scattering
Return type:xarray.DataArray
calc_average()[source]

Calculates the scattering from the avarage structure

Returns:DataArray containing calculated average scattering
Return type:xarray.DataArray
grid = None

The grid attribute defines the reciprocal volume from which the scattering will be calculated. Must of type javelin.grid.Grid And check javelin.grid.Grid for details on how to change the grid.

lots

The size of lots

Getter:Returns the lots size
Setter:Sets the lots size
Type:list of 3 integers or None
magnetic

This sets the options of calculating the magnetic scattering instead of nuclear. This assume neutrons are being used.

Getter:Returns bool of magnetic scattering option
Setter:Sets whether magnetic sacttering is calculated
Type:bool
number_of_lots

The number of lots to use

Getter:Returns the number of lots
Setter:Sets the number of lots
Type:int
radiation

The radiation used

Getter:Returns the radiation selected
Setter:Sets the radiation
Type:str (‘xray’ or ‘neutron’)
structure

The structure from which fourier transform is calculated. The calculation work with any of the following types of structures javelin.structure.Structure, ase.Atoms or diffpy.Structure.structure.Structure but if you are using average stucture subtraction or the lots option it needs to be javelin.structure.Structure type.

Getter:Returns the structure
Setter:Sets the structure
Type:javelin.structure.Structure, ase.Atoms, diffpy.Structure.structure.Structure
javelin.fourier.create_xarray_dataarray(values, grid)[source]

Create a xarry DataArray from the input numpy array and grid object.

Parameters:
  • values (numpy.ndarray) – Input array containing the scattering intensities
  • numbers (javelin.grid.Grid) – Grid object describing the array properties
Returns:

DataArray produced from the values and grid object

Return type:

xarray.DataArray

javelin.fourier.get_ff(atomic_number, radiation, q=None)[source]

Returns the form factor for a given atomic number, radiation and q values

Parameters:
  • atomic_number (int) – atomic number
  • radiation (str) – type of radiation (‘xray’ or ‘neutron’)
  • q (float, list, numpy.ndarray) – value or values of q for which to get form factors
Returns:

form factors for given q

Return type:

float, numpy.ndarray

Examples:
>>> get_ff(8, 'neutron')
5.805
>>> get_ff(8, 'xray', q=2.0)
6.31826029176493
>>> get_ff(8, 'xray', q=[0.0, 3.5, 7.0])
array([ 7.999706  ,  4.38417867,  2.08928068])
javelin.fourier.get_mag_ff(atomic_number, q, ion=0, j=0)[source]

Returns the j0 magnetic form factor for a given atomic number, radiation and q values

Parameters:
  • atomic_number (int) – atomic number
  • q (float, list, numpy.ndarray) – value or values of q for which to get form factors
  • ion (int) – charge of selected atom
  • j (int) – order of spherical Bessel function (0, 2, 4 or 6)
Returns:

magnetic form factor for given q

Return type:

float, numpy.ndarray

Examples:
>>> get_mag_ff(8, q=2, ion=1)
0.58510426376585045
>>> get_mag_ff(26, q=[0.0, 3.5, 7.0], ion=2)
array([ 1.        ,  0.49729671,  0.09979243])
>>> get_mag_ff(26, q=[0.0, 3.5, 7.0], ion=4)
array([ 0.9997    ,  0.58273549,  0.13948496])
>>> get_mag_ff(26, q=[0.0, 3.5, 7.0], ion=4, j=4)
array([ 0.       ,  0.0149604,  0.0759222])

grid

class javelin.grid.Grid[source]

Grid class to allow the Q-space grid to be definied in different ways. The grid can be defined be either specifying the corners of the volume or by the axis vectors.

The grid is defined by three vectors v1, v2 and v3; the range of these vectors r1, r2 and r3; and the number of bin in each of these directions.

Examples:

Setting grid by defining corners

>>> grid = Grid()
>>> grid.set_corners(ll=[-2,-2,-3], lr=[2,2,-3], ul=[-2,-2,3])
>>> grid.bins = 5, 3
>>> grid.v1
array([ 1.,  1.,  0.])
>>> grid.v2
array([ 0.,  0.,  1.])
>>> grid.v3
array([ 1., -1.,  0.])
>>> grid.r1
array([-2., -1.,  0.,  1.,  2.])
>>> grid.r2
array([-3.,  0.,  3.])
>>> grid.r3
array([ 0.])
>>> print(grid)
lower left  corner :     [-2. -2. -3.]
lower right corner :     [ 2.  2. -3.]
upper left  corner :     [-2. -2.  3.]
top   left  corner :     [-2. -2. -3.]

hor. increment     :     [ 1.  1.  0.]
vert. increment    :     [ 0.  0.  3.]
top   increment    :     [ 1. -1.  0.]

# of points        :     5 x 3 x 1

Setting grid by vectors and ranges

>>> grid = Grid()
>>> grid.v1 = [1, 1, 1]
>>> grid.v2 = [2, -1, -1]
>>> grid.v3 = [0, 1, -1]
>>> grid.r1 = [-1, 1]
>>> grid.r2 = [0, 2]
>>> grid.r3 = [0, 2]
>>> grid.bins = 3, 3, 3
>>> print(grid)
lower left  corner :     [-1 -1 -1]
lower right corner :     [1 1 1]
upper left  corner :     [ 3 -3 -3]
top   left  corner :     [-1  1 -3]

hor. increment     :     [ 1.  1.  1.]
vert. increment    :     [ 2. -1. -1.]
top   increment    :     [ 0.  1. -1.]

# of points        :     3 x 3 x 3
bins

The number of bins in each direction

>>> grid = Grid()
>>> grid.bins
(101, 101, 1)
>>> grid.bins = 5
>>> grid.bins
(5, 1, 1)
>>> grid.bins = 5, 6
>>> grid.bins
(5, 6, 1)
>>> grid.bins = 5, 6, 7
>>> grid.bins
(5, 6, 7)
Getter:Returns the number of bins in each direction
Setter:Sets the number of bins, provide 1, 2 or 3 integers
Type:numpy.ndarray of int
get_axes_names()[source]
>>> grid = Grid()
>>> grid.get_axes_names()
('[1 0 0]', '[0 1 0]', '[0 0 1]')
Returns:Axis names, vector of each direction
Return type:tuple of str
get_q_meshgrid()[source]

Equivalent to numpy.mgrid for this volume

Returns:mesh-grid numpy.ndarray all of the same dimensions
Return type:tuple of numpy.ndarray
get_squashed_q_meshgrid()[source]

Almost equivalent to numpy.ogrid for this volume. It may have more than one dimension not equal to 1. This can be used with numpy broadcasting.

Returns:mesh-grid numpy.ndarray with some dimension equal to 1
Return type:tuple of numpy.ndarray
ll
Returns:Lower-left corner of reciprocal volume
Return type:numpy.ndarray
lr
Returns:Lower-right corner of reciprocal volume
Return type:numpy.ndarray
r1

Set the range of the first axis, two values min and max

Getter:Array of values for each bin in the axis
Setter:Range of first axis, two values
Type:numpy.ndarray
r2

Set the range of the second axis, two values min and max

Getter:Array of values for each bin in the axis
Setter:Range of second axis, two values
Type:numpy.ndarray
r3

Set the range of the third axis, two values min and max

Getter:Array of values for each bin in the axis
Setter:Range of third axis, two values
Type:numpy.ndarray
set_corners(ll=(0, 0, 0), lr=None, ul=None, tl=None)[source]

Define the axis vectors by the corners of the reciprocal volume. The corners values with be converted into axis vectors, see javelin.grid.corners_to_vectors() for details.

Parameters:
  • ll (array-like of length 3) – lower-left corner
  • lr (array-like of length 3) – lower-right corner
  • ul (array-like of length 3) – upper-left corner
  • tl (array-like of length 3) – top-left corner
tl
Returns:Top-left corner of reciprocal volume
Return type:numpy.ndarray
ul
Returns:Upper-left corner of reciprocal volume
Return type:numpy.ndarray
v1
Getter:Set the first axis
Setter:Get the first axis
Type:numpy.ndarray
v1_delta
Returns:Increment vector of first axis
Return type:numpy.ndarray
v2
Getter:Set the second axis
Setter:Get the second axis
Type:numpy.ndarray
v2_delta
Returns:Increment vector of second axis
Return type:numpy.ndarray
v3
Getter:Set the third axis
Setter:Get the third axis
Type:numpy.ndarray
v3_delta
Returns:Increment vector of third axis
Return type:numpy.ndarray
javelin.grid.angle(v1, v2)[source]

Calculates the angle between two vectors

Parameters:
  • v1 (array-like object of numbers) – vector 1
  • v2 (array-like object of numbers) – vector 2
Returns:

angle (radians)

Return type:

float

Examples:
>>> angle([1,0,0], [-1,0,0])
3.1415926535897931
>>> angle([1,0,0], [0,1,0])
1.5707963267948966
javelin.grid.check_parallel(v1, v2)[source]

Checks if two vectors are parallel

Parameters:
  • v1 (array-like object of numbers) – vector1
  • v2 (array-like object of numbers) – vector2
Returns:

if parallel

Return type:

bool

Examples:
>>> check_parallel([1,0,0], [-1,0,0])
True
>>> check_parallel([1,0,0], [0,1,0])
False
javelin.grid.corners_to_vectors(ll=None, lr=None, ul=None, tl=None)[source]

This function converts the provided corners into axes vectors and axes ranges. It will also calculate sensible vectors for any unprovided corners.

You must provide at minimum the lower-left (ll) and lower-right (lr) corners.

Parameters:
  • ll (array-like object of numbers) – lower-left corner (required)
  • lr (array-like object of numbers) – lower-right corner (requried)
  • ul (array-like object of numbers) – upper-left corner
  • tl (array-like object of numbers) – top-left corner
Returns:

three axes vectors and three axes ranges

Return type:

tuple of three numpy.ndarray and three tuple ranges

Examples:

Using only ll and lr, the other two vector are calculated using javelin.grid.find_other_vectors()

>>> v1, v2, v3, r1, r2, r3 = corners_to_vectors(ll=[-3,-3,0], lr=[3, 3, 0])
>>> print(v1, v2, v3)
[ 1.  1.  0.] [ 1.  0.  0.] [ 0.  0.  1.]
>>> print(r1, r2, r3)
(-3.0, 3.0) (0.0, 0.0) (0.0, 0.0)

Using ll, lr and ul, the other vector is the javelin.grid.norm1() of the cross product of the first two vectors defined by the corners.

>>> v1, v2, v3, r1, r2, r3 = corners_to_vectors(ll=[-3,-3,-2], lr=[3, 3, -2], ul=[-3, -3, 2])
>>> print(v1, v2, v3)
[ 1.  1.  0.] [ 0.  0.  1.] [ 1. -1.  0.]
>>> print(r1, r2, r3)
(-3.0, 3.0) (-2.0, 2.0) (0.0, 0.0)

Finally defining all corners

>>> v1,v2,v3,r1,r2,r3 = corners_to_vectors(ll=[-5,-6,-7],lr=[-5,-6,7],ul=[-5,6,-7],tl=[5,-6,-7])
>>> print(v1, v2, v3)
[ 0.  0.  1.] [ 0.  1.  0.] [ 1.  0.  0.]
>>> print(r1, r2, r3)
(-7.0, 7.0) (-6.0, 6.0) (-5.0, 5.0)

If you provided corners which will create parallel vectors you will get a ValueError

>>> corners_to_vectors(ll=[0, 0, 0], lr=[1, 0, 0], ul=[1, 0, 0])
Traceback (most recent call last):
    ...
ValueError: Vector from ll to lr is parallel with vector from ll to ul
javelin.grid.find_other_vectors(v)[source]

This will find two new vectors which in combination with the provided vector (v) will form a basis for a complete space filling set.

Parameters:v (array-like object of numbers) – vector
Returns:two new space filling vectors
Return type:tuple of two numpy.ndarray
Examples:
>>> find_other_vectors([1, 0, 0])
(array([ 0.,  1.,  0.]), array([ 0.,  0.,  1.]))
>>> find_other_vectors([0, 0, 1])
(array([ 1.,  0.,  0.]), array([ 0.,  1.,  0.]))
>>> find_other_vectors([1, 1, 0])
(array([ 1.,  0.,  0.]), array([ 0.,  0.,  1.]))
javelin.grid.get_vector_from_points(p1, p2)[source]

Calculates the vector form two points

Parameters:
  • p1 (array-like object of numbers) – point 1
  • p2 (array-like object of numbers) – point 2
Returns:

vector between points

Return type:

numpy.ndarray

Examples:
>>> get_vector_from_points([-1, -1, 0], [1, 1, 0])
array([ 1.,  1.,  0.])
>>> get_vector_from_points([0, 0, 0], [2, 2, 4])
array([ 1.,  1.,  2.])
javelin.grid.length(v)[source]

Calculates the length of a vector

Parameters:v (array-like object of numbers) – vector
Returns:length of vector
Return type:float
Examples:
>>> length([1,0,0])
1.0
>>> length([1,-1,0])
1.4142135623730951
>>> length([2,2,2])
3.4641016151377544
javelin.grid.norm(v)[source]

Calculates the normalised vector

Parameters:v (array-like object of numbers) – vector
Returns:normalised vector
Return type:numpy.ndarray
Examples:
>>> norm([5, 0, 0])
array([ 1.,  0.,  0.])
>>> norm([1, 1, 0])
array([ 0.70710678,  0.70710678,  0.        ])
javelin.grid.norm1(v)[source]

Calculate the equivalent vector with the smallest non-zero component equal to one.

Parameters:v (array-like object of numbers) – vector
Returns:normalised1 vector
Return type:numpy.ndarray
Examples:
>>> norm1([5, 10, 0])
array([ 1.,  2.,  0.])
>>> norm1([1, 1, 0])
array([ 1.,  1.,  0.])

io

javelin.io.load_HDF5_to_xarray(filename)[source]

Load HDF file into an xarray DataArray using pandas HDFStore

requries pytables

javelin.io.numpy_to_vti(array, origin, spacing, filename)[source]

This function write a VtkImageData vti file from a numpy array.

Parameters:
  • array (numpy.ndarray) – input array
  • origin (array like object of values) – the origin of the array
  • spacing (array like object of values) – the step in each dimension
  • filename (str) – output filename (.vti)
javelin.io.read_mantid_MDHisto(filename)[source]

Read the saved MDHisto from from Mantid and returns an xarray.DataArray object

javelin.io.read_stru(filename, starting_cell=(1, 1, 1))[source]

Read in a .stru file saved from DISCUS into a javelin Structure

If the line ncell is not present in the file all the atoms will be read into a single cell.

javelin.io.read_stru_to_ase(filename)[source]

This function read the legacy DISCUS stru file format into a ASE Atoms object.

Parameters:filename (str) – filename of DISCUS stru file
Returns:ASE Atoms object
Return type:ase.Atoms
javelin.io.save_mantid_MDHisto(dataArray, filename)[source]

Save a file that can be read in using Mantid’s LoadMD

javelin.io.save_xarray_to_HDF5(dataArray, filename, complib=None)[source]

Save the xarray DataArray to HDF file using pandas HDFStore

attrs will be saved as metadata via pickle

requries pytables

complib : {‘zlib’, ‘bzip2’, ‘lzo’, ‘blosc’, None}, default None

structure

class javelin.structure.Structure(symbols=None, numbers=None, unitcell=1, ncells=None, positions=None, rotations=False, translations=False, magnetic_moments=False)[source]

The structure class is made up of a unitcell and a list of atoms

Structure can be initialize using either another javelin.structure.Structure, ase.Atoms or diffpy.Structure.structure.Structure. It is recommended you use javelin.structure.Structure.reindex() after initializing from a foreign type in order to get the correct unitcell structure type.

Parameters:
  • symbols (list) – atoms symbols to initialize structure
  • numbers (list) – atomic numbers to initialize structure
  • unitcell (javelin.unitcell.UnitCell) – unitcell of structure, can be javelin.unitcell.UnitCell or values used to initialize the UnitCell
  • ncells (list) – ncells has four components, (i, j, k, n) where i, j, k are the number of unitcell in each direction and n is the number of site positions in each unitcell. The product of ncells must equal the total number of atoms in the structure.
  • positions (3 x n array-like) – array of atom coordinates
add_atom(i=0, j=0, k=0, site=0, Z=None, symbol=None, position=None)[source]

Adds a single atom to the structure. It the atom exist as provided i, j, k and site it will be replaced.

Parameters:
  • i (int) – unitcell index
  • j (int) – unitcell index
  • k (int) – unitcell index
  • site (int) – site index
  • Z (int) – atomic number
  • symbol (int) – chemical symbol
  • position (vector) – position within the unitcell
>>> stru=Structure(unitcell=5.64)
>>> stru.atoms
Empty DataFrame
Columns: [Z, symbol, x, y, z]
Index: []
>>> stru.add_atom(Z=12, position=[0,0,0])
>>> stru.atoms 
             Z symbol  x  y  z
i j k site
0 0 0 0     12     Mg  0  0  0
>>> stru.add_atom(Z=13, position=[0.5,0,0], i=1)
>>> stru.atoms 
             Z symbol    x  y  z
i j k site
0 0 0 0     12     Mg    0  0  0
1 0 0 0     13     Al  0.5  0  0
atoms = None

Attribute storing list of atom type and positions as a pandas.DataFrame

Example:
>>> stru = Structure(symbols=['Na','Cl','Na'],positions=[[0,0,0],[0.5,0.5,0.5],[0,1,0]])
>>> stru.atoms
             Z symbol    x    y    z
i j k site
0 0 0 0     11     Na  0.0  0.0  0.0
      1     17     Cl  0.5  0.5  0.5
      2     11     Na  0.0  1.0  0.0
element

Array of all elements in the stucture

Returns:array of element symbols
Return type:numpy.ndarray
Example:
>>> stru = Structure(symbols=['Na','Cl','Na'],positions=[[0,0,0],[0.5,0.5,0.5],[0,1,0]])
>>> stru.element
array(['Na', 'Cl', 'Na'], dtype=object)
get_atom_Zs()[source]

Get a list of unique atomic number in structure

Returns:array of Zs
Return type:numpy.ndarray
Example:
>>> stru = Structure(symbols=['Na','Cl','Na'],positions=[[0,0,0],[0.5,0.5,0.5],[0,1,0]])
>>> stru.get_atom_Zs()
array([11, 17])
get_atom_count()[source]

Returns a count of each different type of atom in the structure

Returns:series of atom count
Return type:pandas.Series
Example:
>>> stru = Structure(symbols=['Na','Na'],positions=[[0,0,0],[0.5,0.5,0.5]])
>>> stru.get_atom_count()
Na    2
Name: symbol, dtype: int64
get_atom_symbols()[source]

Get a list of unique atom symbols in structure

Returns:array of atom symbols
Return type:numpy.ndarray
Example:
>>> stru = Structure(symbols=['Na','Cl','Na'],positions=[[0,0,0],[0.5,0.5,0.5],[0,1,0]])
>>> stru.get_atom_symbols()
array(['Na', 'Cl'], dtype=object)
get_atomic_numbers()[source]

Array of all atomic numbers in the stucture

Returns:array of atomic numbers
Return type:numpy.ndarray
Example:
>>> stru = Structure(symbols=['Na','Cl','Na'],positions=[[0,0,0],[0.5,0.5,0.5],[0,1,0]])
>>> stru.get_atomic_numbers()
array([11, 17, 11])
get_chemical_formula()[source]

Returns the chemical formula of the structure

Returns:chemical formula
Return type:str
Example:
>>> stru = Structure(symbols=['Na','Cl','Na'],positions=[[0,0,0],[0.5,0.5,0.5],[0,1,0]])
>>> stru.get_chemical_formula()
'Na2Cl1'
get_chemical_symbols()[source]

Same as javelin.structure.Structure.element

get_magnetic_moments()[source]
get_positions()[source]

Same as javelin.structure.Structure.xyz_cartn

get_scaled_positions()[source]

Array of all xyz positions in fractional lattice units of the atoms in the structure

Returns:3 x n array of atom positions
Return type:numpy.ndarray
Example:
>>> stru = Structure(symbols=['Na', 'Cl'], positions=[[0,0,0],[0.5,0.5,0.5]], unitcell=5.64)
>>> stru.get_scaled_positions()
array([[ 0. ,  0. ,  0. ],
       [ 0.5,  0.5,  0.5]])
number_of_atoms

The total number of atoms in the structure

Returns:number of atoms in structure
Return type:int
Example:
>>> stru = Structure(symbols=['Na','Cl','Na'],positions=[[0,0,0],[0.5,0.5,0.5],[0,1,0]])
>>> stru.number_of_atoms
3
rattle(scale=0.001, seed=None)[source]

Randomly move all atoms by a normal distbution with a standard deviation given by scale.

Parameters:
  • scale (float) – standard deviation
  • seed (int) – seed for random number generator
Example:
>>> stru = Structure(symbols=['Na', 'Cl'], positions=[[0,0,0],[0.5,0.5,0.5]], unitcell=5.64)
>>> print(stru.xyz)
[[ 0.   0.   0. ]
 [ 0.5  0.5  0.5]]
>>> stru.rattle(seed=42)
>>> print(stru.xyz)
[[  4.96714153e-04  -1.38264301e-04   6.47688538e-04]
 [  5.01523030e-01   4.99765847e-01   4.99765863e-01]]
reindex(ncells)[source]

This will reindex the list of atoms into the unitcell framework of this structure

ncells has four components, (i, j, k, n) where i, j, k are the number of unitcell in each direction and n is the number of site positions in each unitcell. The product of ncells must equal the total number of atoms in the structure.

Example:
>>> stru = Structure(symbols=['Na', 'Cl'], positions=[[0,0,0],[0.5,0.5,0.5]], unitcell=5.64)
>>> stru.atoms 
             Z symbol    x    y    z
i j k site
0 0 0 0     11     Na  0.0  0.0  0.0
      1     17     Cl  0.5  0.5  0.5
>>> stru.reindex([2,1,1,1])
>>> stru.atoms 
             Z symbol    x    y    z
i j k site
0 0 0 0     11     Na  0.0  0.0  0.0
1 0 0 0     17     Cl  0.5  0.5  0.5
repeat(rep)[source]

Repeat the cells a number of time along each dimension

rep argument should be either three value like (1,2,3) or a single value r equivalent to (r,r,r).

Parameters:rep (1 or 3 ints) – repeating rate
Examples:
>>> stru = Structure(symbols=['Na', 'Cl'], positions=[[0,0,0],[0.5,0.5,0.5]], unitcell=5.64)
>>> print(stru.element)
['Na' 'Cl']
>>> print(stru.xyz_cartn)
[[ 0.    0.    0.  ]
 [ 2.82  2.82  2.82]]
>>> stru.repeat((2,1,1))
>>> print([str(e) for e in stru.element])
['Na', 'Cl', 'Na', 'Cl']
>>> print(stru.xyz_cartn)
[[  0.00000000e+00   0.00000000e+00   0.00000000e+00]
 [  2.82000000e+00   2.82000000e+00   2.82000000e+00]
 [  5.64000000e+00   9.06981174e-16   9.06981174e-16]
 [  8.46000000e+00   2.82000000e+00   2.82000000e+00]]
>>> stru = Structure(symbols=['Na'], positions=[[0,0,0]], unitcell=5.64)
>>> print(stru.element)
['Na']
>>> print(stru.xyz_cartn)
[[ 0.  0.  0.]]
>>> stru.repeat(2)
>>> print([str(e) for e in stru.element])
['Na', 'Na', 'Na', 'Na', 'Na', 'Na', 'Na', 'Na']
>>> print(stru.xyz_cartn)
[[  0.00000000e+00   0.00000000e+00   0.00000000e+00]
 [  0.00000000e+00   0.00000000e+00   5.64000000e+00]
 [  0.00000000e+00   5.64000000e+00   3.45350397e-16]
 [  0.00000000e+00   5.64000000e+00   5.64000000e+00]
 [  5.64000000e+00   9.06981174e-16   9.06981174e-16]
 [  5.64000000e+00   9.06981174e-16   5.64000000e+00]
 [  5.64000000e+00   5.64000000e+00   1.25233157e-15]
 [  5.64000000e+00   5.64000000e+00   5.64000000e+00]]
unitcell = None

Attribute containing the unitcell of the structure. Must be of type javelin.unitcell.UnitCell

x

Array of all x positions of in fractional lattice units of the atoms within the unitcell of the structure

Returns:array of atom x positions
Return type:numpy.ndarray
Example:
>>> stru = Structure(symbols=['Na', 'Cl'], positions=[[0,0,0],[0.5,0.5,0.5]], unitcell=5.64)
>>> stru.x
array([ 0. ,  0.5])
xyz

Array of all xyz positions in fractional lattice units of the atoms within the unitcell of the structure

Returns:3 x n array of atom positions
Return type:numpy.ndarray
Example:
>>> stru = Structure(symbols=['Na', 'Cl'], positions=[[0,0,0],[0.5,0.5,0.5]], unitcell=5.64)
>>> stru.xyz
array([[ 0. ,  0. ,  0. ],
       [ 0.5,  0.5,  0.5]])
xyz_cartn

Array of all xyz positions in cartesian coordinates of the atoms in the structure

Returns:3 x n array of atom positions
Return type:numpy.ndarray
Example:
>>> stru = Structure(symbols=['Na', 'Cl'], positions=[[0,0,0],[0.5,0.5,0.5]], unitcell=5.64)
>>> stru.xyz_cartn
array([[ 0.  ,  0.  ,  0.  ],
       [ 2.82,  2.82,  2.82]])
y

Array of all y positions of in fractional lattice units of the atoms within the unitcell of the structure

Returns:array of atom y positions
Return type:numpy.ndarray
Example:
>>> stru = Structure(symbols=['Na', 'Cl'], positions=[[0,0,0],[0.5,0.5,0.5]], unitcell=5.64)
>>> stru.y
array([ 0. ,  0.5])
z

Array of all z positions of in fractional lattice units of the atoms within the unitcell of the structure

Returns:array of atom z positions
Return type:numpy.ndarray
Example:
>>> stru = Structure(symbols=['Na', 'Cl'], positions=[[0,0,0],[0.5,0.5,0.5]], unitcell=5.64)
>>> stru.z
array([ 0. ,  0.5])
javelin.structure.axisAngle2Versor(x, y, z, angle, unit='degrees')[source]
javelin.structure.get_miindex(l=0, ncells=None)[source]
javelin.structure.get_rotation_matrix(l, m, n, theta, unit='degrees')[source]
javelin.structure.get_rotation_matrix_from_versor(w, x, y, z)[source]

unitcell

class javelin.unitcell.UnitCell(*args)[source]

The UnitCell ojbect can be set with either 1, 3 or 6 parameters corresponding to cubic a parameters, (a, b, c) or (a, b, c, alpha, beta, gamma), where angles are in degrees.

>>> cubic = UnitCell(5)
>>> cubic.cell
(5.0, 5.0, 5.0, 90.0, 90.0, 90.0)
>>> orthorhombic = UnitCell(5, 6, 7)
>>> orthorhombic.cell
(5.0, 6.0, 7.0, 90.0, 90.0, 90.0)
>>> unitcell = UnitCell(4.0, 3.0, 6.0, 89.0, 90.0, 97.0)
>>> unitcell.cell
(4.0, 3.0, 6.0, 89.0, 90.0, 97.0)

UnitCell objects can be set after being created simply by

>>> unitcell = UnitCell()
>>> unitcell.cell = 6
>>> unitcell.cell
(6.0, 6.0, 6.0, 90.0, 90.0, 90.0)
>>> unitcell.cell = 3, 4, 5
>>> unitcell.cell
(3.0, 4.0, 5.0, 90.0, 90.0, 90.0)
>>> unitcell.cell = 6, 7, 8, 91.0, 90, 89
>>> unitcell.cell
(6.0, 7.0, 8.0, 91.0, 90.0, 89.0)
>>> # or using a list or tuple
>>> unitcell.cell = [8, 7, 6, 89, 90, 90]
>>> unitcell.cell
(8.0, 7.0, 6.0, 89.0, 90.0, 90.0)
B

Returns the B matrix

Binv

Returns the inverse B matrix

G

Returns the metric tensor G

Gstar

Returns the reciprocal metric tensor G*

cartesian(u)[source]

Return Cartesian coordinates of a lattice vector.

>>> unitcell = UnitCell(3,4,5,90,90,120)
>>> unitcell.cartesian([1,0,0])
array([[  2.59807621e+00,  -1.50000000e+00,   3.25954010e-16]])

A array of atoms position can also be passed

>>> positions = [[1,0,0], [0,0,0.5]]
>>> unitcell.cartesian(positions)
array([[  2.59807621e+00,  -1.50000000e+00,   3.25954010e-16],
       [  0.00000000e+00,   0.00000000e+00,   2.50000000e+00]])
cell

Return the unit cell parameters (a, b, c, alpha, beta, gamma) in degrees.

d(h, k, l)[source]

Returns d-spacing for given h,k,l

dstar(h, k, l)[source]

Returns d*=1/d for given h,k,l

fractional(u)[source]

Return Cartesian coordinates of a lattice vector.

>>> unitcell = UnitCell(3,4,5,90,90,120)
>>> unitcell.fractional([0,4,0])
array([[  0.00000000e+00,   1.00000000e+00,  -4.89858720e-17]])

A array of atoms position can also be passed

>>> positions = [[0,2,0], [0,0,5]]
>>> unitcell.fractional(positions)
array([[  0.00000000e+00,   5.00000000e-01,  -2.44929360e-17],
       [  0.00000000e+00,   0.00000000e+00,   1.00000000e+00]])
recAngle(h1, k1, l1, h2, k2, l2, degrees=False)[source]

Calculates the angle between two reciprocal vectors

reciprocalCell

Return the reciprocal unit cell parameters (a*, b*, c*, alpha*, beta*, gamma*) in degrees.

reciprocalVolume

Returns the unit cell reciprocal volume

volume

Returns the unit cell volume

utils

javelin.utils.get_atomic_number_symbol(Z=None, symbol=None)[source]

This function returns a tuple of matching arrays of atomic numbers (Z) and chemical symbols (symbol).

Parameters:
  • Z (int, array like object of int's) – atomic numbers
  • symbol (str, array like object of str) – chemical symbols
Returns:

arrays of atomic numbers and chemical symbols

Return type:

tuple of numpy.ndarray

Note: If both Z and symbol are provided the symbol will win out and change the Z to match.

Examples:
>>> Z, symbol = get_atomic_number_symbol(Z=[12, 24, 26, 48])
>>> print(Z)
[12 24 26 48]
>>> print(symbol[0])
Mg
>>> print(symbol[1])
Cr
>>> print(symbol[2])
Fe
>>> print(symbol[3])
Cd
>>> Z, symbol = get_atomic_number_symbol(symbol=['C', 'H', 'N', 'O'])
>>> print(Z)
[6 1 7 8]
>>> print(symbol)
['C' 'H' 'N' 'O']
javelin.utils.get_atomic_numbers(structure)[source]

Wrapper to get the atomic numbers from different structure classes

javelin.utils.get_positions(structure)[source]

Wrapper to get the positions from different structure classes

javelin.utils.get_unitcell(structure)[source]

Wrapper to get the unit cell from different structure classes

javelin.utils.is_structure(structure)[source]

Check if an object is a stucture that javelin can understand.

ase.atoms with have cell, get_scaled_positions and get_atomic_numbers attributes diffpy.structure with have lattice, xyz, and element attributes