Crystals

Xtals.jl maintains a data structure Crystal that stores information about a crystal structure file.

Reading in a Crystal Structure File

Currently, the crystal structure file reader accepts .cif and .cssr file formats. Xtals.jl looks for the crystal structure files in rc[:paths][:crystals] which is by default ./data/crystals (relative to pwd() at module loading). By typing rc[:paths][:crystals] = "my_crystal_dir", Xtals.jl now looks for the crystal structure file in my_crystal_dir. The files can be read as:

xtal = Crystal("IRMOF-1.cif")       # The crystal reader stores the information in xtal
xtal.name                           # The name of the crystal structure file
xtal.box                            # The unit cell information
xtal.atoms                          # The atom coordinates (in fractional space) and the atom identities
xtal.charges                        # The charge magnitude and coordinates (in fractional space)
xtal.bonds                          # Bonding information in the structure. By default this is an empty graph,
#  but use `read_bonds_from_file=true` argument in `Crystal` to read from crystal structure file
xtal.symmetry                       # Symmetry information of the crystal. By default converts the symmetry to P1 symmetry.
#  Use `convert_to_p1=false` argument in `Crystal` to keep original symmetry

Fixing Atom Species Labels

Often, the atoms species are appended by numbers. This messes with the internal workings of Xtals.jl. To circumvent this problem, the function strip_numbers_from_atom_labels! removes the appended numbers. It is important to use this function prior to GCMC or Henry coefficient calculations and bond inference operations.

xtal = Crystal("IRMOF-1.cif"; species_col=["_atom_site_label"])
xtal.atoms.species[1]

# output

:Zn1
strip_numbers_from_atom_labels!(xtal)
xtal.atoms.species[1]

# output

:Zn

Converting the Coordinates to Cartesian Space

The coordinates of the Crystal's Atoms are stored in Fractional coordinates. If one needs to analyze the Cartesian coordinates of the Crystal, that can be done by using the unit cell (Box) information.

xtal.atoms.coords.xf                                    # array of fractional coordinates
cart_coords = xtal.box.f_to_c * xtal.atoms.coords.xf    # array of cartesian coordinates

Creating a Super-Cell

For many simulations, one needs to replicate the unit cell multiple times to create a bigger super-cell. This is done with replicate:

super_xtal = replicate(xtal, (2, 2, 2))       # Replicates the original unit cell once in each dimension
xtal.atoms.n, super_xtal.atoms.n

# output

(424, 3392)

Finding Other Properties

rho = crystal_density(xtal)         # Crystal density of the crystal in kg/m^2
mw = molecular_weight(xtal)         # The molecular weight of the unit cell in amu
formula = empirical_formula(xtal)   # The irreducible chemical formula of the crystal

Assigning New Charges

If the structure file does not contain partial charges, we provide methods to assign new Charges to the Crystal.

species_to_charge = Dict(:Zn => 2.0, :C => 0.0, :H => 0.0, :O => -0.61538)  # This method assigns a static charge to atom species
charged_xtal = assign_charges(xtal, species_to_charge, 1e-3)                # This function creates a new charged `Crystal` object.
#   The function checks for charge neutrality with a tolerance of 1e-3
new_charges = Charges(zeros(xtal.atoms.n), xtal.atoms.coords)
other_charged_xtal = Crystal(
    xtal.name,
    xtal.box,
    xtal.atoms,               # Here we create a new `Charges` object using an array of new charges.
    new_charges,
    xtal.bonds,
    xtal.symmetry
)        #   The number of charges in the array has to be equal to the number of atoms
#   and finally a new `Crystal` object is manually created

Writing Crystal Files

We provide methods to write both .xyz and .cif files.

write_cif(xtal, "my_new_cif_file.cif")      # Stored in the current directory
write_xyz(xtal, "my_new_xyz_file.xyz")      # stored in the current directory

Detailed Docs

Xtals.CrystalType
crystal = Crystal(filename;
    check_neutrality=true, net_charge_tol=1e-4,
    check_overlap=true, convert_to_p1=true,
    read_bonds_from_file=false, wrap_coords=true,
    include_zero_charges=true,
    remove_duplicates=false,
    species_col=["_atom_site_type_symbol", "_atom_site_label"]
    ) # read from file

crystal = Crystal(name, box, atoms, charges) # construct from matter, no bonds, P1-symmetry assumed

Read a crystal structure file (.cif or .cssr) and populate a Crystal data structure, or construct a Crystal data structure directly.

Arguments

  • filename::AbstractString: the name of the crystal structure file (include ".cif" or ".cssr") read from PATH_TO_CRYSTALS.
  • check_neutrality::Bool: check for charge neutrality
  • net_charge_tol::Float64: when checking for charge neutrality, throw an error if the absolute value of the net charge is larger than this value.
  • check_overlap::Bool: throw an error if overlapping atoms are detected.
  • convert_to_p1::Bool: If the structure is not in P1 it will be converted to P1 symmetry using the symmetry rules from the _symmetry_equiv_pos_as_xyz list in the .cif file. (We do not use the space groups name to look up symmetry rules).
  • read_bonds_from_file::Bool: Whether or not to read bonding information from cif file. If false, the bonds can be inferred later. note that, if the crystal is not in P1 symmetry, we cannot both read bonds and convert to P1 symmetry.
  • wrap_coords::Bool: if true, enforce that fractional coords of atoms and charges are in [0,1]³ by mod(x, 1)
  • include_zero_charges::Bool: if false, do not include in crystal.charges atoms which have zero charges, in order to speed up the electrostatic calculations. If true, include the atoms in crystal.charges that have zero charge, ensuring that the number of atoms is equal to the number of charges and that crystal.charges.coords.xf and crystal.atoms.coords.xf are the same.
  • remove_duplicates::Bool: remove duplicate atoms and charges. an atom is duplicate only if it is the same species.
  • species_col::Array{String}: which column to use for species identification for crystal.atoms.species. we use a priority list: we check for the first entry of species_col in the .cif file; if not present, we then use the second entry, and so on.
  • infer_bonds::Bool: if true, bonds are inferred automatically. If set, must specify periodic_boundaries. By default, bonds are not inferred.
  • periodic_boundaries::Union{Bool, Missing}: use with infer_bonds to specify treatment of the unit cell boundary. Set true to treat the unit cell edge as a periodic boundary (allow bonds across it); set false to restrict bonding to within the local unit cell. By default, no bonds are inferred.
  • absolute_path::Bool: set true to load from an absolute path instead of the default path stored in rc[:paths][:crystals].

Returns

  • crystal::Crystal: A crystal containing the crystal structure information

Attributes

  • name::AbstractString: name of crystal structure
  • box::Box: unit cell (Bravais Lattice)
  • atoms::Atoms: list of Atoms in crystal unit cell
  • charges::Charges: list of point charges in crystal unit cell
  • bonds::MetaGraph: Unweighted, undirected graph showing all of the atoms that are bonded within the crystal
  • symmetry::SymmetryInfo: symmetry inforomation
source
Xtals.SymmetryInfoType
SymmetryInfo(symmetry, space_group, is_p1)

Attributes

  • operations::Array{Function, 2}: 2D array of anonymous functions that represent the symmetry operations. If the structure is in P1 there will be one symmetry operation.
  • space_group::AbstractString: The name of the space group. This is stored so that it can be written out again in the write_cif function. The space group is not used to verify the symmetry rules.
  • is_p1::Bool: Stores whether the crystal is currently in P1 symmetry.
source
Xtals.strip_numbers_from_atom_labels!Function
strip_numbers_from_atom_labels!(crystal)

Strip numbers from labels for crystal.atoms. Precisely, for atom in crystal.atoms, find the first number that appears in atom. Remove this number and all following characters from atom. e.g. C12 –> C Ba12A_3 –> Ba

Arguments

  • crystal::Crystal: The crystal containing the crystal structure information
source
Xtals.replicateFunction
new_box = replicate(original_box, repfactors)

Replicates a Box in positive directions to construct a new Box representing a supercell. The original_box is replicated according to the factors in repfactors. Note replicate(original_box, repfactors=(1, 1, 1)) returns same Box. The new fractional coordinates as described by f_to_c and c_to_f still ∈ [0, 1].

Arguments

  • original_box::Box: The box that you want to replicate
  • repfactors::Tuple{Int, Int, Int}: The factor you want to replicate the box by

Returns

  • box::Box: Fully formed Box object
source
replicated_crystal = replicate(crystal, repfactors)

replicate the atoms and charges in a Crystal in positive directions to construct a new Crystal. Note replicate(crystal, (1, 1, 1)) returns the same Crystal. the fractional coordinates will be rescaled to be in [0, 1].

arguments

  • crystal::Crystal: The crystal to replicate
  • repfactors::Tuple{Int, Int, Int}: The factors by which to replicate the crystal structure in each crystallographic direction (a, b, c).

returns

  • replicated_frame::Crystal: replicated crystal
source
Xtals.molecular_weightFunction
mass_of_crystal = molecular_weight(crystal)

Calculates the molecular weight of a unit cell of the crystal in amu using information stored in data/atomicmasses.csv.

Arguments

  • crystal::Crystal: The crystal containing the crystal structure information

Returns

  • mass_of_crystal::Float64: The molecular weight of a unit cell of the crystal in amu
source
Xtals.crystal_densityFunction
ρ = crystal_density(crystal) # kg/m²

Compute the crystal density of a crystal. Pulls atomic masses from read_atomic_masses.

Arguments

  • crystal::Crystal: The crystal containing the crystal structure information

Returns

  • ρ::Float64: The crystal density of a crystal in kg/m³ # –> kg/m3
source
Xtals.empirical_formulaFunction
formula = empirical_formula(crystal, verbose=false)

Find the irreducible chemical formula of a crystal structure.

Arguments

  • crystal::Crystal: The crystal containing the crystal structure information
  • verbose::Bool: If true, will print the chemical formula as well

Returns

  • formula::Dict{Symbol, Int}: A dictionary with the irreducible chemical formula of a crystal structure
source
Xtals.assign_chargesFunction
crystal_with_charges = assign_charges(crystal, species_to_charge, net_charge_tol=1e-5)

assign charges to the atoms present in the crystal based on atom type. pass a dictionary species_to_charge that maps atomic species to a charge.

if the crystal already has charges, the charges are removed and new charges are added. a warning is thrown if this is the case.

checks for charge neutrality in the end.

returns a new crystal.

Examples

species_to_charge = Dict(:Ca => 2.0, :C => 1.0, :H => -1.0)
crystal_with_charges = assign_charges(crystal, species_to_charge, 1e-7)
crystal_with_charges = assign_charges(crystal, species_to_charge) # tol 1e-5 default

Arguments

  • crystal::Crystal: the crystal
  • species_to_charge::Dict{Symbol, Float64}: a dictionary that maps atomic species to charge
  • net_charge_tol::Float64: the net charge tolerated when asserting charge neutrality of the resulting crystal
source
Xtals.write_cifFunction
write_cif(crystal, filename; fractional_coords=true, number_atoms=true)
write_cif(crystal) # writes to file crystal.name[.cif]

Write a crystal::Crystal to a .cif file.

arguments

  • crystal::Crystal: crystal to write to file
  • filename::AbstractString: the filename of the .cif file. if ".cif" is not included as an extension, it will automatically be appended to the filename string.
  • fractional_coords::Bool=true: write the coordinates of the atoms as fractional coords if true. if false, write Cartesian coords.
  • number_atoms::Bool=true: write the atoms as "C1", "C2", "C3", ..., "N1", "N2", ... etc. to give each atom a unique identifier
source
Xtals.write_xyzFunction
write_xyz(atoms, filename; comment="")
write_xyz(crystal; comment="", center_at_origin=false)
write_xyz(molecules, box, filename; comment="") # fractional
write_xyz(molecules, box, filename; comment="") # Cartesian

write atoms to an .xyz file.

Arguments

  • atoms::Atoms: the set of atoms.
  • filename::AbstractString: the filename (absolute path) of the .xyz file. (".xyz" appended automatically if the extension is not provided.)
  • comment::AbstractString: comment if you'd like to write to the file.
  • center_at_origin::Bool: (for crystal only) if true, translate all coords such that the origin is the center of the unit cell.
source
Xtals.read_xyzFunction
atoms = read_xyz("molecule.xyz")

read a list of atomic species and their corresponding coordinates from an .xyz file.

Arguments

  • filename::AbstractString: the path to and filename of the .xyz file

Returns

  • atoms::Atoms{Cart}: the set of atoms read from the .xyz file.
source
Xtals.read_molFunction
atoms, bonds = read_mol("molecule.mol")

read a .mol file, which contains info about both atoms and bonds. see here for the anatomy of a .mol file.

Arguments

  • filename::AbstractString: the path to and filename of the .mol file (must pass extension)

Returns

  • atoms::Atoms{Cart}: the set of atoms read from the .mol file.
  • bonds::MetaGraph: the bonding graph of the atoms read from the .mol file.
  • bond_types::Array{Int, 1}: the array of bond types.
source