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 Xtals.PATH_TO_CRYSTALS which is by default ./data/crystals. By typing set_path_to_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

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!(xtal) removes the appending numbers. It is important to use this function prior to GCMC or Henry coefficient calculations.

xtal.atoms.species              # [:C1, :C2, :O1, ...]
strip_numbers_from_atom_labels!(xtal)
xtal.atoms.species              # [:C, :C, :O, ...]

converting the coordinates to cartesian space

The coordinates of the crystals 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 information of the crystal.

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.

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

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 = chemical_formula(xtal)    # The irreducible chemical formula of the crystal

assigning new charges

If the crystal structure file does not contains partial charges, we provide methods to assign new charges to the crystal

species_to_charges = Dict(:Ca => 2.0, :C => 1.0, :H => -1.0)                # This method assigns a static charge to atom species
charged_xtal = assign_charges(xtal, species_to_charge, 1e-5)                # This function creates a new charged `Crystal` object.
                                                                            #   The function checks for charge neutrality with a tolerance of 1e-5
new_charges = Charges([2.0, 1.0, -1.0, -1.0, ...], 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, overlap_tol=0.1,
    convert_to_p1=true, read_bonds_from_file=false, wrap_coords=true,
    include_zero_charges=false,
    remove_duplicates=false,
    species_col=["_atom_site_label", "_atom_site_type_symbol"]
    ) # 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::String: 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::Union{Bool, Missing} : if true, bonds are inferred across periodic unit cell boundaries; if false, bonds are only inferred within the local unit cell; if missing (default), bonds are not inferred.

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::SimpleGraph: Unweighted, undirected graph showing all of the atoms that are bonded within the crystal
  • symmetry::SymmetryInfo: symmetry inforomation
source
Missing docstring.

Missing docstring for SymmetryInfo. Check Documenter's build log for details.

Xtals.set_path_to_dataFunction
set_path_to_data("../data/")

set the PATH_TO_DATA or PATH_TO_CRYSTALS variable. by default, sets PATH_TO_DATA to path, and PATH_TO_CRYSTALS to path/crystals.

Arguments

  • path::String The path to use for setting the environment variable.
  • relpath_xtals::Bool Specify true to update path to crystals relative to new data path.
  • print::Bool Specify true to print path variables.
source
Xtals.set_path_to_crystalsFunction
set_path_to_crystals("../other_crystals/")

set Xtals.PATH_TO_CRYSTALS.

Arguments

  • path::String The path to use for setting the environment variable.
  • print::Bool Specify true to print path variables.
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
Missing docstring.

Missing docstring for replicate. Check Documenter's build log for details.

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³
source
Xtals.chemical_formulaFunction
formula = chemical_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 a crystal::Crystal to a .cif file with filename::AbstractString. If filename does not include the .cif extension, it will automatically be added. the fractional_coords flag allows us to write either fractional or Cartesian coordinates.

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