ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
System::System Class Reference

Main system class. More...

#include <System.hpp>

+ Inheritance diagram for System::System:
+ Collaboration diagram for System::System:

Public Member Functions

 System (Private)
 
auto get_time_step () const
 Get time_step.
 
void set_time_step (double value)
 Set time_step.
 
auto get_sim_time () const
 Get sim_time.
 
void set_sim_time (double value)
 Set sim_time.
 
auto get_force_cap () const
 Get force_cap.
 
void set_force_cap (double value)
 Set force_cap.
 
auto get_min_global_cut () const
 Get min_global_cut.
 
void set_min_global_cut (double value)
 Set min_global_cut.
 
void set_box_l (Utils::Vector3d const &box_l)
 Change the box dimensions.
 
void tune_verlet_skin (double min_skin, double max_skin, double tol, int int_steps, bool adjust_max_skin)
 Tune the Verlet skin.
 
void set_cell_structure_topology (CellStructureType topology)
 Change cell structure topology.
 
void rebuild_cell_structure ()
 Rebuild cell lists.
 
void rebuild_aosoa ()
 
double maximal_cutoff () const
 Calculate the maximal cutoff of all interactions.
 
double get_interaction_range () const
 Get the interaction range.
 
unsigned get_global_ghost_flags () const
 Returns the ghost flags required for running pair kernels for the global state, e.g.
 
bool long_range_interactions_sanity_checks () const
 Check electrostatic and magnetostatic methods are properly initialized.
 
std::shared_ptr< Observable_statcalculate_energy ()
 Calculate the total energy.
 
std::shared_ptr< Observable_statcalculate_pressure ()
 Calculate the pressure from a virial expansion.
 
double get_instantaneous_pressure ()
 get the instantaneous pressure with (q(t+dt), p(t+dt/2))
 
double get_instantaneous_pressure_virial ()
 get the instantaneous virial pressure with q(t+dt)
 
void synchronize_npt_state ()
 Synchronize NpT state such as instantaneous and average pressure.
 
void npt_ensemble_init (bool recalc_forces)
 Reinitialize the NpT state.
 
void npt_add_virial_contribution (double energy)
 
bool has_npt_enabled () const
 
Utils::Vector3dget_npt_virial () const
 
void calculate_forces ()
 Calculate all forces.
 
void calculate_long_range_fields ()
 Calculate dipole fields.
 
bool has_collision_detection_enabled () const
 
double particle_short_range_energy_contribution (int pid)
 Compute the short-range energy of a particle.
 
std::optional< doubleparticle_bond_energy (int pid, int bond_id, std::vector< int > partners)
 Compute the energy of a given bond which has to exist on the given particle.
 
int integrate (int n_steps, int reuse_forces)
 Integrate equations of motion.
 
int integrate_with_signal_handler (int n_steps, int reuse_forces, bool update_accumulators)
 
void lb_couple_particles ()
 Calculate particle-lattice interactions.
 
void update_dependent_particles ()
 Update particles with properties depending on other particles, namely virtual sites and ICC charges.
 
void update_used_propagations ()
 Update the global propagation bitmask.
 
void check_kT (double value) const
 Veto temperature change.
 
Hook procedures

These procedures are called if several significant changes to the system happen which may make a reinitialization of subsystems necessary.

void on_boxl_change (bool skip_method_adaption=false)
 Called when the box length has changed.
 
void on_node_grid_change ()
 
void on_periodicity_change ()
 
void on_cell_structure_change ()
 
void on_thermostat_param_change ()
 
void on_temperature_change ()
 
void on_verlet_skin_change ()
 
void on_timestep_change ()
 
void on_integration_start ()
 
void on_short_range_ia_change ()
 
void on_non_bonded_ia_change ()
 
void on_coulomb_change ()
 
void on_dipoles_change ()
 
void on_constraint_change ()
 Called every time a constraint is changed.
 
void on_lb_boundary_conditions_change ()
 Called when the LB boundary conditions change (geometry, slip velocity, or both).
 
void on_particle_local_change ()
 Called every time a particle local property changes.
 
void on_particle_change ()
 Called every time a particle property changes.
 
void on_particle_charge_change ()
 Called every time a particle charge changes.
 
void on_observable_calc ()
 called before calculating observables, i.e.
 
void on_lees_edwards_change ()
 
void veto_boxl_change (bool skip_particle_checks=false) const
 

Static Public Member Functions

static std::shared_ptr< Systemcreate ()
 

Public Attributes

GpuParticleData gpu
 
ResourceCleanup cleanup_queue
 
Coulomb::Solver coulomb
 
Dipoles::Solver dipoles
 
LB::Solver lb
 
EK::Solver ek
 
std::shared_ptr< BoxGeometrybox_geo
 
std::shared_ptr< LocalBoxlocal_geo
 
std::shared_ptr< CellStructurecell_structure
 
std::shared_ptr< Propagationpropagation
 
std::shared_ptr< BondedInteractionsMapbonded_ias
 
std::shared_ptr< InteractionsNonBondednonbonded_ias
 
std::shared_ptr< Thermostat::Thermostatthermostat
 
std::shared_ptr< ComFixedcomfixed
 
std::shared_ptr< Galileigalilei
 
std::shared_ptr< OifGlobaloif_global
 
std::shared_ptr< ImmersedBoundariesimmersed_boundaries
 
std::shared_ptr< CollisionDetection::CollisionDetectioncollision_detection
 
std::shared_ptr< BondBreakage::BondBreakagebond_breakage
 
std::shared_ptr< LeesEdwards::LeesEdwardslees_edwards
 
std::shared_ptr< Accumulators::AutoUpdateAccumulatorsauto_update_accumulators
 
std::shared_ptr< Constraints::Constraintsconstraints
 
std::shared_ptr< SteepestDescentsteepest_descent
 
std::shared_ptr< StokesianDynamicsstokesian_dynamics
 
std::shared_ptr< NptIsoParametersnptiso
 
std::shared_ptr< InstantaneousPressurenpt_inst_pressure
 

Protected Member Functions

void update_local_geo ()
 
void update_icc_particles ()
 
bool has_icc_enabled () const
 
void integrate_magnetodynamics ()
 Run magnetodynamics update for local virtual particles.
 

Protected Attributes

bool reinit_thermo
 Whether the thermostat has to be reinitialized before integration.
 
double time_step
 Molecular dynamics integrator time step.
 
double sim_time
 Molecular dynamics integrator simulation time.
 
double force_cap
 Molecular dynamics integrator force capping.
 
double min_global_cut
 Minimal global interaction cutoff.
 

Detailed Description

Main system class.

Most components follow the composite pattern and the opaque pointer pattern. See System class design for more details.

Definition at line 81 of file core/system/System.hpp.

Constructor & Destructor Documentation

◆ System()

System::System::System ( Private  )

Definition at line 70 of file core/system/System.cpp.

References inactive_cutoff, and kokkos_handle.

Member Function Documentation

◆ calculate_energy()

std::shared_ptr< Observable_stat > System::System::calculate_energy ( )

◆ calculate_forces()

◆ calculate_long_range_fields()

void System::System::calculate_long_range_fields ( )

Calculate dipole fields.

◆ calculate_pressure()

◆ check_kT()

void System::System::check_kT ( double  value) const

Veto temperature change.

Definition at line 152 of file core/system/System.cpp.

◆ create()

std::shared_ptr< System > System::System::create ( )
static

Definition at line 64 of file core/system/System.cpp.

References stream.

Referenced by ScriptInterface::System::System::do_construct().

◆ get_force_cap()

auto System::System::get_force_cap ( ) const
inline

Get force_cap.

Definition at line 109 of file core/system/System.hpp.

◆ get_global_ghost_flags()

unsigned System::System::get_global_ghost_flags ( ) const

Returns the ghost flags required for running pair kernels for the global state, e.g.

the force calculation.

Returns
Required data parts;

Definition at line 540 of file core/system/System.cpp.

References Cells::DATA_PART_BONDS, Cells::DATA_PART_MOMENTUM, Cells::DATA_PART_POSITION, Cells::DATA_PART_PROPERTIES, stream, THERMO_BOND, and THERMO_DPD.

◆ get_instantaneous_pressure()

double System::System::get_instantaneous_pressure ( )

get the instantaneous pressure with (q(t+dt), p(t+dt/2))

◆ get_instantaneous_pressure_virial()

double System::System::get_instantaneous_pressure_virial ( )

get the instantaneous virial pressure with q(t+dt)

◆ get_interaction_range()

double System::System::get_interaction_range ( ) const

Get the interaction range.

Definition at line 472 of file core/system/System.cpp.

References inactive_cutoff, and stream.

◆ get_min_global_cut()

auto System::System::get_min_global_cut ( ) const
inline

Get min_global_cut.

Definition at line 115 of file core/system/System.hpp.

◆ get_npt_virial()

Utils::Vector3d * System::System::get_npt_virial ( ) const

Definition at line 571 of file core/system/System.cpp.

◆ get_sim_time()

auto System::System::get_sim_time ( ) const
inline

Get sim_time.

Definition at line 103 of file core/system/System.hpp.

◆ get_time_step()

auto System::System::get_time_step ( ) const
inline

Get time_step.

Definition at line 97 of file core/system/System.hpp.

◆ has_collision_detection_enabled()

bool System::System::has_collision_detection_enabled ( ) const

Definition at line 581 of file core/system/System.cpp.

References stream.

◆ has_icc_enabled()

bool System::System::has_icc_enabled ( ) const
protected

Definition at line 320 of file icc.cpp.

References stream.

◆ has_npt_enabled()

bool System::System::has_npt_enabled ( ) const

◆ integrate()

int System::System::integrate ( int  n_steps,
int  reuse_forces 
)

Integrate equations of motion.

Parameters
n_stepsNumber of integration steps, can be zero
reuse_forcesDecide when to re-calculate forces

This function calls two hooks for propagation kernels such as velocity verlet, velocity verlet + npt box changes, and steepest_descent. One hook is called before and one after the force calculation. It is up to the propagation kernels to increment the simulation time.

This function propagates the system according to the choice of integrator stored in Propagation::integ_switch. The general structure is:

  • if reuse_forces is zero, recalculate the forces based on the current state of the system
  • Loop over the number of simulation steps:
    1. initialization (e.g., RATTLE)
    2. First hook for propagation kernels
    3. Update dependent particles and properties (RATTLE, virtual sites)
    4. Calculate forces for the current state of the system. This includes forces added by the Langevin thermostat and the Lattice-Boltzmann-particle coupling
    5. Second hook for propagation kernels
    6. Update dependent properties (Virtual sites, RATTLE)
    7. Run single step algorithms (Lattice-Boltzmann propagation, collision detection, NpT update)
  • Final update of dependent properties and statistics/counters

High-level documentation of the integration and thermostatting schemes can be found in doc/sphinx/system_setup.rst and /doc/sphinx/running.rst

Returns
number of steps that have been integrated, or a negative error code

Definition at line 495 of file integrate.cpp.

References check_runtime_errors(), comm_cart, convert_initial_torques(), correct_position_shake(), correct_velocity_shake(), Cells::DATA_PART_PROPERTIES, INTEG_ERROR_RUNTIME, INTEG_ERROR_SIGINT, INTEG_METHOD_BD, INTEG_METHOD_STEEPEST_DESCENT, INTEG_REUSE_FORCES_ALWAYS, INTEG_REUSE_FORCES_NEVER, integrator_step_1(), integrator_step_2(), lb_tracers_add_particle_force_to_fluid(), lb_tracers_propagate(), LEES_EDWARDS, Cells::RESORT_LOCAL, resort_particles_if_needed(), PropagationMode::ROT_VS_INDEPENDENT, PropagationMode::ROT_VS_RELATIVE, runtimeErrorMsg, save_old_position(), stream, PropagationMode::TRANS_LB_MOMENTUM_EXCHANGE, PropagationMode::TRANS_LB_TRACER, PropagationMode::TRANS_VS_CENTER_OF_MASS, PropagationMode::TRANS_VS_RELATIVE, vs_com_update_particles(), and vs_relative_update_particles().

◆ integrate_magnetodynamics()

void System::System::integrate_magnetodynamics ( )
protected

Run magnetodynamics update for local virtual particles.

Iterate over local particles and update the dipole moment of virtual particles according to the thermal Stoner-Wohlfarth model. Collect active homogeneous external magnetic fields from constraints and add the per-particle dipolar contribution before performing either the simplified no-field update or the full thermal Stoner-Wohlfarth update.

Definition at line 257 of file stoner_wohlfarth_thermal.cpp.

References Particle::dip_fld(), get_external_field(), get_reference_particle(), Particle::id(), Particle::is_virtual(), Particle::stoner_wohlfarth_is_enabled(), stoner_wohlfarth_main(), stoner_wohlfarth_no_field(), stream, THERMO_LANGEVIN, and Utils::uniform().

◆ integrate_with_signal_handler()

int System::System::integrate_with_signal_handler ( int  n_steps,
int  reuse_forces,
bool  update_accumulators 
)

◆ lb_couple_particles()

void System::System::lb_couple_particles ( )

Calculate particle-lattice interactions.

Definition at line 350 of file particle_coupling.cpp.

References LB::is_tracer(), LB::lb_coupling_sanity_checks(), runtimeErrorMsg, and stream.

◆ long_range_interactions_sanity_checks()

bool System::System::long_range_interactions_sanity_checks ( ) const

Check electrostatic and magnetostatic methods are properly initialized.

Returns
true if sanity checks failed.

Definition at line 457 of file core/system/System.cpp.

References runtimeErrorMsg, and stream.

◆ maximal_cutoff()

double System::System::maximal_cutoff ( ) const

Calculate the maximal cutoff of all interactions.

Definition at line 439 of file core/system/System.cpp.

References communicator, and inactive_cutoff.

◆ npt_add_virial_contribution()

◆ npt_ensemble_init()

void System::System::npt_ensemble_init ( bool  recalc_forces)

Reinitialize the NpT state.

Definition at line 100 of file npt.cpp.

References comm_cart, Utils::Mpi::gather_buffer(), Utils::Vector< T, N >::size(), and this_node.

◆ on_boxl_change()

void System::System::on_boxl_change ( bool  skip_method_adaption = false)

Called when the box length has changed.

This routine is relatively fast, and changing the box length every time step as for example necessary for NpT is more or less ok.

Parameters
skip_method_adaptionskip the long-range methods adaptions

Definition at line 201 of file core/system/System.cpp.

References stream.

◆ on_cell_structure_change()

void System::System::on_cell_structure_change ( )

Definition at line 265 of file core/system/System.cpp.

References clear_particle_node().

◆ on_constraint_change()

void System::System::on_constraint_change ( )

Called every time a constraint is changed.

Definition at line 327 of file core/system/System.cpp.

◆ on_coulomb_change()

void System::System::on_coulomb_change ( )

◆ on_dipoles_change()

void System::System::on_dipoles_change ( )

◆ on_integration_start()

◆ on_lb_boundary_conditions_change()

void System::System::on_lb_boundary_conditions_change ( )

Called when the LB boundary conditions change (geometry, slip velocity, or both).

Definition at line 329 of file core/system/System.cpp.

Referenced by ScriptInterface::walberla::LBFluid::do_call_method().

◆ on_lees_edwards_change()

void System::System::on_lees_edwards_change ( )

Definition at line 431 of file core/system/System.cpp.

◆ on_node_grid_change()

void System::System::on_node_grid_change ( )

Definition at line 233 of file core/system/System.cpp.

◆ on_non_bonded_ia_change()

void System::System::on_non_bonded_ia_change ( )

Definition at line 306 of file core/system/System.cpp.

◆ on_observable_calc()

void System::System::on_observable_calc ( )

called before calculating observables, i.e.

energy, pressure or the integrator (forces). Initialize any methods here which are not initialized immediately (P3M etc.).

Definition at line 391 of file core/system/System.cpp.

References clear_particle_node().

◆ on_particle_change()

◆ on_particle_charge_change()

void System::System::on_particle_charge_change ( )

Called every time a particle charge changes.

Definition at line 359 of file core/system/System.cpp.

◆ on_particle_local_change()

void System::System::on_particle_local_change ( )

Called every time a particle local property changes.

Definition at line 333 of file core/system/System.cpp.

Referenced by ReactionMethods::ReactionAlgorithm::make_reaction_attempt().

◆ on_periodicity_change()

void System::System::on_periodicity_change ( )

Definition at line 246 of file core/system/System.cpp.

References INTEG_METHOD_SD, runtimeErrorMsg, and stream.

◆ on_short_range_ia_change()

void System::System::on_short_range_ia_change ( )

Definition at line 301 of file core/system/System.cpp.

◆ on_temperature_change()

void System::System::on_temperature_change ( )

Definition at line 290 of file core/system/System.cpp.

◆ on_thermostat_param_change()

void System::System::on_thermostat_param_change ( )

◆ on_timestep_change()

void System::System::on_timestep_change ( )

Definition at line 295 of file core/system/System.cpp.

◆ on_verlet_skin_change()

void System::System::on_verlet_skin_change ( )

Definition at line 279 of file core/system/System.cpp.

◆ particle_bond_energy()

std::optional< double > System::System::particle_bond_energy ( int  pid,
int  bond_id,
std::vector< int partners 
)

Compute the energy of a given bond which has to exist on the given particle.

Requires that bond partners are visible on the same MPI rank as the primary particle. Returns nothing if the primary particle is not owned by this MPI rank.

Parameters
pidParticle id
bond_idBond id
partnersParticle ids of the bond partners
Returns
energy of the bond given the primary particle and bond partners

Definition at line 148 of file energy.cpp.

References bond_broken_error(), calc_bonded_energy(), get_ptr(), Particle::id(), Particle::is_ghost(), and stream.

◆ particle_short_range_energy_contribution()

double System::System::particle_short_range_energy_contribution ( int  pid)

Compute the short-range energy of a particle.

Iterate through particles inside cell and neighboring cells and compute energy contribution for a specific particle.

Parameters
pidParticle id
Returns
Non-bonded energy of the particle.

Definition at line 122 of file energy.cpp.

References calc_non_bonded_pair_energy(), do_nonbonded(), get_ptr(), and stream.

◆ rebuild_aosoa()

void System::System::rebuild_aosoa ( )

Definition at line 412 of file core/system/System.cpp.

References inactive_cutoff, stream, and update_cabana_state().

◆ rebuild_cell_structure()

void System::System::rebuild_cell_structure ( )

Rebuild cell lists.

Use e.g. after a skin change.

Definition at line 197 of file core/system/System.cpp.

◆ set_box_l()

void System::System::set_box_l ( Utils::Vector3d const box_l)

Change the box dimensions.

Definition at line 479 of file core/system/System.cpp.

References stream.

◆ set_cell_structure_topology()

void System::System::set_cell_structure_topology ( CellStructureType  topology)

Change cell structure topology.

Definition at line 171 of file core/system/System.cpp.

References HYBRID, NSQUARE, REGULAR, and stream.

◆ set_force_cap()

void System::System::set_force_cap ( double  value)

Set force_cap.

Definition at line 161 of file core/system/System.cpp.

◆ set_min_global_cut()

void System::System::set_min_global_cut ( double  value)

Set min_global_cut.

Definition at line 166 of file core/system/System.cpp.

◆ set_sim_time()

void System::System::set_sim_time ( double  value)

Set sim_time.

Definition at line 891 of file integrate.cpp.

◆ set_time_step()

void System::System::set_time_step ( double  value)

Set time_step.

Definition at line 139 of file core/system/System.cpp.

◆ synchronize_npt_state()

void System::System::synchronize_npt_state ( )

Synchronize NpT state such as instantaneous and average pressure.

Definition at line 43 of file npt.cpp.

References comm_cart.

◆ tune_verlet_skin()

void System::System::tune_verlet_skin ( double  min_skin,
double  max_skin,
double  tol,
int  int_steps,
bool  adjust_max_skin 
)

Tune the Verlet skin.

Choose the optimal Verlet list skin between min_skin and max_skin by bisection to tolerance tol.

Definition at line 126 of file tuning.cpp.

References stream, and time_calc().

◆ update_dependent_particles()

void System::System::update_dependent_particles ( )

Update particles with properties depending on other particles, namely virtual sites and ICC charges.

Definition at line 365 of file core/system/System.cpp.

References vs_com_update_particles(), and vs_relative_update_particles().

◆ update_icc_particles()

void System::System::update_icc_particles ( )
protected

Definition at line 326 of file icc.cpp.

References get_ptr(), and stream.

◆ update_local_geo()

void System::System::update_local_geo ( )
protected

◆ update_used_propagations()

void System::System::update_used_propagations ( )

Update the global propagation bitmask.

Definition at line 185 of file integrate.cpp.

References comm_cart, PropagationMode::NONE, and PropagationMode::SYSTEM_DEFAULT.

◆ veto_boxl_change()

void System::System::veto_boxl_change ( bool  skip_particle_checks = false) const

Definition at line 219 of file core/system/System.cpp.

References comm_cart, and stream.

Member Data Documentation

◆ auto_update_accumulators

std::shared_ptr<Accumulators::AutoUpdateAccumulators> System::System::auto_update_accumulators

Definition at line 336 of file core/system/System.hpp.

◆ bond_breakage

std::shared_ptr<BondBreakage::BondBreakage> System::System::bond_breakage

Definition at line 333 of file core/system/System.hpp.

◆ bonded_ias

std::shared_ptr<BondedInteractionsMap> System::System::bonded_ias

◆ box_geo

std::shared_ptr<BoxGeometry> System::System::box_geo

Definition at line 319 of file core/system/System.hpp.

Referenced by Constraints::Constraints::add_energy(), Constraints::Constraints::add_forces(), DipolarP3MHeffte< FloatType, Architecture, FFTConfig >::calc_average_self_energy_k_space(), DipolarP3MHeffte< FloatType, Architecture, FFTConfig >::calc_energy_correction(), CoulombP3MHeffte< FloatType, Architecture, FFTConfig >::calc_influence_function_energy(), DipolarP3MHeffte< FloatType, Architecture, FFTConfig >::calc_influence_function_energy(), CoulombP3MHeffte< FloatType, Architecture, FFTConfig >::calc_influence_function_force(), DipolarP3MHeffte< FloatType, Architecture, FFTConfig >::calc_influence_function_force(), DipolarTuningAlgorithm< FloatType, Architecture, FFTConfig >::calculate_accuracy(), TuningAlgorithm::commit(), PairCriteria::DistanceCriterion::decide(), PairCriteria::EnergyCriterion::decide(), TuningAlgorithm::determine_r_cut_limits(), ScriptInterface::walberla::LBFluid::do_call_method(), ScriptInterface::walberla::LatticeWalberla::do_construct(), Observables::BondAngles::evaluate(), Observables::BondDihedrals::evaluate(), Observables::CosPersistenceAngles::evaluate(), Observables::CylindricalDensityProfile::evaluate(), Observables::CylindricalFluxDensityProfile::evaluate(), Observables::CylindricalVelocityProfile::evaluate(), Observables::DensityProfile::evaluate(), Observables::FluxDensityProfile::evaluate(), Observables::ForceDensityProfile::evaluate(), Observables::ParticleDistances::evaluate(), LB::Solver::get_interpolated_densities(), LB::Solver::get_interpolated_density(), LB::Solver::get_interpolated_velocity(), TuningAlgorithm::get_mc_time(), ReactionMethods::ReactionAlgorithm::get_random_position_in_box(), CoulombP3MHeffte< FloatType, Architecture, FFTConfig >::long_range_pressure(), LB::Solver::make_lattice_position_checker(), maybe_insert_particle(), pack_particles(), ScriptInterface::Particles::ParticleHandle::ParticleHandle(), ParticleObservables::traits< Particle >::position(), CoulombP3MHeffte< FloatType, Architecture, FFTConfig >::scaleby_box_l(), DipolarP3MHeffte< FloatType, Architecture, FFTConfig >::scaleby_box_l(), ReactionMethods::ReactionAlgorithm::set_cyl_constraint(), ReactionMethods::ReactionAlgorithm::set_slab_constraint(), CoulombTuningAlgorithm< FloatType, Architecture, FFTConfig >::setup_logger(), DipolarTuningAlgorithm< FloatType, Architecture, FFTConfig >::setup_logger(), and ReactionMethods::ReactionAlgorithm::update_volume().

◆ cell_structure

◆ cleanup_queue

ResourceCleanup System::System::cleanup_queue

Definition at line 94 of file core/system/System.hpp.

◆ collision_detection

std::shared_ptr<CollisionDetection::CollisionDetection> System::System::collision_detection

Definition at line 331 of file core/system/System.hpp.

◆ comfixed

std::shared_ptr<ComFixed> System::System::comfixed

Definition at line 326 of file core/system/System.hpp.

◆ constraints

std::shared_ptr<Constraints::Constraints> System::System::constraints

Definition at line 337 of file core/system/System.hpp.

◆ coulomb

◆ dipoles

Dipoles::Solver System::System::dipoles

Definition at line 316 of file core/system/System.hpp.

◆ ek

◆ force_cap

double System::System::force_cap
protected

Molecular dynamics integrator force capping.

Definition at line 356 of file core/system/System.hpp.

◆ galilei

std::shared_ptr<Galilei> System::System::galilei

Definition at line 327 of file core/system/System.hpp.

◆ gpu

◆ immersed_boundaries

std::shared_ptr<ImmersedBoundaries> System::System::immersed_boundaries

Definition at line 329 of file core/system/System.hpp.

◆ lb

◆ lees_edwards

std::shared_ptr<LeesEdwards::LeesEdwards> System::System::lees_edwards

Definition at line 334 of file core/system/System.hpp.

◆ local_geo

std::shared_ptr<LocalBox> System::System::local_geo

◆ min_global_cut

double System::System::min_global_cut
protected

Minimal global interaction cutoff.

Particles with a distance smaller than this are guaranteed to be available on the same node (through ghosts).

Definition at line 362 of file core/system/System.hpp.

◆ nonbonded_ias

◆ npt_inst_pressure

std::shared_ptr<InstantaneousPressure> System::System::npt_inst_pressure

Definition at line 344 of file core/system/System.hpp.

◆ nptiso

std::shared_ptr<NptIsoParameters> System::System::nptiso

◆ oif_global

std::shared_ptr<OifGlobal> System::System::oif_global

Definition at line 328 of file core/system/System.hpp.

◆ propagation

◆ reinit_thermo

bool System::System::reinit_thermo
protected

Whether the thermostat has to be reinitialized before integration.

Definition at line 350 of file core/system/System.hpp.

◆ sim_time

double System::System::sim_time
protected

Molecular dynamics integrator simulation time.

Definition at line 354 of file core/system/System.hpp.

◆ steepest_descent

std::shared_ptr<SteepestDescent> System::System::steepest_descent

◆ stokesian_dynamics

std::shared_ptr<StokesianDynamics> System::System::stokesian_dynamics

Definition at line 340 of file core/system/System.hpp.

◆ thermostat

std::shared_ptr<Thermostat::Thermostat> System::System::thermostat

Definition at line 325 of file core/system/System.hpp.

◆ time_step

double System::System::time_step
protected

Molecular dynamics integrator time step.

Definition at line 352 of file core/system/System.hpp.


The documentation for this class was generated from the following files: