def get_density(molecule_name, temperature=273.15, pressure=101325, cycles=5000, init_cycles="auto", forcefield="CrystalGenerator"): """Calculates the density of a gas through an NPT ensemble. Args: molecule_name: The molecule to test for adsorption. A file of the same name must exist in `$RASPA_DIR/share/raspa/molecules/TraPPE`. temperature: (Optional) The temperature of the simulation, in Kelvin. pressure: (Optional) The pressure of the simulation, in Pascals. cycles: (Optional) The number of simulation cycles to run. init_cycles: (Optional) The number of initialization cycles to run. Defaults to the minimum of cycles / 2 and 10,000. forcefield: (Optional) The forcefield to use. Name must match a folder in `$RASPA_DIR/share/raspa/forcefield`, which contains the properly named `.def` files. Returns: The density, as a float, in kg/m^3. """ print_every = cycles // 10 if init_cycles == "auto": init_cycles = min(cycles // 2, 10000) script = dedent(""" SimulationType {simulation_type} NumberOfCycles {cycles} NumberOfInitializationCycles {init_cycles} PrintEvery {print_every} Forcefield {forcefield} Box 0 BoxLengths 30 30 30 ExternalTemperature {temperature} ExternalPressure {pressure} VolumeChangeProbability 0.25 Component 0 MoleculeName {molecule_name} MoleculeDefinition TraPPE TranslationProbability 0.5 ReinsertionProbability 0.5 CreateNumberOfMolecules 256 """.format(**locals())).strip() output = parse(run_script(script)) return output["Average Density"]["[kg/m^3]"][0]
def get_geometric_surface_area(structure, unit_cells=(1, 1, 1), cycles=500, input_file_type="cif", units="m^2/g", forcefield="CrystalGenerator"): """Calculates the geometric surface area of an inputted structure. Args: structure: The structure to use, as a string of type `input_file_type` (default is "cif"). input_file_type: (Optional) The type of input structure. Assumes cif. unit_cells: (Optional) The number of unit cells to use, by dimension. cycles: (Optional) The number of simulation cycles to run. units: (Optional) The units in which to return the surface area. Can be "m^2/g", "A^2", or "m^2/cm^3". forcefield: (Optional) The forcefield to use. Name must match a folder in `$RASPA_DIR/share/raspa/forcefield`, which contains the properly named `.def` files. Returns: The geometric surface area, as a float. """ print_every = cycles // 10 a, b, c = unit_cells script = dedent(""" SimulationType MonteCarlo NumberOfCycles {cycles} PrintEvery {print_every} PrintPropertiesEvery {print_every} Forcefield {forcefield} CutOff 12.8 Framework 0 FrameworkName streamed InputFileType {input_file_type} UnitCells {a} {b} {c} SurfaceAreaProbeDistance Sigma Component 0 MoleculeName N2 StartingBead 0 MoleculeDefinition TraPPE SurfaceAreaProbability 1.0 CreateNumberOfMolecules 0 """.format(**locals())).strip() output = parse(run_script(script, structure)) return output["Average Surface Area"]["[{}]".format(units)][0]
def run(structure, molecule_name, temperature=273.15, pressure=101325, helium_void_fraction=1.0, unit_cells=(1, 1, 1), framework_name="streamed", simulation_type="MonteCarlo", cycles=2000, init_cycles="auto", forcefield="CrystalGenerator", input_file_type="cif"): """Runs a simulation with the specified parameters. Args: structure: The structure to test for adsorption, as a string of type `input_file_type` (default is "cif"). molecule_name: The molecule to test for adsorption. A file of the same name must exist in `$RASPA_DIR/share/raspa/molecules/TraPPE`. temperature: (Optional) The temperature of the simulation, in Kelvin. pressure: (Optional) The pressure of the simulation, in Pascals. helium_void_fraction: (Optional) The helium void fraction of the input structure. Required for excess adsorption back-calculation. unit_cells: (Optional) The number of unit cells to use, by dimension. framework_name: (Optional) If not streaming, this will load the structure at `$RASPA_DIR/share/raspa/structures`. Ignored if streaming. simulation_type: (Optional) The type of simulation to run, defaults to "MonteCarlo". cycles: (Optional) The number of simulation cycles to run. init_cycles: (Optional) The number of initialization cycles to run. Defaults to the minimum of cycles / 2 and 10,000. forcefield: (Optional) The forcefield to use. Name must match a folder in `$RASPA_DIR/share/raspa/forcefield`, which contains the properly named `.def` files. input_file_type: (Optional) The type of input structure. Assumes cif. Returns: A string representing the contents of a simulation input file. The goal of this function is to mask the complexity of RASPA by limiting parameters and assuming sensible defaults. This should streamline common- case usage, but it means that this function won't work in all use cases. In these cases, look into loading your own simulation input file and passing it to `RASPA.run_script`. """ return parse(run_script(create_script(**locals()), structure))
def get_helium_void_fraction(structure, unit_cells=(1, 1, 1), cycles=2000, input_file_type="cif", forcefield="CrystalGenerator"): """Calculates the helium void fraction of the inputted structure. Args: structure: The structure to test for helium void fraction, as a string of type 'input_file_type` (default is "cif"). unit_cells: (Optional) The number of unit cells to use, by dimension. cycles: (Optional) The number of simulation cycles to run. input_file_type: (Optional) The type of input structure. Assumes cif. forcefield: (Optional) The forcefield to use. Name must match a folder in `$RASPA_DIR/share/raspa/forcefield`, which contains the properly named `.def` files. Returns: The helium void fraction of the structure, as a float. """ print_every = cycles // 10 a, b, c = unit_cells script = dedent(""" SimulationType MonteCarlo NumberOfCycles {cycles} PrintEvery {print_every} PrintPropertiesEvery {print_every} Forcefield {forcefield} CutOff 12.8 Framework 0 FrameworkName streamed InputFileType {input_file_type} UnitCells {a} {b} {c} ExternalTemperature 298.0 Component 0 MoleculeName helium MoleculeDefinition TraPPE WidomProbability 1.0 CreateNumberOfMolecules 0 """.format(**locals())).strip() output = parse(run_script(script, structure)) return output["Average Widom Rosenbluth factor"]["Widom"][0]
def run_mixture(structure, molecules, mol_fractions, temperature=273.15, pressure=101325, helium_void_fraction=1.0, unit_cells=(1, 1, 1), simulation_type="MonteCarlo", cycles=2000, init_cycles="auto", forcefield="CrystalGenerator", input_file_type="cif"): """Runs a simulation with mixture of gases. Args: structure: The structure to test for adsorption, as a string of type `input_file_type` (default is "cif"). molecules: The molecules to test for adsorption. Files of the same names must exist in `$RASPA_DIR/share/raspa/molecules/TraPPE`. mol_fractions: The mol fractions of each gas that you want to separate. Corresponds to the `molecules` list. temperature: (Optional) The temperature of the simulation, in Kelvin. pressure: (Optional) The pressure of the simulation, in Pascals. helium_void_fraction: (Optional) The helium void fraction of the input structure. Required for excess adsorption back-calculation. unit_cells: (Optional) The number of unit cells to use, by dimension. simulation_type: (Optional) The type of simulation to run, defaults to "MonteCarlo". cycles: (Optional) The number of simulation cycles to run. init_cycles: (Optional) The number of initialization cycles to run. Defaults to the minimum of cycles / 2 and 10,000. forcefield: (Optional) The forcefield to use. Name must match a folder in `$RASPA_DIR/share/raspa/forcefield`, which contains the properly named `.def` files. input_file_type: (Optional) The type of input structure. Assumes cif. Returns: A string representing the contents of a simulation input file. The goal of this function is to mask the complexity of RASPA by limiting parameters and assuming sensible defaults. This should streamline common- case usage, but it means that this function won't work in all use cases. In these cases, look into loading your own simulation input file and passing it to `RASPA.run_script`. """ is_mol = "yes" if input_file_type.lower() == "mol" else "no" print_every = cycles // 10 a, b, c = unit_cells molecule_list = " ".join(str(n) for n in range(len(molecules))) molecule_count = len(molecules) if init_cycles == "auto": init_cycles = min(cycles // 2, 10000) script = dedent(""" SimulationType {simulation_type} NumberOfCycles {cycles} NumberOfInitializationCycles {init_cycles} PrintEvery {print_every} RestartFile no Forcefield {forcefield} CutOff 12.8 ChargeMethod Ewald EwaldPrecision 1e-6 UseChargesFromMOLFile {is_mol} Framework 0 FrameworkName streamed InputFileType {input_file_type} UnitCells {a} {b} {c} HeliumVoidFraction {helium_void_fraction} ExternalTemperature {temperature} ExternalPressure {pressure} Movies no WriteMoviesEvery 100 """.format(**locals())).strip() for i, (molecule, fraction) in enumerate(zip(molecules, mol_fractions)): script += dedent(""" Component {i} MoleculeName {molecule} StartingBead 0 MoleculeDefinition TraPPE MolFraction {fraction} IdentityChangeProbability 1.0 NumberOfIdentityChanges {molecule_count} IdentityChangesList {molecule_list} IdealGasRosenbluthWeight 1.0 TranslationProbability 1.0 RotationProbability 1.0 ReinsertionProbability 1.0 SwapProbability 1.0 CreateNumberOfMolecules 0 """.format(**locals())) return parse(run_script(script, structure))