ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
LB::Solver Struct Reference

#include <Solver.hpp>

+ Inheritance diagram for LB::Solver:
+ Collaboration diagram for LB::Solver:

Classes

struct  Conversions
 
struct  Implementation
 

Public Member Functions

 Solver ()
 
bool is_solver_set () const
 Return true if a LB solver is active.
 
void reset ()
 Remove the LB solver.
 
template<typename LB , class... Args>
void set (Args... args)
 Set the LB solver.
 
template<class Connector >
void connect (Connector &&connector) const
 Connector to the implementation internal state.
 
void propagate ()
 Propagate the LB fluid.
 
void ghost_communication ()
 Perform a full ghost communication.
 
void ghost_communication_pdf ()
 Perform a ghost communication of the PDF field.
 
void ghost_communication_vel ()
 Perform a ghost communication of the velocity field.
 
void init () const
 Perform a full initialization of the lattice-Boltzmann system.
 
void sanity_checks () const
 Perform LB parameter and boundary velocity checks.
 
void veto_time_step (double time_step) const
 Check if a MD time step is compatible with the LB tau.
 
void veto_kT (double kT) const
 Check if a thermostat is compatible with the LB temperature.
 
void lebc_sanity_checks (unsigned int shear_direction, unsigned int shear_plane_normal) const
 Perform LB LEbc parameter checks.
 
bool is_gpu () const
 
double get_tau () const
 Get the LB time step.
 
double get_agrid () const
 Get the LB grid spacing.
 
double get_kT () const
 Get the thermal energy parameter of the LB fluid.
 
auto get_lattice_speed () const
 Get the lattice speed (agrid/tau).
 
Utils::VectorXd< 9 > get_pressure_tensor () const
 
Utils::Vector3d get_momentum () const
 
std::optional< Utils::Vector3dget_interpolated_velocity (Utils::Vector3d const &pos) const
 Calculate the interpolated fluid velocity in LB units.
 
std::optional< double > get_interpolated_density (Utils::Vector3d const &pos) const
 Calculate the interpolated fluid density in LB units.
 
Utils::Vector3d get_coupling_interpolated_velocity (Utils::Vector3d const &pos) const
 Calculate the interpolated fluid velocity in MD units.
 
std::vector< Utils::Vector3dget_coupling_interpolated_velocities (std::vector< Utils::Vector3d > const &pos) const
 
void add_forces_at_pos (std::vector< Utils::Vector3d > const &pos, std::vector< Utils::Vector3d > const &forces)
 
void add_force_density (Utils::Vector3d const &pos, Utils::Vector3d const &force_density)
 Add a force density to the fluid at the given position.
 
void on_boxl_change ()
 
void on_node_grid_change ()
 
void on_cell_structure_change ()
 
void on_timestep_change ()
 
void on_temperature_change ()
 
void on_lees_edwards_change ()
 
void update_collision_model ()
 
void veto_boxl_change () const
 
template<>
void set (std::shared_ptr< LBNone > lb_instance)
 
template<>
void set (std::shared_ptr< LBWalberlaBase > lb_fluid, std::shared_ptr< LBWalberlaParams > lb_params)
 
- Public Member Functions inherited from System::Leaf< Solver >
void bind_system (std::shared_ptr< System > const &system)
 
void detach_system (std::shared_ptr< System > const &system)
 

Additional Inherited Members

- Protected Member Functions inherited from System::Leaf< Solver >
auto & get_system ()
 
auto & get_system () const
 
- Protected Attributes inherited from System::Leaf< Solver >
std::weak_ptr< Systemm_system
 

Detailed Description

Definition at line 35 of file lb/Solver.hpp.

Constructor & Destructor Documentation

◆ Solver()

LB::Solver::Solver ( )

Definition at line 52 of file lb/Solver.cpp.

Member Function Documentation

◆ add_force_density()

void LB::Solver::add_force_density ( Utils::Vector3d const &  pos,
Utils::Vector3d const &  force_density 
)

Add a force density to the fluid at the given position.

Parameters
posPosition at which the force density is to be applied.
force_densityForce density to apply.

Definition at line 269 of file lb/Solver.cpp.

References LB::Solver::Conversions::force_to_lb, and LB::Solver::Conversions::pos_to_lb.

Referenced by lb_tracers_add_particle_force_to_fluid().

◆ add_forces_at_pos()

void LB::Solver::add_forces_at_pos ( std::vector< Utils::Vector3d > const &  pos,
std::vector< Utils::Vector3d > const &  forces 
)

◆ connect()

template<class Connector >
void LB::Solver::connect ( Connector &&  connector) const
inline

Connector to the implementation internal state.

For developers: use this mechanism to access the underlying variant.

Definition at line 61 of file lb/Solver.hpp.

Referenced by EK::EKWalberla::propagate().

◆ get_agrid()

double LB::Solver::get_agrid ( ) const

Get the LB grid spacing.

Definition at line 175 of file lb/Solver.cpp.

References LB::check_solver().

Referenced by get_lattice_speed(), LB::ParticleCoupling::kernel(), and lb_tracers_add_particle_force_to_fluid().

◆ get_coupling_interpolated_velocities()

std::vector< Utils::Vector3d > LB::Solver::get_coupling_interpolated_velocities ( std::vector< Utils::Vector3d > const &  pos) const

◆ get_coupling_interpolated_velocity()

Utils::Vector3d LB::Solver::get_coupling_interpolated_velocity ( Utils::Vector3d const &  pos) const

Calculate the interpolated fluid velocity in MD units.

Special method used only for particle coupling. Uses the LB ghost layer. Achieved by linear interpolation (eq. 11 in [2]).

Parameters
posPosition in MD units at which the velocity is to be calculated.
Return values
interpolatedfluid velocity.

Definition at line 222 of file lb/Solver.cpp.

References LB::Solver::Conversions::pos_to_lb, and LB::Solver::Conversions::vel_to_md.

Referenced by lb_drag_force(), and lb_tracers_propagate().

◆ get_interpolated_density()

std::optional< double > LB::Solver::get_interpolated_density ( Utils::Vector3d const &  pos) const

Calculate the interpolated fluid density in LB units.

Use this function in MPI-parallel code. The LB ghost layer is ignored.

Parameters
posPosition in MD units at which the density is to be calculated.
Return values
interpolatedfluid density.

Definition at line 211 of file lb/Solver.cpp.

References System::System::box_geo, System::get_system(), and LB::Solver::Conversions::pos_to_lb.

◆ get_interpolated_velocity()

std::optional< Utils::Vector3d > LB::Solver::get_interpolated_velocity ( Utils::Vector3d const &  pos) const

Calculate the interpolated fluid velocity in LB units.

Use this function in MPI-parallel code. The LB ghost layer is ignored.

Parameters
posPosition in MD units at which the velocity is to be calculated.
Return values
interpolatedfluid velocity.

Definition at line 197 of file lb/Solver.cpp.

References System::System::box_geo, System::get_system(), and LB::Solver::Conversions::pos_to_lb.

◆ get_kT()

double LB::Solver::get_kT ( ) const

Get the thermal energy parameter of the LB fluid.

Definition at line 185 of file lb/Solver.cpp.

References LB::check_solver().

Referenced by LB::ParticleCoupling::ParticleCoupling().

◆ get_lattice_speed()

auto LB::Solver::get_lattice_speed ( ) const
inline

◆ get_momentum()

Utils::Vector3d LB::Solver::get_momentum ( ) const

Definition at line 281 of file lb/Solver.cpp.

References LB::check_solver().

Referenced by calc_linear_momentum().

◆ get_pressure_tensor()

Utils::VectorXd< 9 > LB::Solver::get_pressure_tensor ( ) const

Definition at line 190 of file lb/Solver.cpp.

References LB::check_solver().

Referenced by Observables::LBFluidPressureTensor::operator()().

◆ get_tau()

double LB::Solver::get_tau ( ) const

Get the LB time step.

Definition at line 180 of file lb/Solver.cpp.

References LB::check_solver().

Referenced by get_lattice_speed(), and LB::ParticleCoupling::ParticleCoupling().

◆ ghost_communication()

void LB::Solver::ghost_communication ( )

Perform a full ghost communication.

Definition at line 76 of file lb/Solver.cpp.

References LB::check_solver().

◆ ghost_communication_pdf()

void LB::Solver::ghost_communication_pdf ( )

Perform a ghost communication of the PDF field.

Definition at line 81 of file lb/Solver.cpp.

References LB::check_solver().

◆ ghost_communication_vel()

void LB::Solver::ghost_communication_vel ( )

Perform a ghost communication of the velocity field.

Definition at line 86 of file lb/Solver.cpp.

References LB::check_solver().

◆ init()

void LB::Solver::init ( ) const
inline

Perform a full initialization of the lattice-Boltzmann system.

All derived parameters and the fluid are reset to their default values.

Definition at line 90 of file lb/Solver.hpp.

◆ is_gpu()

bool LB::Solver::is_gpu ( ) const

Definition at line 170 of file lb/Solver.cpp.

References LB::check_solver().

◆ is_solver_set()

bool LB::Solver::is_solver_set ( ) const

Return true if a LB solver is active.

Definition at line 64 of file lb/Solver.cpp.

References LB::is_solver_set().

Referenced by calc_linear_momentum(), and lb_sanity_checks().

◆ lebc_sanity_checks()

void LB::Solver::lebc_sanity_checks ( unsigned int  shear_direction,
unsigned int  shear_plane_normal 
) const

Perform LB LEbc parameter checks.

Definition at line 111 of file lb/Solver.cpp.

◆ on_boxl_change()

void LB::Solver::on_boxl_change ( )

Definition at line 134 of file lb/Solver.cpp.

◆ on_cell_structure_change()

void LB::Solver::on_cell_structure_change ( )

Definition at line 121 of file lb/Solver.cpp.

◆ on_lees_edwards_change()

void LB::Solver::on_lees_edwards_change ( )

Definition at line 158 of file lb/Solver.cpp.

◆ on_node_grid_change()

void LB::Solver::on_node_grid_change ( )

Definition at line 140 of file lb/Solver.cpp.

◆ on_temperature_change()

void LB::Solver::on_temperature_change ( )

Definition at line 152 of file lb/Solver.cpp.

◆ on_timestep_change()

void LB::Solver::on_timestep_change ( )

Definition at line 146 of file lb/Solver.cpp.

◆ propagate()

void LB::Solver::propagate ( )

Propagate the LB fluid.

Definition at line 71 of file lb/Solver.cpp.

References LB::check_solver().

◆ reset()

void LB::Solver::reset ( )

Remove the LB solver.

Definition at line 66 of file lb/Solver.cpp.

References System::get_system(), and System::System::lb.

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

◆ sanity_checks()

void LB::Solver::sanity_checks ( ) const

Perform LB parameter and boundary velocity checks.

Definition at line 91 of file lb/Solver.cpp.

References System::Leaf< Solver >::get_system().

◆ set() [1/3]

template<typename LB , class... Args>
void LB::Solver::set ( Args...  args)

Set the LB solver.

For developers: a specialization must exist for every new LB type.

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

◆ set() [2/3]

template<>
void LB::Solver::set ( std::shared_ptr< LBNone lb_instance)

Definition at line 287 of file lb/Solver.cpp.

◆ set() [3/3]

template<>
void LB::Solver::set ( std::shared_ptr< LBWalberlaBase lb_fluid,
std::shared_ptr< LBWalberlaParams lb_params 
)

Definition at line 295 of file lb/Solver.cpp.

References System::Leaf< Solver >::get_system().

◆ update_collision_model()

void LB::Solver::update_collision_model ( )

Definition at line 164 of file lb/Solver.cpp.

◆ veto_boxl_change()

void LB::Solver::veto_boxl_change ( ) const

Definition at line 128 of file lb/Solver.cpp.

◆ veto_kT()

void LB::Solver::veto_kT ( double  kT) const

Check if a thermostat is compatible with the LB temperature.

Definition at line 105 of file lb/Solver.cpp.

◆ veto_time_step()

void LB::Solver::veto_time_step ( double  time_step) const

Check if a MD time step is compatible with the LB tau.

Definition at line 98 of file lb/Solver.cpp.


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