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

Represents the central system for simulations. More...

Collaboration diagram for System:

Public Member Functions

 System ()=default
 Default constructor for the System class.
 
 System (std::size_t id, ForceField forcefield, std::optional< SimulationBox > box, double T, std::optional< double > P, double heliumVoidFraction, std::optional< Framework > framework, std::vector< Component > components, std::vector< std::vector< double3 > > initialPositions, std::vector< std::size_t > initialNumberOfMolecules, std::size_t numberOfBlocks, const MCMoveProbabilities &systemProbabilities=MCMoveProbabilities(), std::optional< std::size_t > sampleMoviesEvery=std::nullopt)
 Constructs a System object programmatically.
 
 System (std::size_t id, double T, std::optional< double > P, double heliumVoidFraction, std::optional< Framework > framework, std::vector< Component > components)
 
std::optional< double > frameworkMass () const
 
std::size_t indexOfGCFractionalMoleculesPerComponent_CFCMC (std::size_t selectedComponent)
 The fractional molecule for grand-canonical is stored first.
 
std::size_t indexOfPairGCFractionalMoleculesPerComponent_CFCMC (std::size_t selectedComponent)
 The fractional molecule for grand-canonical pair-insertion is stored second.
 
std::size_t indexOfGibbsFractionalMoleculesPerComponent_CFCMC (std::size_t selectedComponent)
 The fractional molecule for Gibbs is stored third.
 
void addComponent (const Component &&component) noexcept(false)
 
void createFrameworks ()
 
void createInitialMolecules (const std::vector< std::vector< double3 > > &initialPositions)
 
void checkCartesianPositions ()
 
void precomputeTotalRigidEnergy () noexcept
 
void precomputeTotalGradients () noexcept
 
RunningEnergy computeTotalEnergies () noexcept
 
RunningEnergy computePolarizationEnergy () noexcept
 
RunningEnergy computeTotalGradients () noexcept
 
void computeTotalElectrostaticPotential () noexcept
 
void computeTotalElectricField () noexcept
 
std::size_t randomFramework (RandomNumber &random)
 
std::size_t randomComponent (RandomNumber &random)
 
std::size_t numerOfAdsorbateComponents ()
 
std::size_t randomMoleculeOfComponent (RandomNumber &random, std::size_t selectedComponent)
 
std::size_t randomIntegerMoleculeOfComponent (RandomNumber &random, std::size_t selectedComponent)
 
std::size_t indexOfFirstMolecule (std::size_t selectedComponent)
 
std::vector< Atom >::iterator iteratorForMolecule (std::size_t selectedComponent, std::size_t selectedMolecule)
 
std::vector< double3 >::iterator iteratorForElectricField (std::size_t selectedComponent, std::size_t selectedMolecule)
 
std::vector< Molecule >::iterator indexForMolecule (std::size_t selectedComponent, std::size_t selectedMolecule)
 
std::size_t moleculeIndexOfComponent (std::size_t selectedComponent, std::size_t selectedMolecule)
 
std::span< AtomspanOfMolecule (std::size_t selectedComponent, std::size_t selectedMolecule)
 
const std::span< const AtomspanOfMolecule (std::size_t selectedComponent, std::size_t selectedMolecule) const
 
const std::span< const AtomspanOfIntegerAtomsOfComponent (std::size_t selectedComponent) const
 
std::span< const AtomspanOfFrameworkAtoms () const
 
std::span< AtomspanOfFrameworkAtoms ()
 
std::span< const AtomspanOfRigidFrameworkAtoms () const
 
std::span< const AtomspanOfFlexibleAtoms () const
 
std::span< const AtomspanOfMoleculeAtoms () const
 
std::span< AtomspanOfMoleculeAtoms ()
 
std::span< double > spanOfMoleculeElectrostaticPotential ()
 
std::span< double3 > spanOfMoleculeElectricField ()
 
std::span< double3 > spanOfMoleculeElectricFieldNew ()
 
std::span< double3 > spanElectricFieldNew (std::size_t selectedComponent, std::size_t selectedMolecule)
 
const std::span< const double3 > spanElectricFieldNew (std::size_t selectedComponent, std::size_t selectedMolecule) const
 
std::span< double3 > spanElectricFieldOld (std::size_t selectedComponent, std::size_t selectedMolecule)
 
const std::span< const double3 > spanElectricFieldOld (std::size_t selectedComponent, std::size_t selectedMolecule) const
 
std::size_t numberOfMolecules () const
 
std::size_t numberOfIntegerMolecules () const
 
double weight () const
 
void removeRedundantMoves ()
 
void rescaleMoveProbabilities ()
 
void determineSwappableComponents ()
 
void determineFractionalComponents ()
 
void rescaleMolarFractions ()
 
void computeComponentFluidProperties ()
 
void computeNumberOfPseudoAtoms ()
 
void optimizeMCMoves ()
 
std::string writeOutputHeader () const
 
std::string writeNumberOfPseudoAtoms () const
 
std::string writeInitializationStatusReport (std::size_t currentCycle, std::size_t numberOfCycles) const
 
std::string writeEquilibrationStatusReportMC (std::size_t currentCycle, std::size_t numberOfCycles) const
 
std::string writeEquilibrationStatusReportMD (std::size_t currentCycle, std::size_t numberOfCycles) const
 
std::string writeProductionStatusReportMC (const std::string &statusLine) const
 
std::string writeProductionStatusReportMD (std::size_t currentCycle, std::size_t numberOfCycles) const
 
std::string writeSystemStatus () const
 
std::string writeComponentStatus () const
 
std::string writeMCMoveStatistics () const
 
nlohmann::json jsonSystemStatus () const
 
nlohmann::json jsonComponentStatus () const
 
nlohmann::json jsonMCMoveStatistics () const
 
void insertMolecule (std::size_t selectedComponent, const Molecule &molecule, std::vector< Atom > atoms)
 
void insertMoleculePolarization (std::size_t selectedComponent, const Molecule &molecule, std::vector< Atom > atoms, std::span< double3 > electricField)
 
void insertFractionalMolecule (std::size_t selectedComponent, const Molecule &molecule, std::vector< Atom > atoms, std::size_t moleculeId)
 
void deleteMolecule (std::size_t selectedComponent, std::size_t selectedMolecule, const std::span< Atom > atoms)
 
void updateMoleculeAtomInformation ()
 
void checkMoleculeIds ()
 
std::vector< AtomrandomConfiguration (RandomNumber &random, std::size_t selectedComponent, const std::span< const Atom > atoms)
 
bool insideBlockedPockets (const Component &component, std::span< const Atom > molecule_atoms) const
 
void sampleProperties (std::size_t currentBlock, std::size_t currentCycle)
 
void writeCPUTimeStatistics (std::ostream &stream) const
 
std::pair< EnergyStatus, double3x3 > computeMolecularPressure () noexcept
 
std::pair< std::vector< Molecule >, std::vector< Atom > > scaledCenterOfMassPositions (double scale) const
 
void writeComponentFittingStatus (std::ostream &stream, const std::vector< std::pair< double, double > > &rawData) const
 
void createInterpolationGrids (std::ostream &stream)
 
void writeRestartFile ()
 
std::string repr () const
 

Public Attributes

std::uint64_t versionNumber {1}
 
std::size_t systemId {}
 
double temperature {300.0}
 
double pressure {1e4}
 
double input_pressure {1e4}
 
double beta {1.0 / (Units::KB * 300.0)}
 
double heliumVoidFraction {0.29}
 
std::size_t numberOfFrameworks {0}
 
std::size_t numberOfFrameworkAtoms {0}
 
std::size_t numberOfRigidFrameworkAtoms {0}
 
std::optional< Frameworkframework
 
std::vector< Componentcomponents
 
EquationOfState equationOfState
 
std::optional< Thermostatthermostat
 
Loadings loadings
 
std::vector< std::size_t > swappableComponents {}
 
std::vector< std::size_t > initialNumberOfMolecules {}
 
std::vector< std::size_t > numberOfMoleculesPerComponent {}
 
std::vector< std::size_t > numberOfIntegerMoleculesPerComponent {}
 
std::vector< std::size_t > numberOfFractionalMoleculesPerComponent {}
 
std::vector< std::size_t > numberOfGCFractionalMoleculesPerComponent_CFCMC {}
 
std::vector< std::size_t > numberOfPairGCFractionalMoleculesPerComponent_CFCMC {}
 
std::vector< std::size_t > numberOfGibbsFractionalMoleculesPerComponent_CFCMC {}
 
std::vector< std::vector< std::size_t > > numberOfReactionFractionalMoleculesPerComponent_CFCMC {}
 
std::vector< double > idealGasEnergiesPerComponent {}
 
ForceField forceField
 
bool hasExternalField
 
std::vector< std::vector< std::size_t > > numberOfPseudoAtoms
 
std::vector< std::size_t > totalNumberOfPseudoAtoms
 
std::vector< std::pair< double, std::size_t > > numberOfIntegerPseudoAtoms
 
std::vector< std::pair< double, std::size_t > > numberOfFractionalPseudoAtoms
 
std::size_t translationalCenterOfMassConstraint {}
 
std::size_t translationalDegreesOfFreedom {}
 
std::size_t rotationalDegreesOfFreedom {}
 
double timeStep {0.0005}
 
SimulationBox simulationBox
 
bool containsTheFractionalMolecule {true}
 
std::vector< AtomatomData
 
std::vector< MoleculemoleculeData
 
std::vector< double > electricPotential
 
std::vector< double3 > electricField
 
std::vector< double3 > electricFieldNew
 
double conservedEnergy {}
 
double referenceEnergy {}
 
double accumulatedDrift {}
 
RunningEnergy rigidEnergies
 
RunningEnergy runningEnergies
 
double3x3 currentExcessPressureTensor
 
EnergyStatus currentEnergyStatus
 
std::size_t numberOfHybridMCSteps {10}
 
std::vector< std::complex< double > > eik_xy {}
 
std::vector< std::complex< double > > eik_x {}
 
std::vector< std::complex< double > > eik_y {}
 
std::vector< std::complex< double > > eik_z {}
 
std::vector< std::pair< std::complex< double >, std::complex< double > > > storedEik {}
 
std::vector< std::pair< std::complex< double >, std::complex< double > > > fixedFrameworkStoredEik {}
 
std::vector< std::pair< std::complex< double >, std::complex< double > > > totalEik {}
 
double CoulombicFourierEnergySingleIon {0.0}
 
double netCharge {0.0}
 
double netChargeFramework {0.0}
 
double netChargeAdsorbates {0.0}
 
std::vector< double > netChargePerComponent
 
MCMoveProbabilities mc_moves_probabilities
 
MCMoveStatistics mc_moves_statistics
 
MCMoveCpuTime mc_moves_cputime
 
Reactions reactions
 
TransitionMatrix tmmc
 
PropertyEnergy averageEnergies
 
PropertyLoading averageLoadings
 
PropertyEnthalpy averageEnthalpiesOfAdsorption
 
PropertyTemperature averageTemperature
 
PropertyTemperature averageTranslationalTemperature
 
PropertyTemperature averageRotationalTemperature
 
PropertyPressure averagePressure
 
PropertySimulationBox averageSimulationBox
 
std::optional< SampleMoviesamplePDBMovie
 
std::optional< PropertyConventionalRadialDistributionFunctionpropertyConventionalRadialDistributionFunction
 
std::optional< PropertyRadialDistributionFunctionpropertyRadialDistributionFunction
 
std::optional< PropertyDensityGridpropertyDensityGrid
 
std::optional< PropertyEnergyHistogramaverageEnergyHistogram
 
std::optional< PropertyNumberOfMoleculesHistogramaverageNumberOfMoleculesHistogram
 
std::optional< PropertyMeanSquaredDisplacementpropertyMSD
 
std::optional< PropertyVelocityAutoCorrelationFunctionpropertyVACF
 
std::optional< WriteLammpsDatawriteLammpsData
 
std::size_t columnNumberOfGridPoints {100}
 
double columnTotalPressure {1e5}
 
double columnPressureGradient {0.0}
 
double columnVoidFraction {0.4}
 
double columnParticleDensity {1000}
 
double columnEntranceVelocity {0.1}
 
double columnLength {0.3}
 
double columnTimeStep {0.0005}
 
std::size_t columnNumberOfTimeSteps {0}
 
bool columnAutoNumberOfTimeSteps {true}
 
MultiSiteIsotherm::PredictionMethod mixturePredictionMethod {MultiSiteIsotherm::PredictionMethod::IAST}
 
PressureRange pressure_range
 
std::size_t numberOfCarrierGases {0}
 
std::size_t carrierGasComponent {0}
 
std::size_t maxIsothermTerms {0}
 
std::vector< std::optional< InterpolationEnergyGrid > > interpolationGrids
 

Friends

Archive< std::ofstream > & operator<< (Archive< std::ofstream > &archive, const System &s)
 
Archive< std::ifstream > & operator>> (Archive< std::ifstream > &archive, System &s)
 

Detailed Description

Represents the central system for simulations.

The System struct holds all objects required for simulations, such as atom lists, frameworks, and simulation boxes. It is passed to various simulation methods, like MonteCarlo, to conduct simulations. The system can contain frameworks (zeolite/MOF) or a simulation box, and it manages various simulation parameters and states.

Constructor & Destructor Documentation

◆ System() [1/2]

System::System ( )
default

Default constructor for the System class.

Initializes a System object with default values.

◆ System() [2/2]

System::System ( std::size_t  id,
ForceField  forcefield,
std::optional< SimulationBox box,
double  T,
std::optional< double >  P,
double  heliumVoidFraction,
std::optional< Framework framework,
std::vector< Component components,
std::vector< std::vector< double3 > >  initialPositions,
std::vector< std::size_t >  initialNumberOfMolecules,
std::size_t  numberOfBlocks,
const MCMoveProbabilities systemProbabilities = MCMoveProbabilities(),
std::optional< std::size_t >  sampleMoviesEvery = std::nullopt 
)

Constructs a System object programmatically.

Initializes a System with the provided parameters, including system ID, simulation box, temperature, pressure, force field, framework components, and other simulation-related configurations.

Parameters
idThe identifier for the system.
boxThe optional simulation box for the system.
TThe temperature of the system.
PThe optional pressure of the system.
forcefieldThe force field used in the simulation.
frameworkComponentThe list of framework components in the system.
componentsThe list of components in the system.
initialNumberOfMoleculesThe initial number of molecules per component.
numberOfBlocksThe number of blocks in the simulation.
systemProbabilitiesThe move probabilities for the Monte Carlo simulation.
sampleMoviesEveryInterval in which movies are written to PDB.

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