Loading in Molecule Files#
Molecule input files are stored in PorousMaterials.PATH_TO_MOLECULES. Each molecule possesses its own directory and contains two files: point_charges.csv and lennard_jones_spheres.csv, comma-separated-value files describing the point charges and Lennard Jones spheres, respectively, comprising the molecule. Only rigid molecules are currently supported. Units of length are in Angstrom; units of charges are electrons.
using PorousMaterials
m = Molecule("CO2")
PorousMaterials will then output information about the molecule you just loaded:
Molecule species: CO2
Center of mass (fractional coords): [0.0, 0.0, 0.0]
Atoms:
atom = C_CO2, xf = [0.000, 0.000, 0.000]
atom = O_CO2, xf = [-1.160, 0.000, 0.000]
atom = O_CO2, xf = [1.160, 0.000, 0.000]
Point charges:
charge = 0.700000, xf = [0.000, 0.000, 0.000]
charge = -0.350000, xf = [-1.160, 0.000, 0.000]
charge = -0.350000, xf = [1.160, 0.000, 0.000]
Building Blocks of PorousMaterials: Molecules#
molecule = Molecule("CO2") # fractional coords in terms of unit cube box
# access Lennard-Jones spheres & point charges that comprise molecule
molecule.atoms
molecule.charges
# translate to [1.0, 2.0, 3.0] fractional coordinates
translate_to!(molecule, [1.0, 2.0, 3.0])
# translate by [0.1, 0.0, 0.0] fractional coordinates
translate_by!(molecule, [0.1, 0.0, 0.0])
# conduct a uniform random rotation
rotate!(molecule, UnitCube()) # b/c now fractional coords defined in context of a unit cube
Molecules#
#
PorousMaterials.Molecule — Type.
Data structure for a molecule/adsorbate.
Attributes
species::Symbol: Species of molecule, e.g.:CO2atoms::Atoms: array of Lennard-Jones spheres comprising the moleculecharges::Charges: array of point charges comprising the moleculexf_com::Array{Float64, 1}: center of mass of the molecule in fractional coordinates
#
PorousMaterials.n_atoms — Function.
num_atoms = n_atoms(molecules)
calculates the total number of atoms in an array of molecules
Arguments
molecule::Array{Molecule, 1}: The molecules to count the number of atoms in
Returns
- The total number of atoms in the molecules passed in
#
PorousMaterials.translate_to! — Function.
translate_to!(molecule, xf)
translate_to!(molecule, x, box)
Translate a molecule a molecule to point xf in fractional coordinate space or to x in Cartesian coordinate space. For the latter, a unit cell box is required for context. The molecule is translated such that its center of mass is at xf/x`.
Arguments
molecule::Molecule: The molecule which will be translated toxfxf::Array{Float64, 1}: A vector containing the coordinates of the final destination of the molecule
#
PorousMaterials.rotate! — Function.
rotate!(molecule, box)
Conduct a random rotation of the molecule about its center of mass. The box is needed because the molecule contains only its fractional coordinates.
Arguments
molecule::Molecule: The molecule which will be subject to a random rotationbox::Box: The molecule only contains fractional coordinates, so the box is needed for a correct rotation
#
PorousMaterials.rotation_matrix — Method.
r = rotation_matrix()
Generate a 3x3 random rotation matrix r such that when a point x is rotated using this rotation matrix via r * x, this point x is placed at a uniform random distributed position on the surface of a sphere of radius norm(x). See James Arvo. Fast Random Rotation Matrices.
https://pdfs.semanticscholar.org/04f3/beeee1ce89b9adf17a6fabde1221a328dbad.pdf
Returns
r::Array{Float64, 2}: A 3x3 random rotation matrix
#
PorousMaterials.rotation_matrix — Method.
R = rotation_matrix(θ, u, assume_unit_vector=false) # 3 by 3 rotation matrix, angle θ about vector u
R = rotation_matrix(θ, dim) # 3 by 3 rotation matrix, angle θ about axis `dim`
Determine the 3D rotation matrix to rotate an angle θ (radians) about axis u.
See Wikipedia.
Arguments
θ::Float64: angle to rotate about an axis, in radiansu::Array{Float64, 1}: axis about which to rotatedim::Int: 1, 2, 3 for rotation about x-, y-, or z-axis, respectively.assume_unit_vector::Bool: assumeuis a unit vector; otherwise,uwill be normalized
internal to this function.
Returns
R::Array{Float64, 2}: 3D rotation matrix. soR * xwill rotate vectorxas desired.
#
PorousMaterials.rand_point_on_unit_sphere — Function.
u = rand_point_on_unit_sphere()
Generate a unit vector with a random orientation.
Returns
u::Array{Float64, 1}: A unit vector with a random orientation
#
PorousMaterials.charged — Method.
charged_flag = charged(molecule, verbose=false)
Determine if a molecule has point charges
Arguments
molecule::Molecule: The molecule which will be checked for chargesverbose::Bool: Will print result iftrue
Returns
charged_flag::Bool:trueif molecule is charged,falseotherwise
Molecular Movement#
#
PorousMaterials.insert_molecule! — Function.
insert_molecule!(molecules, box, template)
Inserts an additional adsorbate molecule into the simulation box using the template provided. The center of mass of the molecule is chosen at a uniform random position in the simulation box. A uniformly random orientation of the molecule is chosen by rotating about the center of mass.
Arguments
molecules::Array{Molecule, 1}: An array of Molecule objectsbox::Box: The simulation boxtemplate::Molecule: A template molecule used as reference when inserting molecules
#
PorousMaterials.delete_molecule! — Function.
delete_molecule!(molecule_id, molecules)
Removes a random molecule from the current molecules in the framework. molecule_id decides which molecule will be deleted, for a simulation, it must be a randomly generated value
Arguments
molecule_id::Int: The molecule ID is used to determine which molecule inmoleculesshould be removedmolecules::Array{Molecule, 1}: An array of Molecule objects
#
PorousMaterials.translate_molecule! — Function.
translate_molecule!(molecule, box)
Perturbs the Cartesian coordinates of a molecule about its center of mass by a random vector of max length δ. Applies periodic boundary conditions to keep the molecule inside the simulation box. Returns a deep copy of the old molecule in case it needs replaced if the Monte Carlo proposal is rejected.
Arguments
molecule::Molecule: The molecule we want to perturbbox::Box: The simulation box
Returns
old_molecule::Molecule: The old molecule in case the MC proposal is rejected
#
PorousMaterials.reinsert_molecule! — Function.
reinsert_molecule(molecule, box)
Move molecule to a new center of mass randomly distrubted in the unit cell and choose a random orientation for it. Return a deep copy of the starting molecule for possible restoration. This MC move can be viewed as a more aggressive translate_molecule!.
Arguments
molecule::Molecule: The molecule we want to perturbbox::Box: The simulation box
#
PorousMaterials.rotatable — Function.
need_to_rotate = rotatable(molecule)
Determines whether or not a given molecule needs to be rotated. For example, rotating a single atom isn't necessary.
Arguments
molecule::Molecule: The molecule being tested. This function determines if a rotation of this molecule will do anything.
Returns
is_rotatable::Bool: A boolean describing whether or not rotating the molecule will alter its interactions with other molecules