Molecules
Loading Molecule Files
Molecule input files are stored in PorousMaterials.PATH_TO_MOLECULES
. Each molecule possesses its own directory containing two files: charges.csv
and atoms.csv
, comma-separated-value files, which describe the point charges and Lennard Jones spheres, respectively, that compose the molecule. Only rigid molecules are currently supported. Units of length are in Angstroms ($\AA$); units of charges are electrons.
molecule = Molecule("CO2")
# output
Molecule species: CO2
Center of mass (fractional coords): Cart([0.0; 0.0; 0.0;;])
Atoms:
atom = C_CO2, x = [0.000, 0.000, 0.000]
atom = O_CO2, x = [-1.160, 0.000, 0.000]
atom = O_CO2, x = [1.160, 0.000, 0.000]
Point charges:
charge = 0.700000, x = [0.000, 0.000, 0.000]
charge = -0.350000, x = [-1.160, 0.000, 0.000]
charge = -0.350000, x = [1.160, 0.000, 0.000]
Building Blocks of PorousMaterials: Molecules
# access the attributes that comprise the molecule object
molecule.species # molecule species
molecule.com # center-of-mass
molecule.atoms # Lennard-Jones spheres
molecule.charges # point charges
To see specific information about the atoms and charges attributes of the molecule see Atoms
and Charges
.
Moving Molecules
We can translate and roatate a molecule:
# convert to Molecule{Frac}
molecule = Frac(molecule, unit_cube())
# translate center-of-mass to [1.0, 2.0, 3.0] in fractional coordinates
x = Cart([1.0, 2.0, 3.0])
translate_to!(molecule, x, unit_cube())
# translate by [0.1, 0.0, 0.0] in fractional coordinates
dx = Cart([0.1, 0.0, 0.0])
translate_by!(molecule, dx, unit_cube())
# conduct a uniform random rotation about the center-of-mass
random_rotation!(molecule, unit_cube())
detailed docs
Molecules
PorousMaterials.Molecule
— TypeData structure for a molecule/adsorbate.
Attributes
species::Symbol
: Species of molecule, e.g.:CO2
atoms::Atoms
: array of Lennard-Jones spheres comprising the moleculecharges::Charges
: array of point charges comprising the moleculecom::Coords
: center of mass
PorousMaterials.translate_to!
— Functiontranslate_to!(molecule, xf)
translate_to!(molecule, x)
translate_to!(molecule, xf, box)
translate_to!(molecule, x, box)
Translate a molecule so that its center of masss is at a point xf
in fractional coordinate space or at x
in Cartesian coordinate space. For the latter, a unit cell box is required for context.
PorousMaterials.random_rotation!
— Functionrandom_rotation!(molecule{Frac}, box)
random_rotation!(molecule{Cart})
randomly rotate a molecule about its center of mass.
PorousMaterials.random_rotation_matrix
— Methodr = random_rotation_matrix() # rotation matrix in cartesian coords
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)
. the point x
is in Cartesian coordinates here. 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.ion
— Functionmolecule = ion(q, coords)
Facilitate constructing a point charge by constructing a molecule: Molecule(:ion, Atoms{Frac}(0), Charges(q, coords), coords)
Arguments
q::Float64
: value of point charge, units: electronscoords::Frac
: fractional coordinates of the charge
Returns
molecule::Molecule{Frac}
: the ion as a molecule with Fractional coordinates
PorousMaterials.distortion
— Functionis_distorted = distortion(molecule, ref_molecule, box;
atol=1e-12, throw_warning=true)
Determine whether a molecule has distortion w.r.t. a reference molecule via pairwise distance comparison of the atoms and charges coordinates.
Arguments
molecule::Molecule{Frac}
: molecule you want to compareref_molecule::Molecule{Frac}
: reference moleculebox::Box
: box used for the fractional coordinatesatol::Float64=1e-12
: absolute tolerance for distance comparisonthrow_warning::Bool=true
: issue a warning if there is distortion
Returns
is_distorted::Bool
: true if there is distortion w.r.t. reference molecule
Molecular Movement
PorousMaterials.apply_periodic_boundary_condition!
— Functionapply_periodic_boundary_condition!(molecule::Molecule{Frac})
Check if each of a molecule's center-of-mass coordinates is within the bounds (0.0, 1.0) in fractional coordinates, and translate if needed.
Arguments
molecule::Molecule{Frac}
: the molecule to be checked
Returns
nothing
, if the molecule is within the boundary; ohterwise, the coordinates of the input molecule will be modified
PorousMaterials.random_insertion!
— Functionrandom_insertion!(molecules::Array{Molecule{Frac}, 1}, box::Box, template::Molecule{Cart})
Insert a molecule into the simulation box and perform a random rotation if needed. this function calls (insert_w_random_orientation!
)[@ref], supplying it with a random position vector.
Arguments
molecules::Array{Molecule{Frac}, 1}
: array containing the molecules in the simulationbox::Box
: the box used for fractional coordinatstemplate::Molecule{Cart}
: reference molecule of the type inserted
PorousMaterials.remove_molecule!
— Functionremove_molecule!(molecule_id, molecules)
Remove a molecule from the array of molecules.
Arguments
molecule_id::Int
: the ID of the molecule to be removedmolecules::Array{<:Molecule, 1}
: array of molecules to be modified
PorousMaterials.random_translation!
— Functionrandom_reinsertion!(molecule, box)
Perform a translational perturbation in Cartesian coordinates on a molecule, apply the periodic boundry conditions, and keep a copy of the original in case it needs to be restored.
Arguements
molecule::Molecule{Frac}
: molecule to be translatedbox::Box
: the box used in rotation
Returns
old_molecule::Molecule{Frac}
: a copy of the original molecule
PorousMaterials.random_reinsertion!
— Functionrandom_reinsertion!(molecule, box)
Perform a translation and rotated (if needed) on a molecule, and keep a copy of the original in case it needs to be restored.
Arguements
molecule::Molecule{Frac}
: molecule to be translated and rotated (if needed)box::Box
: the box used in rotation
Returns
old_molecule::Molecule{Frac}
: a copy of the original molecule
PorousMaterials.needs_rotations
— Functionneeds_rotations(molecule) # true or false
Determine whether a molecule needs to undergo rotation.
`true` if molecule.atoms.n + molecule.charges.n > 1