ROOT TTree objects passing between Narval Actors¶
Thanks for this documentation to Gennaro Tortone. You can find the original article here
Introduction¶
Narval framework is a powerful tool to design a Data Acquisition Software as a composition of modules (actors) where the transport layer of data is achieved using common network protocols or other techniques (e.g. shared memory, custom bus, ...).
Data format is a key feature in this scenario and each actor must be able to manage event structures in proper way in order to carry out processing or storage. This tutorial will show how to transmit and receive ROOT TTree data objects between Narval C++ actors. The Narval topology of the example is a simple Producer → Filter → Consumer pattern that is very common in DAQ design. GRU will be also used in Filter actor to visually show transmitted data.
Producer code¶
The “core” of each Narval C++ actor is the process_block method that has the following signature:
void pbuffer::process_block (void *output_buffer,
unsigned int size_of_output_buffer,
unsigned int *used_size_of_output_buffer,
unsigned int *error_code)
output_buffer is a pointer to pre-allocated memory area - size definition of output_buffer area is contained in topology definition file and available in C++ code trough size_of_output_buffer parameter.
The Producer code fills a TTree branch with random values then serialise the tree for network transmission.
TTree *t1 = new TTree();
TBufferFile *tbuf = new TBufferFile(TBuffer::kWrite,32000);
Double_t random;
char *source;
cout << "PROD::process_block : libproduce.cc start\n";
t1->Branch("random",&random,"random/D");
for (Int_t i=0;i<100;i++) {
random = gRandom->Rndm();
t1->Fill();
}
t1->Streamer(*tbuf);
source = tbuf->Buffer();
*used_size_of_output_buffer = tbuf->BufferSize();
memcpy((char *) output_buffer, source, tbuf->BufferSize());
cout <<"PROD::process_block : libproduce.cc done\n";
some_algo_dependent_data++;
*error_code = 0;
tbuf->Delete();
t1->Delete();
At lines 1 and 2 a new TTree and TBufferFile are initialized. tbuf object is needed for data serialisation and has a size of 32000 bytes.
Starting from line 8 the TTree object is populated with a branch called random and with 100 random values.
Line 14 invokes Streamer method of TTree object and, as result, it is serialized in tbuf.
At the end (lines from 16 to 19) output_buffer is updated with serialised data (array of char) with memcpy function and lines 26 - 27 provide deallocation of temporary objects.
Filter code¶
The process_block of Narval C++ Filter has the following signature:
void fbuffer::process_block (void *input_buffer,
unsigned int size_of_input_buffer,
void *output_buffer,
unsigned int size_of_output_buffer,
unsigned int *used_size_of_output_buffer,
unsigned int *error_code)
input_buffer and output_buffer are pointers to pre-allocated memory area - size definition of input_buffer and output_buffer area are contained in topology definition file and available in C++ code trough size_of_input_buffer and size_of_output_buffer parameters.
The Filter receives TTree serialised objects in its input_buffer parameter. As next step object deserialisation is performed and all TTree random values are used to update a TH1I ROOT histogram available on network socket thanks to GRU.
GRU socket histogram is initialised in main object constructor:
fbuffer::fbuffer ()
{
Int_t port = 9090;
MyDB = new GSpectra();
serv = new GNetServerRoot(MyDB);
serv->SetPort(port);
serv->StartGNetServer();
spe = new TH1I ("MySpectra","MySpectra",1000,0,1000);
MyDB->AddSpectrum(spe,"SpectraFamily");
}
The filtering routine is the process_block method:
TTree *t2 = new TTree();
TBranch *bpoint;
TBufferFile *tbuf2 = new TBufferFile(TBuffer::kWrite, 32000);
Double_t value;
cout << "FILT::process_block : libfilter.cc start\n";
tbuf22->Reset();
tbuf2->SetWriteMode();
tbuf2->SetBuffer((char *) input_buffer, size_of_input_buffer);
tbuf2->SetReadMode();
t2->Streamer(*tbuf2);
bpoint = t2->GetBranch("random");
bpoint->SetAddress(&value);
Int_t nentries = (Int_t) t2->GetEntries();
for (Int_t i=0;i<nentries;i++) {
t2->GetEntry(i);
bpoint->GetEntry(i);
spe->Fill(value*1000);
}
memcpy((char *) output_buffer, source, tbuf->BufferSize());
*used_size_of_output_buffer = size_of_input_buffer;
*error_code = 0;
tbuf2->Delete();
t2->Delete();
cout << "FILT::process_block : libfilter.cc done\n";
At lines 1 – 3 a new TTree, TBufferFile and TBranch are initialized. tbuf object is needed for data serialization and has a size of 32000 bytes, bpoint is an object needed to point to TTree branch data.
On lines 8 – 10 the tbuf object is initialized with input_buffer (serialized data from Producer) and on line 13 the local TTree object is initialized with serialised data. Lines from 15 to 22 browse TTree object and fill ROOT TH1I spe histogram.
At the end (lines from 25 to 27) output_buffer is updated with original serialised data (array of char) with memcpy function and lines 29 - 30 provide deallocation of temporary objects.
Consumer actor is not shown but it is possible to transfer TTree object following the above method. With VIGRU it is possible to visualize in real-time histogram content and apply fitting functions.
Complete source code for Producer, Filter and Consumer of this tutorial is available below. Topology file and start/stop scripts are also included.
IMPORTANT: The order to load ROOT libraries in Narval can change depending on your ROOT installation.