# Marshaling for VTK
The `vtkObjectManager` synchronizes VTK objects that may or may not reside in the same process. For an object to be synchronized, these conditions must be met:
1. The type of an object must be derived from `vtkObjectBase`.
2. A serialization handler `std::function<vtknlohmann::json(vtkObjectBase*, vtkSerializer*)>` must be registered with the serializer. See [Marshaling with VTK CMake Module API](#Marshaling-with-VTK-CMake-Module-API)
3. A deserialization handler of type `std::function<void(const vtknlohmann::json&, vtkObjectBase*, vtkDeserializer*)>` must be registered with the deserializer. See [Marshaling with VTK CMake Module API](#Marshaling-with-VTK-CMake-Module-API)
4. A constructor of type `std::function<vtkObjectBase*()>` must be registered with the deserializer.
## Object manager
`vtkObjectManager` maintains internal instances of `vtkSerializer` and a `vtkDeserializer` to serialize and deserialize VTK objects respectively.
### Serialization
The vtkObjectManager facilitates serialization of objects by registering them, updating their state, and providing methods to retrieve both the serialized data (blobs) and object states based on their unique identifiers.
```python
manager: vtkObjectManager
manager.Initialize()
vtkObject0--|
|
|---> manager.RegisterObject(vtkObject0) -> id0
vtkObject1--|
|
|---> manager.RegisterObject(vtkObject1) -> id1
......
manager.Update() -> [id0, id1, id2, ...]---|
| |
| |
| |---> manager.GetBlobHashes(ids) -> ["hash0", "hash1", ...]
| |
| |
<foreach> <foreach>
| |
| |---> blob = manager.GetBlob(hash)
|
|---> state = manager.GetState(id)
```
### Deserialization
The vtkObjectManager also facilitates deserialization of objects by registering their states and data (blobs) and constructing or updating VTK objects based on MTime.
```python
manager: vtkObjectManager
manager.Initialize()
states : [str]
|
<foreach>
|---> manager.RegisterState(state) -> id
hashesAndBlobs : dict<str, str>
|
<foreach>
|---> manager.RegisterBlob(hash, blob) -> id
manager.Update() -> [id0, id1, id2, ...]
```
The object manager derives from `vtkObject`. It has the convenience of autogenerated python wrappers.
### vtkObjectManager API
#### Initialization
- `bool Initialize()`:
Loads the default (de)serialization handlers and constructors for VTK classes. See [Marshaling with VTK CMake Module API](#Marshaling-with-VTK-CMake-Module-API)
- Returns true on success, false otherwise.
- `bool InitializeExtensionModuleHandlers(vector<RegistrarType>)`:
Loads user provided handlers. See [Marshaling with VTK CMake Module API](#Marshaling-with-VTK-CMake-Module-API)
- Returns true on success, false otherwise.
#### Objects and states
- `vtkIdType RegisterObject(vtkObjectBase* object)`
Adds `object` into an internal container and returns a unique identifier.
- The identifier can be used in any of the methods that accept `id` or a vector of `id`.
- Returns an integer > 0 if object was successfully registered.
- Returns -1 in all other cases.
**Note:** *Objects registered here are immune to `PruneUnusedObjects()`. The only way to remove these objects is to invoke `UnRegisterObject(id)`*
- `vtkIdType RegisterState(const std::string& state)`
Adds `state` into an internal container and returns a unique identifier.
- The state
1. must be valid json.
2. must have a key-value pair `{'Id': n}` where n is an integer of type `vtkIdType`.
- The identifier can be used in any of the methods that accept `id` or a vector of `id`.
- Returns an integer > 0 if state was successfully parsed into json.
- Returns -1 in all other cases.
- `bool UnRegisterObject(vtkIdType id)`
Removes the `object` found at `id` and it's `state`.
- Returns `true` if the `object` and `state` was found and removed successfully.
- Returns `false` otherwise.
- `void Clear()`
Reset the object manager back to initial state.
- All registered objects are removed and no longer tracked.
- All registered states are also removed.
- `void PruneUnusedObjects()`
Unregisters objects whose identifier is not found in the state of any other objects.
- `vtkIdType GetId(vtkObjectBase* object)`
Get the identifier for `object`.
- Returns an integer > 0 if `object` was previously registered directly or indirectly i.e, as a dependency of another registered object.
- Returns -1 in all other cases.
- `std::string GetState(vtkIdType id)`
Get state of the object at `id`.
- Returns a non empty json valid string if an object registered directly or indirectly at `id` has a state.
- Returns an empty string in all other cases.
- `vtkObjectBase* GetObjectWithId(vtkIdType id)`
Get object at `id`.
- Returns a non-null pointer to a `vtkObjectBase` instance registered directly or indirectly at `id`.
- Returns `nullptr` in all other cases.
- `std::vector<vtkIdType> GetDirectDependencies(vtkIdType id)`
Get the direct dependencies of an object at `id`.
- Returns a non-empty vector of identifiers that correspond to the direct descendants of an object registered directly or indirectly at `id`.
- Returns an empty vector in all other cases.
- `std::vector<vtkIdType> GetAllDependencies(vtkIdType id) const`
Get all identifiers of objects that are dependencies of the given object
- Returns a non-empty vector of identifiers of all objects that depend on an object with the given identifier.
- Returns an empty vector if there are no dependents.
#### Data array management
- `std::vector<std::string> GetBlobHashes(const std::vector<vtkIdType>& ids)`
Get the hashes of blobs referenced by the objects at specified identifiers.
- Returns a non-empty vector of hash strings that correspond to blobs used by the registered objects at each identifier in `ids`.
- Returns an empty vector in all other cases.
- `vtkSmartPointer<vtkTypeUInt8Array> GetBlob(const std::string& hash) const`
Get a blob for a give hash.
- Returns a non-null pointer to `vtkTypeUInt8Array` that corresponds to the `hash`.
- Returns `nullptr` in all other cases.
- `bool RegisterBlob(const std::string& hash, vtkSmartPointer<vtkTypeUInt8Array> blob)`
Stores a `blob` that corresponds to `hash`. All object states that reference `hash` will resolve to `blob`.
- Returns `true` if the `blob` is valid and successfully registered
- Returns`false` in all other cases.
- `bool UnRegisterBlob(const std::string& hash)`
Removes a `blob` found at `hash`.
- Returns `true` if the `blob` was found and removed successfully.
- Returns `false` in all other cases.
- `void PruneUnusedBlobs()`
Removes all `blob`(s) whose `hash` is not found in the state of any object registered directly or indirectly.
#### Synchronization
- `std::vector<vtkIdType> GetIdsOfObjectsGreaterThanGivenMTime(vtkMTimeType mtime = 0)`
Returns identifiers of all registered objects whose MTime (modification time) is greater than `mtime`.
- `vtkMTimeType GetLatestMTimeFromState()`
Returns the latest MTime from the states of all objects.
- `vtkMTimeType GetLatestMTimeFromObjects()`
Returns the latest MTime from all objects.
- `std::vector<vtkIdType> Update()`
Synchronizes objects with the state based on their MTime.
- Returns a non-empty vector of identifiers that correspond to objects which were either serialized or deserialized.
- Returns an empty vector in all other cases.
## Marshaling with VTK CMake Module API
The new vtk-cmake function `vtk_module_marshal` is similar to `vtk_module_wrap_python` and auto-generates code for the handlers and class registration.
Serialization, deserialization and construction handlers can either be auto-generated or hand-written for VTK classes. It all boils down to defining two C++ functions and registering them with the serializer and deserializer inside one function with C linkage.
### Handlers
The handlers of a class must be registered inside a C function.
Here is an example of the autogenerated registrar for `vtkActor`.
```cpp
// vtkObjectManager::Initialize() calls handler registrar function for a library,
// which in turn invokes handlers for all classes.
// The object manager passes in pointers to a serializer and a deserializer.
// -> vtkObjectManager::Initialize()
// -> RegisterClasses_vtkModuleMarshaling(ser, deser)
// -> RegisterHandlers_vtkClassName1Marshalling(ser, deser)
// -> RegisterHandlers_vtkClassName2Marshalling(ser, deser)
// -> RegisterHandlers_vtkClassName3Marshalling(ser, deser)
// -> RegisterHandlers_vtkClassName4Marshalling(ser, deser)
// ...
extern "C" { int RegisterHandlers_vtkActorMarshaling(void* serializer, void* desererializer); }
int RegisterHandlers_vtkActorMarshaling(void* ser, void* deser)
{
int success = 0;
if (auto* asObjectBase = static_cast<vtkObjectBase*>(ser))
{
if (auto* serializer = vtkSerializer::SafeDownCast(asObjectBase))
{
serializer->RegisterHandler(typeid(vtkActor), Serialize_vtkActor);
success++;
}
}
if (auto* asObjectBase = static_cast<vtkObjectBase*>(deser))
{
if (auto* deserializer = vtkDeserializer::SafeDownCast(asObjectBase))
{
deserializer->RegisterHandler(typeid(vtkActor), Deserialize_vtkActor);
deserializer->RegisterConstructor("vtkActor", []() { return vtkActor::New(); });
success++;
}
}
return success;
}
```
#### Construction
This handler is a C++ function which returns a new `vtkObjectBase`. Here is an example of the auto-generated construction handler for `vtkActor`.
```cpp
deserializer->RegisterConstructor("vtkActor", []() { return vtkActor::New(); });
```
#### Serialization
This handler is a C++ function. It populates a `json` object with entries that make up the state of an object.
A typical serialization handler must initialize the state by calling the superclass's serialization handler.
The hierarchy information of the class must be encoded in two lists.
1. *SuperClasses* is a list of strings that correspond to the super class names. The first entry in this list is always `vtkObjectBase`. Handlers must append the superclass name into this list.
2. *SuperModules* is a list of strings that correspond to the modules in which the super classes exist. The first entry in this list is always "(null)". Handlers must append the supermodule name into this list. "(null)" is used if the super class is in the same module.
Here is a commented example of the auto-generated serialization handler for `vtkActor`
```cpp
static vtknlohmann::json Serialize_vtkActor(vtkObjectBase* objectBase, vtkSerializer* serializer)
{
using json = vtknlohmann::json;
json state;
auto object = vtkActor::SafeDownCast(objectBase);
// find the serialization handler of the Superclass.
if (auto f = serializer->GetHandler(typeid(vtkActor::Superclass)))
{
// Serialize all properties of vtkProp3D, which in turn
// will serialize all properties of vtkProp, vtkObject and vtkObjectBase.
state = f(object, serializer);
}
// Fill in hierarchy information obtained from vtkParseHierarchy.
// Our superclass is vtkProp3D.
state["SuperClassNames"].push_back("vtkProp3D");
// Our superclass is in the same module i.e, vtkRenderingCore.
state["SuperModuleNames"].push_back("(null)");
// Types that are acceptable in the constructor of basic_json<> can be copy-assigned.
state["RedrawMTime"] = object->GetRedrawMTime();
// When encountering a vtkObject, the auto generated code takes additional steps.
// First, check if the mapper is non-null.
if (object->GetMapper())
{
// serialize the mapper using `serializer`. it returns a simple identifier like '{"Id" : id}
state["Mapper"] = serializer->SerializeJSON(reinterpret_cast<void*>(object->GetMapper()));
}
// First, check if the property is non-null.
if (object->GetProperty())
{
// serialize the property using `serializer`. it returns a simple identifier like '{"Id" : id}
state["Property"] = serializer->SerializeJSON(reinterpret_cast<void*>(object->GetProperty()));
}
// First, check if the backface property is non-null.
if (object->GetBackfaceProperty())
{
// serialize the property using `serializer`. it returns a simple identifier like '{"Id" : id}
state["BackfaceProperty"] = serializer->SerializeJSON(reinterpret_cast<void*>(object->GetBackfaceProperty()));
}
// First, check if the texture is non-null.
if (object->GetTexture())
{
// serialize the texture using `serializer`. it returns a simple identifier like '{"Id" : id}
state["Texture"] = serializer->SerializeJSON(reinterpret_cast<void*>(object->GetTexture()));
}
return state;
}
```
#### Deserialization
This handler is a C++ function. It applies the state on a `vtkObjectBase`. It is always assumed that the correct type of `vtkObjectBase` was constructed using the construction handler.
A typical handler of this kind must invoke the super class's deserialization handler on the state and object before deserializing it's own properties.
Here is the auto-generated code for `vtkActor`.
```cpp
static void Deserialize_vtkActor(const vtknlohmann::json& state, vtkObjectBase* objectBase,vtkDeserializer* deserializer)
{
using json = vtknlohmann::json;
auto object = vtkActor::SafeDownCast(objectBase);
if (auto f = deserializer->GetHandler(typeid(vtkActor::Superclass)))
{
f(state, object, deserializer);
}
// deserialize object from state[id]
{
auto iter = state.find("Mapper");
if ((iter != state.end()) && !iter->is_null())
{
const auto indirection = deserializer->GetIndirection();
vtkObjectManagerObjectData data{iter->at("Id").get<vtkIdType>(), nullptr, {}};
indirection->GetObjectWithId(indirection->Context, data);
indirection->GetRegisteredState(indirection->Context, data);
if (data.Object == nullptr)
{
data.Object = deserializer->ConstructObject(data.State);
object->SetMapper(static_cast<vtkMapper*>(static_cast<void*>(data.Object)));
}
deserializer->DeserializeJSON(data.State, data.Object);
}
}
// deserialize object from state[id]
{
auto iter = state.find("Property");
if ((iter != state.end()) && !iter->is_null())
{
const auto indirection = deserializer->GetIndirection();
vtkObjectManagerObjectData data{iter->at("Id").get<vtkIdType>(), nullptr, {}};
indirection->GetObjectWithId(indirection->Context, data);
indirection->GetRegisteredState(indirection->Context, data);
if (data.Object == nullptr)
{
data.Object = deserializer->ConstructObject(data.State);
object->SetProperty(static_cast<vtkProperty*>(static_cast<void*>(data.Object)));
}
deserializer->DeserializeJSON(data.State, data.Object);
}
}
// deserialize object from state[id]
{
auto iter = state.find("BackfaceProperty");
if ((iter != state.end()) && !iter->is_null())
{
const auto indirection = deserializer->GetIndirection();
vtkObjectManagerObjectData data{iter->at("Id").get<vtkIdType>(), nullptr, {}};
indirection->GetObjectWithId(indirection->Context, data);
indirection->GetRegisteredState(indirection->Context, data);
if (data.Object == nullptr)
{
data.Object = deserializer->ConstructObject(data.State);
object->SetBackfaceProperty(static_cast<vtkProperty*>(static_cast<void*>(data.Object)));
}
deserializer->DeserializeJSON(data.State, data.Object);
}
}
// deserialize object from state[id]
{
auto iter = state.find("Texture");
if ((iter != state.end()) && !iter->is_null())
{
const auto indirection = deserializer->GetIndirection();
vtkObjectManagerObjectData data{iter->at("Id").get<vtkIdType>(), nullptr, {}};
indirection->GetObjectWithId(indirection->Context, data);
indirection->GetRegisteredState(indirection->Context, data);
if (data.Object == nullptr)
{
data.Object = deserializer->ConstructObject(data.State);
object->SetTexture(static_cast<vtkTexture*>(static_cast<void*>(data.Object)));
}
deserializer->DeserializeJSON(data.State, data.Object);
}
}
}
```