Henry Coefficient

PorousMaterials.jl allows for Henry coefficient calculations using Widom insertions.

Preparing the Henry coefficient simulation

The simulation requires the following PorousMaterials.jl objects:

  • Crystal structure
  • Molecule adsorbate
  • LJForceField forcefield

In addition the the list above, one has to specify the temperature (in K) and the number of Widom insertions per unit volume (in Angstrom).

xtal = Crystal("SBMOF-1.cif")           # The crystal structure we are interested in
strip_numbers_from_atom_labels!(xtal)   # We have to make sure the atom species have no numbers appended to them
methane = Molecule("CH4")               # Here we choose to use methane as the adsorbate
ljff = LJForceField("UFF")              # We will use the Universal Force Field (UFF) to calculate the interatomic interactions
temp = 298.0                            # Standard temperature (K)
widom_insertions = 2000                 # Number of insertions per unit volume

results = henry_coefficient(xtal, methane, temp, ljff, insertions_per_volume=widom_insertions)

The results are also saved to rc[:paths][:simulations] as a .jld2 file that can be read using the JLD2 package.

The output (and saved file) is a dictionary:

results
# output
Dict{String, Any} with 16 entries:
  "⟨U⟩ (K)"                              => -3623.51
  "err Qst (kJ/mol)"                     => 0.0917643
  "⟨U, vdw⟩ (kJ/mol)"                    => -30.1275
  "⟨U, es⟩ (kJ/mol)"                     => 0.0
  "elapsed time (min)"                   => 0.0978277
  "Qst (kJ/mol)"                         => 32.6052
  "err henry coefficient [mmol/(g-bar)]" => 6.12256
  "xtal"                                 => "SBMOF-1.cif"
  "henry coefficient [mmol/(g-bar)]"     => 46.1068
  "err ⟨U, es⟩ (kJ/mol)"                 => 0.0
  "⟨U, vdw⟩ (K)"                         => -3623.51
  "err ⟨U, vdw⟩ (kJ/mol)"                => 0.0917643
  "⟨U, es⟩ (K)"                          => 0.0
  "⟨U⟩ (kJ/mol)"                         => -30.1275
  "henry coefficient [mol/(kg-Pa)]"      => 0.000461068
  "henry coefficient [mol/(m³-bar)]"     => 72406.2

locating the saved results

The name of the result filenames follow a convention outlined in henry_result_savename.

using JLD2
# determine the canonical filename for the simulation
result_filename = henry_result_savename(xtal, methane, temp, ljff, widom_insertions)
# load the results dictionary
@load joinpath(rc[:paths][:simulations], result_filename) results

detailed docs

PorousMaterials.henry_coefficientFunction
result = henry_coefficient(crystal, molecule, temperature, ljforcefield,
                            insertions_per_volume=200, verbose=true, ewald_precision=1e-6,
                            autosave=true)

Conduct particle insertions to compute the Henry coefficient Kₕ of a molecule in a crystal. Also, for free, the heat of adsorption and ensemble average energy of adsorption is computed. The Henry coefficient is a model for adsorption at infinite dilution (low coverage): ⟨N⟩ = Kₕ P, where P is pressure and Kₕ is the Henry coefficient.

Kₕ = β ⟨e^{-β U}⟩, where the average is over positions and orientations of the molecule in the crystal.

Arguments

  • crystal::Crystal: the porous crystal in which we seek to simulate adsorption
  • molecule::Molecule: the adsorbate molecule
  • temperature::Float64: temperature of bulk gas phase in equilibrium with adsorbed phase in the porous material. units: Kelvin (K)
  • ljforcefield::LJForceField: the molecular model used to describe the energetics of the adsorbate-adsorbate and adsorbate-host van der Waals interactions.
  • insertions_per_volume::Int: number of Widom insertions to perform for computing the

average, per unit cell volume (ų)

  • verbose::Bool: whether or not to print off information during the simulation.
  • ewald_precision::Float64: desired precision for Ewald summations; used to determine

the replication factors in reciprocal space.

  • autosave::Bool: save results file as a .jld2 in rc[:paths][:simulations]
  • filename_comment::AbstractString: An optional comment that will be appended to the name of the saved file.

Returns

  • result::Dict{String, Float64}: A dictionary containing all the results from the Henry coefficient simulation
source
PorousMaterials.henry_result_savenameFunction
save_name = henry_result_savename(crystal, molecule, temperature,
                               ljforcefield, insertions_per_volume;
                               comment="")

Determine the name of files saved while calculating the henry coefficient. It uses many pieces of information from the simulation to ensure the file name accurately describes what it holds.

Arguments

  • crystal::Crystal: The porous crystal being tested
  • molecule::Molecule: The molecule being tested inside the crystal
  • temperature::Float64: The temperature used in the simulation units: Kelvin (K)
  • ljforcefield::LJForceField: The molecular model being used in the simulation to describe the intermolecular Van der Waals forces
  • insertions_per_volume::Union{Int, Float64}: The number of widom insertions per unit volume. Will be scaled according to the crystal we're working with
  • comment::AbstractString: An optional comment that will be appended to the filename
source