RASPA3 3.0.12
A molecular simulation code for computing adsorption and diffusion in nanoporous materials
Loading...
Searching...
No Matches
Component Struct Referenceexport
module component

Represents a component within the simulation system. More...

Collaboration diagram for Component:

Public Types

enum class  Type : std::size_t { Adsorbate = 0 , Cation = 1 }
 Enumeration of component types. More...
 
enum class  GrowType : std::size_t { Rigid = 0 , Flexible = 1 }
 Enumeration of growth types for the component. More...
 
enum class  Shape : std::size_t { NonLinear = 0 , Linear = 1 , Point = 2 }
 Enumeration of molecular shapes. More...
 
enum class  Chirality : std::size_t { NoCharility = 0 , S_Chiral = 1 , R_Chiral = 2 }
 
enum class  PressureScale { Log = 0 , Normal = 1 }
 Enumeration of pressure scaling types. More...
 

Public Member Functions

 Component ()
 Default constructor for the Component struct.
 
 Component (Component::Type type, std::size_t currentComponent, const ForceField &forceField, const std::string &componentName, std::optional< const std::string > fileName, std::size_t numberOfBlocks, std::size_t numberOfLambdaBins, const MCMoveProbabilities &systemProbabilities=MCMoveProbabilities(), std::optional< double > fugacityCoefficient=std::nullopt, bool thermodynamicIntegration=false) noexcept(false)
 Constructs a Component from a file.
 
 Component (std::size_t componentId, const ForceField &forceField, std::string componentName, double T_c, double P_c, double w, std::vector< Atom > definedAtoms, const ConnectivityTable &connectivityTable, const Potentials::IntraMolecularPotentials &intraMolecularPotentials, std::size_t numberOfBlocks, std::size_t numberOfLambdaBins, const MCMoveProbabilities &systemProbabilities=MCMoveProbabilities(), std::optional< double > fugacityCoefficient=std::nullopt, bool thermodynamicIntegration=false) noexcept(false)
 Constructs a Component programmatically.
 
void readComponent (const ForceField &forceField, const std::string &fileName)
 Reads component data from a file.
 
std::string printStatus (const ForceField &forceField, double inputPressure) const
 Generates a string representing the component's current status.
 
std::string printBreakthroughStatus () const
 Generates a string representing the breakthrough status of the component.
 
nlohmann::json jsonStatus () const
 Serializes the component's status to JSON.
 
void computeRigidProperties ()
 Computes rigid body properties of the component.
 
double3 computeCenterOfMass (std::vector< Atom > atom_list) const
 Computes the center of mass for a given list of atoms.
 
std::vector< AtomrotatePositions (const simd_quatd &q) const
 Rotates the positions of all atoms using a given quaternion.
 
std::vector< AtomcopiedAtoms (std::span< Atom > molecule) const
 Creates a copy of the component's atoms based on a molecule span.
 
std::pair< Molecule, std::vector< Atom > > equilibratedMoleculeRandomInBox (RandomNumber &random, const SimulationBox &simulationBox) const
 Generates an equilibrated molecule randomly placed within a simulation box.
 
std::pair< Molecule, std::vector< Atom > > translate (const Molecule &molecule, std::span< Atom > molecule_atoms, double3 displacement) const
 Translates a molecule by a specified displacement.
 
std::pair< Molecule, std::vector< Atom > > rotate (const Molecule &molecule, std::span< Atom > molecule_atoms, simd_quatd rotation) const
 Rotates a molecule using a specified rotation quaternion.
 
ConnectivityTable readConnectivityTable (std::size_t size, const nlohmann::basic_json< nlohmann::raspa_map > &parsed_data)
 
std::vector< BondPotentialreadBondPotentials (const ForceField &forceField, const nlohmann::basic_json< nlohmann::raspa_map > &parsed_data)
 
std::vector< BendPotentialreadBendPotentials (const ForceField &forceField, const nlohmann::basic_json< nlohmann::raspa_map > &parsed_data)
 
std::vector< TorsionPotentialreadTorsionPotentials (const ForceField &forceField, const nlohmann::basic_json< nlohmann::raspa_map > &parsed_data)
 
std::vector< VanDerWaalsPotentialreadVanDerWaalsPotentials (const ForceField &forceField, const nlohmann::basic_json< nlohmann::raspa_map > &parsed_data)
 
std::vector< std::vector< std::size_t > > readPartialReinsertionFixedAtoms (const nlohmann::basic_json< nlohmann::raspa_map > &parsed_data)
 
std::string repr () const
 Returns a string representation of the Component.
 

Public Attributes

std::uint64_t versionNumber {1}
 Version number for serialization.
 
Type type {0}
 Type of the component (Adsorbate or Cation).
 
GrowType growType {0}
 Growth type of the component.
 
std::size_t componentId {0}
 Unique identifier for the component.
 
std::string name {}
 Name of the component.
 
std::optional< std::string > filenameData {}
 Optional filename containing component data.
 
std::string filename {}
 Filename associated with the component.
 
std::vector< double4 > blockingPockets {}
 List of blocking pockets defined by position and radius.
 
bool rigid {true}
 Flag indicating if the component is rigid.
 
std::size_t translationalDegreesOfFreedom {}
 Number of translational degrees of freedom.
 
std::size_t rotationalDegreesOfFreedom {}
 Number of rotational degrees of freedom.
 
double criticalTemperature {0.0}
 Critical temperature of the component [K].
 
double criticalPressure {0.0}
 Critical pressure of the component [Pa].
 
double acentricFactor {0.0}
 Acentric factor of the component [-].
 
double molFraction {1.0}
 Mole fraction of the component [-].
 
bool swappable {false}
 Flag indicating if the component is swappable.
 
double partialPressure {0.0}
 Partial pressure of the component [Pa].
 
double totalMass {0.0}
 Total mass of the component [kg].
 
std::optional< double > fugacityCoefficient {}
 Optional fugacity coefficient [-].
 
double amountOfExcessMolecules {0.0}
 Amount of excess molecules [-].
 
double bulkFluidDensity {0.0}
 Bulk fluid density [kg/m³].
 
double compressibility {0.0}
 Compressibility of the component [-].
 
std::optional< double > idealGasRosenbluthWeight {}
 Optional ideal gas Rosenbluth weight [-].
 
std::optional< double > idealGasEnergy {}
 Optional ideal gas energy [J].
 
double netCharge {0.0}
 Net charge of the component [e].
 
std::size_t startingBead {0}
 Starting bead index for simulations.
 
std::vector< std::pair< Atom, double > > definedAtoms {}
 List of defined atoms and their masses.
 
double3 inertiaVector {}
 Inertia vector of the component.
 
double3 inverseInertiaVector {}
 Inverse of the inertia vector.
 
Shape shapeType
 Shape type of the molecule.
 
std::vector< Atomatoms {}
 List of atoms in the component.
 
ConnectivityTable connectivityTable {}
 Connectivity table for the component.
 
Potentials::IntraMolecularPotentials intraMolecularPotentials {}
 List of internal potentials.
 
std::vector< AtomgrownAtoms {}
 
std::vector< std::vector< std::size_t > > partialReinsertionFixedAtoms {}
 
std::size_t initialNumberOfMolecules {0}
 Initial number of molecules in the component.
 
PropertyLambdaProbabilityHistogram lambdaGC
 Lambda probability histogram for Gibbs-Chebyshev integration.
 
PropertyLambdaProbabilityHistogram lambdaGibbs
 Lambda probability histogram for Gibbs integration.
 
bool hasFractionalMolecule {false}
 Flag indicating if the component has fractional molecules.
 
MCMoveProbabilities mc_moves_probabilities
 Move probabilities for Monte Carlo simulations.
 
MCMoveStatistics mc_moves_statistics
 
MCMoveCpuTime mc_moves_cputime
 CPU time statistics for Monte Carlo moves.
 
std::vector< CBMCMoveStatisticscbmc_moves_statistics
 
PropertyWidom averageRosenbluthWeights
 Average Rosenbluth weights for Widom insertion.
 
double lnPartitionFunction {0}
 Natural logarithm of the partition function [-].
 
MultiSiteIsotherm isotherm {}
 Isotherm information for the component.
 
double massTransferCoefficient {0.0}
 Mass transfer coefficient [1/s].
 
double axialDispersionCoefficient {0.0}
 Axial dispersion coefficient [m²/s].
 
bool isCarrierGas {false}
 Flag indicating if the component is a carrier gas.
 
std::size_t columnPressure {0}
 Column index for pressure data.
 
std::size_t columnLoading {1}
 Column index for loading data.
 
std::size_t columnError {2}
 Column index for error data.
 
PressureScale pressureScale {PressureScale::Log}
 Pressure scaling type.
 

Friends

Archive< std::ofstream > & operator<< (Archive< std::ofstream > &archive, const Component &c)
 
Archive< std::ifstream > & operator>> (Archive< std::ifstream > &archive, Component &c)
 

Detailed Description

Represents a component within the simulation system.

The Component struct encapsulates all properties and behaviors associated with a simulation component, such as adsorbates or cations. It includes various physical properties, molecular structures, and simulation parameters necessary for conducting simulations. The struct provides constructors for initializing components from files or programmatically and includes methods for computing physical properties and manipulating molecular positions.

Member Enumeration Documentation

◆ GrowType

enum class Component::GrowType : std::size_t
strong

Enumeration of growth types for the component.

Enumerator
Rigid 

Rigid growth type.

Flexible 

Flexible growth type.

◆ PressureScale

enum class Component::PressureScale
strong

Enumeration of pressure scaling types.

Enumerator
Log 

Logarithmic pressure scaling.

Normal 

Normal pressure scaling.

◆ Shape

enum class Component::Shape : std::size_t
strong

Enumeration of molecular shapes.

Enumerator
NonLinear 

Non-linear molecular shape.

Linear 

Linear molecular shape.

Point 

Point-particle molecular shape.

◆ Type

enum class Component::Type : std::size_t
strong

Enumeration of component types.

Enumerator
Adsorbate 

Represents an adsorbate component.

Cation 

Represents a cation component.

Constructor & Destructor Documentation

◆ Component() [1/3]

Component::Component ( )

Default constructor for the Component struct.

Initializes a Component object with default values.

◆ Component() [2/3]

Component::Component ( Component::Type  type,
std::size_t  currentComponent,
const ForceField forceField,
const std::string &  componentName,
std::optional< const std::string >  fileName,
std::size_t  numberOfBlocks,
std::size_t  numberOfLambdaBins,
const MCMoveProbabilities systemProbabilities = MCMoveProbabilities(),
std::optional< double >  fugacityCoefficient = std::nullopt,
bool  thermodynamicIntegration = false 
)

Constructs a Component from a file.

Initializes a Component by reading data from a specified file, using the provided force field and simulation parameters.

Parameters
typeThe type of the component (Adsorbate or Cation).
currentComponentThe current component ID.
forceFieldThe force field used for parsing atom types and interactions.
componentNameThe name of the component.
fileNameOptional filename containing component data.
numberOfBlocksThe number of simulation blocks.
numberOfLambdaBinsThe number of lambda bins for scaling.
systemProbabilitiesMove probabilities for the Monte Carlo simulation.
fugacityCoefficientOptional fugacity coefficient.
thermodynamicIntegrationFlag indicating if thermodynamic integration is used.
Exceptions
std::runtime_errorIf the component file cannot be read or parsed.

◆ Component() [3/3]

Component::Component ( std::size_t  componentId,
const ForceField forceField,
std::string  componentName,
double  T_c,
double  P_c,
double  w,
std::vector< Atom definedAtoms,
const ConnectivityTable connectivityTable,
const Potentials::IntraMolecularPotentials intraMolecularPotentials,
std::size_t  numberOfBlocks,
std::size_t  numberOfLambdaBins,
const MCMoveProbabilities systemProbabilities = MCMoveProbabilities(),
std::optional< double >  fugacityCoefficient = std::nullopt,
bool  thermodynamicIntegration = false 
)

Constructs a Component programmatically.

Initializes a Component with specified physical properties and molecular structure.

Parameters
componentIdThe unique identifier for the component.
forceFieldThe force field used for defining atom properties.
componentNameThe name of the component.
T_cThe critical temperature of the component.
P_cThe critical pressure of the component.
wThe acentric factor of the component.
atomListA list of atoms defining the component's molecular structure.
numberOfBlocksThe number of simulation blocks.
numberOfLambdaBinsThe number of lambda bins for scaling.
systemProbabilitiesMove probabilities for the Monte Carlo simulation.
fugacityCoefficientOptional fugacity coefficient.
thermodynamicIntegrationFlag indicating if thermodynamic integration is used.
Exceptions
std::runtime_errorIf pseudo-atoms are not recognized or data is invalid.

Member Function Documentation

◆ computeCenterOfMass()

double3 Component::computeCenterOfMass ( std::vector< Atom atom_list) const

Computes the center of mass for a given list of atoms.

Calculates the center of mass position based on the provided atom list.

Parameters
atom_listThe list of atoms to compute the center of mass for.
Returns
The center of mass position as a double3 vector.

◆ computeRigidProperties()

void Component::computeRigidProperties ( )

Computes rigid body properties of the component.

Calculates properties such as the center of mass, inertia tensor, and degrees of freedom based on the defined atoms.

◆ copiedAtoms()

std::vector< Atom > Component::copiedAtoms ( std::span< Atom molecule) const

Creates a copy of the component's atoms based on a molecule span.

Copies and adjusts the positions of atoms relative to a given molecule.

Parameters
moleculeA span representing the molecule to copy atoms from.
Returns
A vector of copied atoms with adjusted positions.

◆ equilibratedMoleculeRandomInBox()

std::pair< Molecule, std::vector< Atom > > Component::equilibratedMoleculeRandomInBox ( RandomNumber &  random,
const SimulationBox simulationBox 
) const

Generates an equilibrated molecule randomly placed within a simulation box.

Rotates and translates the molecule randomly to ensure equilibration within the simulation box.

Parameters
randomA random number generator instance.
simulationBoxThe simulation box within which to place the molecule.
Returns
A pair containing the equilibrated molecule and its corresponding atoms.

◆ jsonStatus()

nlohmann::json Component::jsonStatus ( ) const

Serializes the component's status to JSON.

Converts the component's current status and properties into a JSON object.

Returns
A JSON object representing the component's status.

◆ printBreakthroughStatus()

std::string Component::printBreakthroughStatus ( ) const

Generates a string representing the breakthrough status of the component.

Compiles breakthrough-related properties and parameters of the component into a formatted string.

Returns
A string detailing the breakthrough status of the component.

◆ printStatus()

std::string Component::printStatus ( const ForceField forceField,
double  inputPressure 
) const

Generates a string representing the component's current status.

Compiles various properties and parameters of the component into a formatted string.

Parameters
forceFieldThe force field used for interpreting atom types.
Returns
A string detailing the component's status.

◆ readComponent()

void Component::readComponent ( const ForceField forceField,
const std::string &  fileName 
)

Reads component data from a file.

Parses component information from the specified file using the provided force field.

Parameters
forceFieldThe force field used for parsing atom types and interactions.
fileNameThe name of the file containing component data.
Exceptions
std::runtime_errorIf the file cannot be found, opened, or parsed correctly.

◆ repr()

std::string Component::repr ( ) const

Returns a string representation of the Component.

Generates a simple string indicating a test representation of the Component.

Returns
A string representing the Component.

◆ rotate()

std::pair< Molecule, std::vector< Atom > > Component::rotate ( const Molecule molecule,
std::span< Atom molecule_atoms,
simd_quatd  rotation 
) const

Rotates a molecule using a specified rotation quaternion.

Applies a rotation to the entire molecule based on the provided quaternion.

Parameters
moleculeThe molecule to rotate.
molecule_atomsA span representing the atoms of the molecule.
rotationThe quaternion representing the rotation to apply.
Returns
A pair containing the rotated molecule and its corresponding atoms.

◆ rotatePositions()

std::vector< Atom > Component::rotatePositions ( const simd_quatd &  q) const

Rotates the positions of all atoms using a given quaternion.

Applies a rotation to the atom positions based on the provided quaternion.

Parameters
qThe quaternion representing the rotation.
Returns
A vector of atoms with updated positions after rotation.

◆ translate()

std::pair< Molecule, std::vector< Atom > > Component::translate ( const Molecule molecule,
std::span< Atom molecule_atoms,
double3  displacement 
) const

Translates a molecule by a specified displacement.

Moves the entire molecule by the given displacement vector.

Parameters
moleculeThe molecule to translate.
molecule_atomsA span representing the atoms of the molecule.
displacementThe displacement vector to apply.
Returns
A pair containing the translated molecule and its corresponding atoms.

The documentation for this struct was generated from the following file: