ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
LatticeModel Class Referenceabstract

Abstract representation of a lattice-based model. More...

#include <LatticeModel.hpp>

+ Inheritance diagram for LatticeModel:
+ Collaboration diagram for LatticeModel:

Public Types

using units_map = std::unordered_map< std::string, double >
 

Public Member Functions

virtual ~LatticeModel ()=default
 
virtual LatticeWalberla const & get_lattice () const noexcept=0
 Get the underlying lattice.
 
std::shared_ptr< VTKHandlecreate_vtk (int delta_N, int initial_count, int flag_observables, units_map const &units_conversion, std::string const &identifier, std::string const &base_folder, std::string const &prefix)
 Create a VTK observable.
 
void write_vtk (std::string const &vtk_uid)
 Write a VTK observable to disk.
 
void switch_vtk (std::string const &vtk_uid, bool status)
 Toggle a VTK observable on/off.
 

Protected Member Functions

virtual void register_vtk_field_writers (walberla::vtk::VTKOutput &vtk_obj, units_map const &units_conversion, int flag_observables)=0
 Register VTK writers.
 
virtual void register_vtk_field_filters (walberla::vtk::VTKOutput &vtk_obj)=0
 
virtual void integrate_vtk_writers ()=0
 

Protected Attributes

std::map< std::string, std::shared_ptr< VTKHandle > > m_vtk_auto
 VTK writers that are executed automatically.
 
std::map< std::string, std::shared_ptr< VTKHandle > > m_vtk_manual
 VTK writers that are executed manually.
 

Detailed Description

Abstract representation of a lattice-based model.

Definition at line 31 of file walberla_bridge/include/walberla_bridge/LatticeModel.hpp.

Member Typedef Documentation

◆ units_map

using LatticeModel::units_map = std::unordered_map<std::string, double>

Constructor & Destructor Documentation

◆ ~LatticeModel()

virtual LatticeModel::~LatticeModel ( )
virtualdefault

Member Function Documentation

◆ create_vtk()

std::shared_ptr< VTKHandle > LatticeModel::create_vtk ( int  delta_N,
int  initial_count,
int  flag_observables,
units_map const &  units_conversion,
std::string const &  identifier,
std::string const &  base_folder,
std::string const &  prefix 
)

Create a VTK observable.

Parameters
delta_NWrite frequency, if 0 write a single frame, otherwise add a callback to write every delta_N EK steps to a new file
initial_countInitial execution count
flag_observablesWhich observables to measure (OR'ing of OutputVTK or EKOutputVTK values)
units_conversionLattice-to-MD units conversion
identifierName of the VTK dataset
base_folderPath to the VTK folder
prefixPrefix of the VTK files

Definition at line 31 of file LatticeModel.cpp.

References LatticeWalberla::get_blocks(), get_lattice(), m_vtk_auto, m_vtk_manual, register_vtk_field_filters(), and register_vtk_field_writers().

◆ get_lattice()

virtual LatticeWalberla const & LatticeModel::get_lattice ( ) const
pure virtualnoexcept

◆ integrate_vtk_writers()

virtual void LatticeModel::integrate_vtk_writers ( )
protectedpure virtual

◆ register_vtk_field_filters()

virtual void LatticeModel::register_vtk_field_filters ( walberla::vtk::VTKOutput &  vtk_obj)
protectedpure virtual

◆ register_vtk_field_writers()

virtual void LatticeModel::register_vtk_field_writers ( walberla::vtk::VTKOutput &  vtk_obj,
units_map const &  units_conversion,
int  flag_observables 
)
protectedpure virtual

Register VTK writers.

Use the multi-piece uniform grid format.

Implemented in walberla::EKinWalberlaImpl< FluxCount, FloatType >, and walberla::LBWalberlaImpl< FloatType, Architecture >.

Referenced by create_vtk().

◆ switch_vtk()

void LatticeModel::switch_vtk ( std::string const &  vtk_uid,
bool  status 
)

Toggle a VTK observable on/off.

Parameters
vtk_uidName of the VTK object
statustrue to switch on, false to switch off

Definition at line 81 of file LatticeModel.cpp.

References m_vtk_auto, and m_vtk_manual.

◆ write_vtk()

void LatticeModel::write_vtk ( std::string const &  vtk_uid)

Write a VTK observable to disk.

Parameters
vtk_uidName of the VTK object

Definition at line 69 of file LatticeModel.cpp.

References m_vtk_auto, and m_vtk_manual.

Member Data Documentation

◆ m_vtk_auto

std::map<std::string, std::shared_ptr<VTKHandle> > LatticeModel::m_vtk_auto
protected

◆ m_vtk_manual

std::map<std::string, std::shared_ptr<VTKHandle> > LatticeModel::m_vtk_manual
protected

VTK writers that are executed manually.

Definition at line 39 of file walberla_bridge/include/walberla_bridge/LatticeModel.hpp.

Referenced by create_vtk(), switch_vtk(), and write_vtk().


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