Skip to content

Software design

Jean-Noël Grad edited this page May 22, 2023 · 15 revisions

This page documents the way ESPResSo is designed and how it interfaces with libraries.

Table of Contents

MPI parallelism

In ESPResSo, most of the program control flow takes place on the process at MPI rank 0 (called the "head node"). The n-1 other MPI processes ("worker nodes") are idle and blocked by an infinite loop that waits for a message from the head node to trigger a function call. This is the reason why the Python package mpi4py doesn't work with ESPResSo.

This MPI callback mechanism is internally implemented using a visitor pattern that controls how the return value of the parallel tasks will be reduced on the head node. The Python interpreter has direct control over the script interface (via Cython bindings) on the head node. How the head node will distribute work to the worker nodes is influenced by the Python class MPI creation policy:

  • "GLOBAL" policy: the head node forwards all user input to the worker nodes immediately
    • execution on the head and worker nodes is almost always synchronous
    • this is the preferred policy when developing new features
  • "LOCAL" policy: the head node decides which tasks to distribute to the worker nodes
    • execution on the head and worker nodes is always asynchronous
    • worker nodes wait for a MPI message that instructs them to either fetch parameters from the head node, or call a function
    • this is the legacy policy
Sometimes, the worker nodes need to return early and end up back in the MPI blocking loop, so that the head node can call an asynchronous callback. This is needed for old code that hasn't been made fully MPI-parallel yet. Since the 4.2 release of ESPResSo, there aren't many places left where this is needed, and they are progressively being replaced by fully synchronous code. Also since 4.2, the script interface exposes the MPI communicator, such that boost::mpi operations can be carried out outside of the core.

Core

Hooks mechanism

The ESPResSo core is designed to be extremely modular. Algorithms implemented in the core will work independently from one another, however it is sometimes necessary to exchange information between two algorithms. For example, the particle propagator must communicate with the lattice-Boltzmann propagator, such that momentum can be exchanged between the particles and fluid, and electrostatic algorithms need to know about changes to the cell system to re-tune their parameters.

This communication is carried out using hooks that are triggered by "events", such as a change in the number of particles, cell system, box geometry or periodicity, or the introduction of numerical solvers for electrostatics, magnetostatics or hydrodynamics. These callbacks are responsible for adapting system-dependent parameters (e.g. long-range method cutoffs), raising an error that blocks an invalid change to the system and rolls back the state of the system (e.g. charge neutrality violation, box periodicity violation) or synchronizing the global state of a method across every MPI rank.

MPI parallel code

Functions in the script interface with "LOCAL" policy that change the global state of the system in the core are implemented using callbacks. Each callback needs to declare its reduction policy using a macro. The callback will typically run synchronously on the head and worker nodes, although there are rare situations where the callback will run synchronously on the worker nodes while the head node runs a different, yet related, callback.

Below is one possible implementation of a callback to calculate the total linear momentum of the system. The momentum calculation is written in the callback calc_linear_momentum_local(). It is registered as a callback and tagged for a reduction operation using the macro REGISTER_CALLBACK_REDUCTION(). The calc_linear_momentum() function is responsible for calling the callbacks on all ranks. This function is exported in the Python interface, where it can be called by the user on the head node.

#include "Particle.hpp"
#include "cells.hpp"
#include "communication.hpp"
#include <utils/Vector.hpp>

/** Callback to calculate the linear momentum of the cell on this node. */
static Utils::Vector3d calc_linear_momentum_local() {
  auto const particles = ::cell_structure.local_particles();
  auto const momentum =
      std::accumulate(particles.begin(), particles.end(), Utils::Vector3d{},
                      [](Utils::Vector3d const &acc, Particle const &p) {
                        return acc + p.mass() * p.v();
                      });
  return momentum;
}

REGISTER_CALLBACK_REDUCTION(calc_linear_momentum_local,
                            std::plus<Utils::Vector3d>())

/** Calculate the total linear momentum of the system. Head node only! */
Utils::Vector3d calc_linear_momentum() {
  return mpi_call(::Communication::Result::reduction,
                  std::plus<Utils::Vector3d>(),
                  calc_linear_momentum_local);
}

There are specialized REGISTER_*() macros and mpi_call_*() overloads to accommodate for a wide range of collective operations. See the doxygen documentation of header files MpiCallbacks.hpp and communication.hpp for more details.

Exception mechanism

Parallel functions

Functions called by MPI callbacks have two options. The first one is to throw an exception. In that case, the script interface caller should wrap the throwing callee within a context()->parallel_try_catch([&]() { /*...*/ } to capture the exception and communicate it to the head node. Example:

#include "script_interface/auto_parameters/AutoParameters.hpp"
#include "script_interface/get_value.hpp"

#include <memory>
#include <string>

namespace ScriptInterface::Feature {
class MyFeature : public AutoParameters<MyFeature> {
  std::shared_ptr<::MyFeatureInTheCore> m_actor;
public:
  void do_construct(VariantMap const &params) override {
    context()->parallel_try_catch([this, &params]() {
      m_actor = std::make_shared<::MyFeatureInTheCore>(
          get_value<double>(params, "prefactor"),
          get_value<bool>(params, "verbose"));
    });
  }
};
} // namespace ScriptInterface::Feature

The second option is to log a runtime error in the global runtime error queue. In that case, the Python caller must add a call to handle_errors() after the throwing callee, which will halt the flow of the program if the queue contains error messages. Example:

/* p3m.hpp */
void p3m_tune();

/* p3m.cpp */
#include "p3m.hpp"
#include "errorhandling.hpp"
#include "communication.hpp"

void p3m_tune_local() {
  auto const error = sanity_checks();
  if (error) {
    runtimeErrorMsg() << "sanity checks failed";
    return;
  }
  /* ... */
}

REGISTER_CALLBACK(p3m_tune_local)

void p3m_tune() {
  mpi_call_all(p3m_tune_local);
}
# electrostatics.pxd
cdef extern from "p3m.hpp":
    void p3m_tune()

# electrostatics.pyx
from . cimport electrostatics  # get symbols from the pxd file
from .utils cimport handle_errors

def tune():
    p3m_tune()
    # convert the logged error messages into a Python exception
    handle_errors("P3M tuning failed")

Non-parallel functions

Functions called exclusively on the head node should throw an exception. The corresponding import in the Python interface should decorate the symbol with except +, so that the C++ exception gets transformed into a Python exception. Example:

# sd.pxd
cdef extern from "sd.hpp":
    void set_sd_viscosity(double eta) except +
/* sd.hpp */
void set_sd_viscosity(double eta);

/* sd.cpp */
#include <stdexcept>
#include <string>

double sd_viscosity;
void set_sd_viscosity(double eta) {
  if (eta < 0.0)
    throw std::runtime_error("invalid viscosity: " + std::to_string(eta));
  sd_viscosity = eta;
}

CUDA functions

CUDA functions return an error code upon failure. This code can be converted to a runtime error using throw cuda_runtime_error_cuda(error_code), or alternatively via the convenience macro CUDA_CHECK().

/* cuda_init.cu */
#include "cuda_utils.cuh"
#include <cuda.h>

int cuda_get_n_gpus() {
  int deviceCount;
  CUDA_CHECK(cudaGetDeviceCount(&deviceCount))
  return deviceCount;
}

In an MPI environment, the error can be caught via its parent class cuda_runtime_error.

/* gpu_callbacks.cpp */
#include "cuda_init.hpp"
#include "cuda_utils.hpp"

int n_devices_local_node() {
  int n_devices;
  try {
    n_devices = cuda_get_n_gpus();
  } catch (cuda_runtime_error const &err) {
    n_devices = 0;
  }
  return n_devices;
}

Cython interface

Whenever possible, Cython should only check variable types and the core should check the variable values. This avoids unnecessary code duplication between the core and interface for range checking.

Exceptions need to be checked in the testsuite using with self.assertRaises(RuntimeError). The test case should also check the system isn't left in an undefined state, i.e. the exceptions are raised before values in the core are changed. See for example commit 8977e50.

Submodules

Internal libraries

Functionality that doesn't depend on the state of the system can be extracted from the EspressoCore target and stored in a submodule of src/. This is used for example by the EspressoUtils module, which is a dependency of EspressoCore.

Such submodules can adopt the Pitchfork layout sketched below. With this layout, passing the compiler flag -Isrc/library/include makes the public header file library/feature.hpp visible to consumers of that library while the implementation details in src/library/src/feature_impl.hpp remain hidden.

src/
 \--- library/
       |--- include/
       |     \--- library/
       |           \--- feature.hpp
       |--- src/
       |     |--- feature_impl.hpp
       |     \--- feature_impl.cpp
       \--- tests/
             \--- feature_test.cpp

External libraries

External libraries can be added to the project after approval from the core team.

CMake

Loading packages via environment variables

Many libraries and executables can be loaded from the system via environment variables, such as $PKG_CONFIG_PATH or $LD_LIBRARY_PATH. Several inclusion methods are available: find_package, find_program, pkg_check_modules.

Examples:

find_package(NumPy REQUIRED)
find_program(IPYTHON_EXECUTABLE NAMES jupyter ipython3 ipython)
pkg_check_modules(SCAFACOS scafacos)
Including subprojects via FetchContent

CMake-based projects that provide functionality to ESPResSo can be included using the FetchContent method. This approach enables pinning the library version and patching the source code.

Here is the minimum set of commands needed to populate the waLBerla library in ESPResSo:

if(ESPRESSO_BUILD_WITH_WALBERLA)
  include(FetchContent)
  FetchContent_Declare(
    walberla
    GIT_REPOSITORY https://i10git.cs.fau.de/walberla/walberla
    GIT_TAG        64a4d68
    PATCH_COMMAND  patch -p0 --ignore-whitespace -i
                   ${PROJECT_SOURCE_DIR}/cmake/waLBerlaFunctions.patch
  )
  FetchContent_Populate(walberla)
  set(WALBERLA_BUILD_TESTS off CACHE BOOL "")
  set(CMAKE_POSITION_INDEPENDENT_CODE on CACHE BOOL "")
  add_subdirectory("${walberla_SOURCE_DIR}" "${walberla_BINARY_DIR}")
  set(WALBERLA_LIBS walberla::core walberla::blockforest walberla::boundary walberla::field walberla::lbm)
endif()

The consumer targets in ESPResSo can then be updated with target_include_directories() to get visibility of the waLBerla header files during compilation, and with target_link_libraries() to get visibility of the waLBerla symbols during linking.

Here is a minimum working example that includes the NLopt library, links ESPResSo against it, introduces an external feature macro in config.hpp and adds an extra CMake option to activate the feature at configure time:

diff --git a/CMakeLists.txt b/CMakeLists.txt
index cab2060..3d3d2a8 100644
--- a/CMakeLists.txt
+++ b/CMakeLists.txt
@@ -93,4 +93,5 @@ option(ESPRESSO_BUILD_WITH_SCAFACOS "Build with ScaFaCoS support" OFF)
 option(ESPRESSO_BUILD_WITH_STOKESIAN_DYNAMICS "Build with Stokesian Dynamics"
        OFF)
+option(ESPRESSO_BUILD_WITH_NLOPT "Build with NLopt support" OFF)
 option(ESPRESSO_BUILD_BENCHMARKS "Enable benchmarks" OFF)
 option(ESPRESSO_BUILD_WITH_VALGRIND_MARKERS
@@ -335,4 +336,24 @@ if(ESPRESSO_BUILD_WITH_STOKESIAN_DYNAMICS)
 endif()
 
+if(ESPRESSO_BUILD_WITH_NLOPT)
+  include(FetchContent)
+  FetchContent_Declare(
+    nlopt
+    GIT_REPOSITORY https://github.com/stevengj/nlopt.git
+    GIT_TAG        f8692f86b5269112b3760704093d0600adb033e7
+  )
+  if(NOT nlopt_POPULATED)
+    FetchContent_Populate(nlopt)
+    set(BUILD_SHARED_LIBS on CACHE BOOL "")
+    set(NLOPT_FORTRAN off CACHE BOOL "")
+    set(NLOPT_PYTHON off CACHE BOOL "")
+    set(NLOPT_OCTAVE off CACHE BOOL "")
+    set(NLOPT_MATLAB off CACHE BOOL "")
+    set(NLOPT_GUILE off CACHE BOOL "")
+    set(NLOPT_SWIG off CACHE BOOL "")
+    add_subdirectory("${nlopt_SOURCE_DIR}" "${nlopt_BINARY_DIR}")
+  endif()
+endif()
+
 if(ESPRESSO_BUILD_WITH_VALGRIND_MARKERS)
   find_package(PkgConfig REQUIRED)
diff --git a/cmake/espresso_cmake_config.cmakein b/cmake/espresso_cmake_config.cmakein
index 884c323..78fe040 100644
--- a/cmake/espresso_cmake_config.cmakein
+++ b/cmake/espresso_cmake_config.cmakein
@@ -14,4 +14,6 @@
 #cmakedefine ESPRESSO_BUILD_WITH_STOKESIAN_DYNAMICS
 
+#cmakedefine ESPRESSO_BUILD_WITH_NLOPT
+
 #cmakedefine ESPRESSO_BUILD_WITH_VALGRIND_MARKERS
 
diff --git a/src/config/features.def b/src/config/features.def
index 10e462e..46925c0 100644
--- a/src/config/features.def
+++ b/src/config/features.def
@@ -121,3 +121,4 @@ SCAFACOS external
 GSL external
 STOKESIAN_DYNAMICS external
+NLOPT external
 VALGRIND_MARKERS external
diff --git a/src/core/CMakeLists.txt b/src/core/CMakeLists.txt
index 2d12bcb..d67c6ec 100644
--- a/src/core/CMakeLists.txt
+++ b/src/core/CMakeLists.txt
@@ -61,4 +61,8 @@ add_library(espresso::core ALIAS espresso_core)
 set_target_properties(espresso_core PROPERTIES CXX_CLANG_TIDY
                                                "${ESPRESSO_CXX_CLANG_TIDY}")
+if(ESPRESSO_BUILD_WITH_NLOPT)
+  target_include_directories(espresso_core PRIVATE "${nlopt_SOURCE_DIR}/src")
+  target_link_libraries(espresso_core PRIVATE nlopt)
+endif()
 
 if(ESPRESSO_BUILD_WITH_CUDA)
diff --git a/src/core/integrate.cpp b/src/core/integrate.cpp
index 1fddf73..c8714fd 100644
--- a/src/core/integrate.cpp
+++ b/src/core/integrate.cpp
@@ -57,4 +57,8 @@
 #include "virtual_sites.hpp"
 
+#ifdef NLOPT
+#include <api/nlopt.h>
+#endif
+
 #include <profiler/profiler.hpp>
 
@@ -74,4 +78,13 @@
 #endif
 
+#ifdef NLOPT
+void* munge_on_copy(void*) { return nullptr; }
+void* munge_on_destroy(void*) { return nullptr; }
+
+void foo() {
+  nlopt_set_munge(nlopt_opt{}, munge_on_destroy, munge_on_copy);
+}
+#endif
+
 int integ_switch = INTEG_METHOD_NVT;
 

The project is then patched and built as follows:

git apply nlopt_mwe.patch
mkdir build
cd build
cmake .. -D ESPRESSO_BUILD_WITH_NLOPT=ON
make -j
./pypresso ../testsuite/python/integrate.py

C++ bridge

For better separation of concerns, the external library should not expose its public interface to the EspressoCore directly, but instead rely on a private submodule. With this approach, the external library is a private dependency of a submodule, and that submodule is a private dependency of EspressoCore. The submodule should use a C++ design pattern to hide the details of the external library.

Using the Pitchfork layout outlined in section Internal libraries as an example, one could declare an abstract base class in feature.hpp and derive a concrete class in feature_impl.hpp. An instance of the class would be stored in EspressoCore, however EspressoCore wouldn't have access to the header files of the external library.

/* /usr/include/external_library/ExternalFeature.hpp */
class ExternalFeature {
public:
  ExternalFeature(int argument): m_argument(argument) {}
  double energy(int position) const;
private:
  int m_argument;
};

/* src/library/include/library/feature.hpp */
class FeatureBase {
public:
  virtual double get_energy(int position) const = 0;
};
FeatureBase *new_feature(int argument);

/* src/library/src/feature_impl.hpp */
#include <library/feature.hpp>
#include <ExternalFeature.hpp>

class Feature: public FeatureBase {
public:
  Feature(int argument): m_feature(ExternalFeature(argument)) {}
  double get_energy(int position) const override;
private:
  ExternalFeature m_feature;
};

/* src/library/src/feature_impl.cpp */
#include <library/feature.hpp>
#include "feature_impl.hpp"

FeatureBase *new_feature(int argument) {
  FeatureBase *instance = new Feature(argument);
  return instance;
}

double Feature::get_energy(int position) const {
  return m_feature.energy(position);
}

/* src/core/feature_instance.hpp */
#include <library/feature.hpp>

FeatureBase *feature();

/* src/core/feature_instance.cpp */
#include "communication.hpp"
#include "feature_instance.hpp"
#include <library/feature.hpp>
#include <stdexcept>

namespace {
FeatureBase *feature_instance = nullptr;
}

FeatureBase *feature() {
  if (!feature_instance)
    throw std::runtime_error("feature is uninitialized.");
  return feature_instance;
}

void init_feature_local(int argument) {
  feature_instance = new_feature(argument);
}

REGISTER_CALLBACK(init_feature_local)

void init_feature(int argument) {
  Communication::mpiCallbacks().call_all(init_feature_local, argument);
}

/* src/core/analysis.hpp */
double feature_energy(int position);

/* src/core/analysis.cpp */
#include "analysis.hpp"
#include "feature_instance.hpp"

double feature_energy(int position) {
  return feature()->get_energy(position);
}

Script interface

Using handles to core and script interface objects

Classes declared in the core can be instantiated in the core and managed by the script interface, or be directly instantiated in the script interface. This is achieved using "glue" classes in src/script_interface, which derive from AutoParameters. These classes set up getters and setters to members of the managed object (via add_parameters in the constructor) and provide a mechanism to forward function calls (via do_call_method()). Values are automatically converted from the python types (int, float, string, etc.) to the corresponding C++ types (int, double, std::string, etc.) and passed between the Python interface and script interface via a boost::variant.

Python classes derived from ScriptInterfaceHelper provide the Python interface. These classes name the AutoParameters class to instantiate (field _so_name), which will in turn instantiate the managed object. The Python class must also explicitly list which methods of the managed object are visible to the user (field _so_bind_methods). Python classes also control on which MPI nodes the AutoParameters object is instantiated (field _so_creation_policy): either on all nodes (GLOBAL creation policy, default) or only on the head node (LOCAL creation policy).

Here is a simplified UML diagram depicting how the ESPResSo core (and external libraries like waLBerla) can exchange data with the Python interpreter via Cython bindings (through the script interface):

classDiagram

class `espressomd.script_interface.ScriptInterfaceHelper` {
str _so_name
str _so_creation_policy
}
class `espressomd.script_interface.ScriptObjectList` {
str _so_name
str _so_creation_policy
}
class `espressomd.lb.LBFluidWalberla` {
str _so_name = "ScriptInterface::LBFluid"
str _so_creation_policy = "GLOBAL"
}

`espressomd.script_interface.ScriptInterfaceHelper` <|-- `espressomd.script_interface.ScriptObjectList`
`espressomd.script_interface.ScriptInterfaceHelper` <|-- `espressomd.lb.LBFluidWalberla`

class `ScriptInterface::ObjectHandle` {
do_construct(VariantMap parameters): void
do_call_method(std::string name, VariantMap parameters): Variant
}
class `ScriptInterface::ObjectList` {
std::vector[std::shared_ptr[Object]] objects
}
class `ScriptInterface::AutoParameters` {
std::unordered_map[std::string, AutoParameter] parameters
}

class `ScriptInterface::LBFluid` {
std::shared_ptr[espresso::LBFluid] lb_fluid
}

`ScriptInterface::AutoParameters` <|-- `ScriptInterface::LBFluid`
`espressomd.script_interface.ScriptObjectList` <--> `ScriptInterface::ObjectList` : cython
`espressomd.script_interface.ScriptInterfaceHelper` <--> `ScriptInterface::ObjectHandle` : cython
`espressomd.lb.LBFluidWalberla` <--> `ScriptInterface::LBFluid` : cython
`ScriptInterface::ObjectHandle` <|-- `ScriptInterface::AutoParameters`
`ScriptInterface::AutoParameters` <|-- `ScriptInterface::ObjectList`

class `espresso::LBFluid`
class `walberla::LatticeModel`

`walberla::LatticeModel` <|--  `espresso::LBFluid`
`ScriptInterface::LBFluid` *-- "1" `espresso::LBFluid`

Below is an example with observable classes. The Observable parent class has a shape() method that can be called directly from an instance of the observable. The AutoParameters class has an additional calculate() method that is not meant to be used directly by the user, but can be called from within the class with call_method("calculate"). This approach is used when the result of a function needs to be post-processed, in this case to reshape the flat array into an N-dimensional array. The script interface class instantiates the managed object and provides a method observable() that returns a shared pointer to it. This is used by Accumulator classes to get a handle to the observable they are accumulating data from, and optionally forward it to the core global variable auto_update_accumulators (a vector of shared pointers to script interface observables) when the user calls system.auto_update_accumulators.add(my_accumulator).

# observables.py
from .script_interface import ScriptInterfaceHelper, script_interface_register

@script_interface_register
class Observable(ScriptInterfaceHelper):
    """
    Base class for all observables.

    Methods
    -------
    shape()
        Return the shape of the observable.
    """
    _so_name = "Observables::Observable"
    _so_bind_methods = ("shape",)
    _so_creation_policy = "LOCAL"
    def calculate(self):
        return np.array(self.call_method("calculate")).reshape(self.shape())

@script_interface_register
class ComVelocity(Observable):
    """Calculates the center of mass velocity for particles with given ids."""
    _so_name = "Observables::ComVelocity"
/* src/script_interface/observables/initialize.hpp */
#include "script_interface/ObjectHandle.hpp"
#include <utils/Factory.hpp>

namespace ScriptInterface {
namespace Observables {
void initialize(Utils::Factory<ObjectHandle> *om);
}
}

/* src/script_interface/observables/initialize.cpp */
#include "script_interface/ScriptInterface.hpp"
#include "core/observables/Observable.hpp"
#include "core/observables/ComVelocity.hpp"
#include "initialize.hpp"
#include <utils/Factory.hpp>
#include <memory>

namespace ScriptInterface {
namespace Observables {

/** Base class for script interfaces to core Observables classes. */
class Observable : public ObjectHandle {
public:
  virtual std::shared_ptr<::Observables::Observable> observable() const = 0;
  Variant do_call_method(std::string const &method,
                         VariantMap const &parameters) override {
    if (method == "calculate")
      return observable()->operator()();
    if (method == "shape")
      return observable()->shape();
    return {};
  }
};

/** Base class for script interfaces to particle-based observables.
 *  @tparam CorePidObs Any particle-based observable core class.
 */
template <typename CorePidObs>
class PidObservable : public AutoParameters<PidObservable<CorePidObs>, Observable> {
public:
  PidObservable() {
    this->add_parameters({
      // AutoParameter object for CorePidObs::ids in the form {name, setter lambda, getter lambda}
      {"ids",
       AutoParameter::read_only,
       [this]() { return m_observable->ids(); }
      }
    });
  }
  void do_construct(VariantMap const &params) override {
    // instantiate the managed object in the script interface
    m_obs = make_shared_from_args<CorePidObs, std::vector<int>>(params, "ids");
  }
  std::shared_ptr<::Observables::Observable> observable() const override {
    return m_obs;
  }
private:
  // pointer to the managed object
  std::shared_ptr<CorePidObs> m_obs;
};

/** Register script interface bridge to core observable. */
void initialize(Utils::Factory<ObjectHandle> *om) {
  om->register_new<PidObservable<::Observables::ComVelocity>>("Observables::ComVelocity");
}

} /* namespace Observables */
} /* namespace ScriptInterface */

The AutoParameters framework

The AutoParameters framework provides tools to register members of the wrapped core class and expose them to the python script interface. It is possible to configure access by passing pointers to the core class setter/getter functions or by passing custom lambdas. Values can also be marked as read-only. This is achieved by constructing a list of AutoParameter objects in the add_parameters() function. The code sample below illustrates some of the AutoParameter constructor overloads (see the full list).

/* file src/core/MyClass.hpp */
class MyClass {
private:
  int m_id;
  double m_radius;
  double m_height;
  double m_phi;

public:
  int &get_id_as_reference() { return m_id; }
  void set_radius(double const &radius) { m_radius = radius; }
  double get_radius() const { return m_radius; }
  void set_phi(double phi) { m_phi = phi; }
  double get_phi() const { return m_phi; }
  double const &get_height_const_ref() const { return m_height; }
};

/* file src/script_interface/MyClass.hpp */
#include "core/MyClass.hpp"
#include "script_interface/auto_parameters/AutoParameters.hpp"
#include <cmath>
#include <memory>

namespace ScriptInterface {
class MyClass : public AutoParameters<MyClass> {
private:
  std::shared_ptr<::MyClass> m_obj;
  constexpr static double deg2rad = M_PI / 180.;

public:
  MyClass() {
    add_parameters({
        // class member 'id' has a getter that returns a reference
        {
           "id",
           m_obj->get_id_as_reference()
        },
        // class member 'radius' has a getter and setter
        {
           "radius",
            m_obj,
            &::MyClass::set_radius,
            &::MyClass::get_radius
        },
        // class member access by lambda (e.g. for unit conversion)
        {
           "phi", 
           [this](Variant const &v) { m_obj->set_phi(get_value<double>(v) * deg2rad); },
           [this]() { return m_obj->get_phi() / deg2rad; }
        },
        // immutable class member
        {
           "height",
           AutoParameter::read_only,
           [this]() { return m_obj->get_height_const_ref(); }
        }
     });
  }
};
} // namespace ScriptInterface

Calling core methods and functions

Script interface classes typically contain a do_call_method() method that maps strings to core functions. When core functions take parameters, they are converted from a variant type to the required type. If conversion fails, for example the user provided a string when an integer was expected, the script interface automatically throws an exception that gets converted to a Python exception. However, the script interface class (or the core class) must do range checks on these converted values.

Sometimes, the return value of the core function needs to be reduced, typically via a sum (e.g. for extensive properties of the system), or via an optional reduction (i.e. when only one MPI rank returns a value and all other ranks return an empty object).

#include "script_interface/ScriptInterface.hpp"
#include "script_interface/auto_parameters/AutoParameters.hpp"
#include <memory>
#include <stdexcept>

namespace ScriptInterface {
class MyClass : public AutoParameters<MyClass> {
  std::shared_ptr<Protocol> m_protocol;
public:
  Variant do_call_method(std::string const &name,
                         VariantMap const &params) override {
    if (name == "set_velocity") {
      if (not m_protocol) {
        throw std::invalid_argument("Set a protocol first");
      }
      m_protocol->protocol().set_velocity = get_value<int>(params, "velocity");
    } else if (name == "set_protocol") {
      context()->parallel_try_catch([this, &params]() {
        auto const protocol = params.at("protocol");
        if (is_none(protocol)) {
          do_set_parameter("protocol", protocol);
          return;
        }
        // check input arguments
        m_protocol = get_value<std::shared_ptr<Protocol>>(protocol);
        auto const direction = get_value<int>(params, "direction");
        auto const normal = get_value<int>(params, "normal");
        if (normal == direction) {
          throw std::invalid_argument("Parameters 'direction' and 'normal' must differ");
        }
        ::core_set_protocol(m_protocol->protocol());
        m_protocol->protocol().set(direction, normal);
        }
      });
    } else if (name == "get_total_force") {
      auto const force_local = m_protocol->get_force();
      return mpi_reduce_sum(context()->get_comm(), force_local);
    }
    return {};
  }
};
} // namespace ScriptInterface

Checkpointing

To enable checkpointing, Python classes need to define special methods to serialize and deserialize themselves. This is usually achieved by implementing adequate __getstate__(self) and __setstate__(self, params) methods. The result of the __getstate__ can be any data type, usually a list or dict, and will be written to a checkpoint file by the pickle module. When loading the system from a checkpoint file, pickle will read the value and pass it to the class as the params argument to the __setstate__ function to instantiate an object. During this instantiation, the __init__(self) method is not called.

Usually, the AutoParameters framework takes care of checkpointing directly. However, if the Python class holds non-regenerable data, it becomes necessary to implement the aforementioned serialization methods. When in doubt, look at how the other Python classes handle checkpointing.

Exporting symbols from the core directly in Cython

For historical reasons, some C++ functions and classes from the ESPResSo core are exported in the Python interface using src/python/espressomd/*.pxd files. These exported symbols can then be used within cdef functions in src/python/espressomd/*.pyx files. Example:

# interface.pxd
cdef extern from "cuda_init.hpp":
    int cuda_get_device()

# interface.pyx
from . cimport interface  # get symbols from the pxd file

cdef device(self):
    """Get device id."""
    dev = cuda_get_device()
    if dev == -1:
        raise Exception("cuda device get error")
    return dev

Likewise, struct and class declarations can be copy-pasted in .pxd files to allow returning a copy of a core object in a Cython function.

This is not the preferred way of exposing core functionality to the Python interface. In ESPResSo 4.2, most of this code was removed, and the remaining Cython code will progressively be replaced by regular Python code.