Belle II Software development
TRGECLUnpackerModule Class Reference

A module of TRG ECL Unpacker. More...

#include <trgeclUnpackerModule.h>

Inheritance diagram for TRGECLUnpackerModule:
Module PathElement

Public Types

enum  EModulePropFlags {
  c_Input = 1 ,
  c_Output = 2 ,
  c_ParallelProcessingCertified = 4 ,
  c_HistogramManager = 8 ,
  c_InternalSerializer = 16 ,
  c_TerminateInAllProcesses = 32 ,
  c_DontCollectStatistics = 64
}
 Each module can be tagged with property flags, which indicate certain features of the module. More...
 
typedef ModuleCondition::EAfterConditionPath EAfterConditionPath
 Forward the EAfterConditionPath definition from the ModuleCondition.
 

Public Member Functions

 TRGECLUnpackerModule ()
 Constructor.
 
virtual ~TRGECLUnpackerModule ()
 Destructor.
 
void initialize () override
 Initilizes TRGECLUnpackerModuel.
 
void event () override
 Called event by event.
 
void terminate () override
 Called when processing ended.
 
void beginRun () override
 Called when new run started.
 
void endRun () override
 Called when run ended.
 
std::string version () const
 returns version of TRGECLUnpackerModule.
 
virtual void readCOPPEREvent (RawTRG *, int, int, int)
 Read data from TRG copper.
 
virtual void checkBuffer (int *, int)
 Unpacker main function.
 
virtual void checkBuffer_v136 (int *, int)
 Unpacker main function for upto version 136.
 
virtual std::vector< std::string > getFileNames (bool outputFiles)
 Return a list of output filenames for this modules.
 
const std::string & getName () const
 Returns the name of the module.
 
const std::string & getType () const
 Returns the type of the module (i.e.
 
const std::string & getPackage () const
 Returns the package this module is in.
 
const std::string & getDescription () const
 Returns the description of the module.
 
void setName (const std::string &name)
 Set the name of the module.
 
void setPropertyFlags (unsigned int propertyFlags)
 Sets the flags for the module properties.
 
LogConfiggetLogConfig ()
 Returns the log system configuration.
 
void setLogConfig (const LogConfig &logConfig)
 Set the log system configuration.
 
void setLogLevel (int logLevel)
 Configure the log level.
 
void setDebugLevel (int debugLevel)
 Configure the debug messaging level.
 
void setAbortLevel (int abortLevel)
 Configure the abort log level.
 
void setLogInfo (int logLevel, unsigned int logInfo)
 Configure the printed log information for the given level.
 
void if_value (const std::string &expression, const std::shared_ptr< Path > &path, EAfterConditionPath afterConditionPath=EAfterConditionPath::c_End)
 Add a condition to the module.
 
void if_false (const std::shared_ptr< Path > &path, EAfterConditionPath afterConditionPath=EAfterConditionPath::c_End)
 A simplified version to add a condition to the module.
 
void if_true (const std::shared_ptr< Path > &path, EAfterConditionPath afterConditionPath=EAfterConditionPath::c_End)
 A simplified version to set the condition of the module.
 
bool hasCondition () const
 Returns true if at least one condition was set for the module.
 
const ModuleConditiongetCondition () const
 Return a pointer to the first condition (or nullptr, if none was set)
 
const std::vector< ModuleCondition > & getAllConditions () const
 Return all set conditions for this module.
 
bool evalCondition () const
 If at least one condition was set, it is evaluated and true returned if at least one condition returns true.
 
std::shared_ptr< PathgetConditionPath () const
 Returns the path of the last true condition (if there is at least one, else reaturn a null pointer).
 
Module::EAfterConditionPath getAfterConditionPath () const
 What to do after the conditional path is finished.
 
std::vector< std::shared_ptr< Path > > getAllConditionPaths () const
 Return all condition paths currently set (no matter if the condition is true or not).
 
bool hasProperties (unsigned int propertyFlags) const
 Returns true if all specified property flags are available in this module.
 
bool hasUnsetForcedParams () const
 Returns true and prints error message if the module has unset parameters which the user has to set in the steering file.
 
const ModuleParamListgetParamList () const
 Return module param list.
 
template<typename T >
ModuleParam< T > & getParam (const std::string &name) const
 Returns a reference to a parameter.
 
bool hasReturnValue () const
 Return true if this module has a valid return value set.
 
int getReturnValue () const
 Return the return value set by this module.
 
std::shared_ptr< PathElementclone () const override
 Create an independent copy of this module.
 
std::shared_ptr< boost::python::list > getParamInfoListPython () const
 Returns a python list of all parameters.
 

Static Public Member Functions

static void exposePythonAPI ()
 Exposes methods of the Module class to Python.
 

Protected Member Functions

virtual void def_initialize ()
 Wrappers to make the methods without "def_" prefix callable from Python.
 
virtual void def_beginRun ()
 Wrapper method for the virtual function beginRun() that has the implementation to be used in a call from Python.
 
virtual void def_event ()
 Wrapper method for the virtual function event() that has the implementation to be used in a call from Python.
 
virtual void def_endRun ()
 This method can receive that the current run ends as a call from the Python side.
 
virtual void def_terminate ()
 Wrapper method for the virtual function terminate() that has the implementation to be used in a call from Python.
 
void setDescription (const std::string &description)
 Sets the description of the module.
 
void setType (const std::string &type)
 Set the module type.
 
template<typename T >
void addParam (const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
 Adds a new parameter to the module.
 
template<typename T >
void addParam (const std::string &name, T &paramVariable, const std::string &description)
 Adds a new enforced parameter to the module.
 
void setReturnValue (int value)
 Sets the return value for this module as integer.
 
void setReturnValue (bool value)
 Sets the return value for this module as bool.
 
void setParamList (const ModuleParamList &params)
 Replace existing parameter list.
 

Protected Attributes

int n_basf2evt
 Event number.
 
int etm_version = 0
 ETM Version.
 
unsigned int nodeid = 0
 Node Id.
 
int nwords = 0
 N Word.
 
int iFiness = 0
 Finess.
 
int trgtype = 0
 Trigger Type.
 

Private Member Functions

std::list< ModulePtrgetModules () const override
 no submodules, return empty list
 
std::string getPathString () const override
 return the module name.
 
void setParamPython (const std::string &name, const boost::python::object &pyObj)
 Implements a method for setting boost::python objects.
 
void setParamPythonDict (const boost::python::dict &dictionary)
 Implements a method for reading the parameter values from a boost::python dictionary.
 

Private Attributes

StoreArray< TRGECLUnpackerStorem_TRGECLTCArray
 ECL Trigger Unpacker TC output.
 
StoreArray< TRGECLUnpackerSumStorem_TRGECLSumArray
 ECL Trigger Unpacker Summary output.
 
StoreArray< TRGECLUnpackerEvtStorem_TRGECLEvtArray
 ECL Trigger Unpacker Event output.
 
StoreArray< TRGECLClusterm_TRGECLClusterArray
 ECL Trigger Cluster output.
 
StoreObjPtr< EventLevelClusteringInfom_eventLevelClusteringInfo
 EventLevelClusteringInfo.
 
std::string m_name
 The name of the module, saved as a string (user-modifiable)
 
std::string m_type
 The type of the module, saved as a string.
 
std::string m_package
 Package this module is found in (may be empty).
 
std::string m_description
 The description of the module.
 
unsigned int m_propertyFlags
 The properties of the module as bitwise or (with |) of EModulePropFlags.
 
LogConfig m_logConfig
 The log system configuration of the module.
 
ModuleParamList m_moduleParamList
 List storing and managing all parameter of the module.
 
bool m_hasReturnValue
 True, if the return value is set.
 
int m_returnValue
 The return value.
 
std::vector< ModuleConditionm_conditions
 Module condition, only non-null if set.
 

Detailed Description

A module of TRG ECL Unpacker.

Definition at line 33 of file trgeclUnpackerModule.h.

Member Typedef Documentation

◆ EAfterConditionPath

Forward the EAfterConditionPath definition from the ModuleCondition.

Definition at line 88 of file Module.h.

Member Enumeration Documentation

◆ EModulePropFlags

enum EModulePropFlags
inherited

Each module can be tagged with property flags, which indicate certain features of the module.

Enumerator
c_Input 

This module is an input module (reads data).

c_Output 

This module is an output module (writes data).

c_ParallelProcessingCertified 

This module can be run in parallel processing mode safely (All I/O must be done through the data store, in particular, the module must not write any files.)

c_HistogramManager 

This module is used to manage histograms accumulated by other modules.

c_InternalSerializer 

This module is an internal serializer/deserializer for parallel processing.

c_TerminateInAllProcesses 

When using parallel processing, call this module's terminate() function in all processes().

This will also ensure that there is exactly one process (single-core if no parallel modules found) or at least one input, one main and one output process.

c_DontCollectStatistics 

No statistics is collected for this module.

Definition at line 77 of file Module.h.

77 {
78 c_Input = 1,
79 c_Output = 2,
85 };
@ c_HistogramManager
This module is used to manage histograms accumulated by other modules.
Definition: Module.h:81
@ c_Input
This module is an input module (reads data).
Definition: Module.h:78
@ c_DontCollectStatistics
No statistics is collected for this module.
Definition: Module.h:84
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
@ c_InternalSerializer
This module is an internal serializer/deserializer for parallel processing.
Definition: Module.h:82
@ c_Output
This module is an output module (writes data).
Definition: Module.h:79
@ c_TerminateInAllProcesses
When using parallel processing, call this module's terminate() function in all processes().
Definition: Module.h:83

Constructor & Destructor Documentation

◆ TRGECLUnpackerModule()

Constructor.

Definition at line 25 of file trgeclUnpackerModule.cc.

27{
28
29 string desc = "TRGECLUnpackerModule(" + version() + ")";
30 setDescription(desc);
32
33 B2DEBUG(20, "trgeclunpacker: Constructor done.");
34}
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
Module()
Constructor.
Definition: Module.cc:30
std::string version() const
returns version of TRGECLUnpackerModule.

◆ ~TRGECLUnpackerModule()

~TRGECLUnpackerModule ( )
virtual

Destructor.

Definition at line 36 of file trgeclUnpackerModule.cc.

36{}

Member Function Documentation

◆ beginRun()

void beginRun ( void  )
overridevirtual

Called when new run started.

Reimplemented from Module.

Definition at line 54 of file trgeclUnpackerModule.cc.

54{}

◆ checkBuffer()

void checkBuffer ( int *  rdat,
int  nnn 
)
virtual

Unpacker main function.

Definition at line 122 of file trgeclUnpackerModule.cc.

123{
124
125 int version_check = (rdat[0] >> 12) & 0xf;
126 if (version_check != 15) return;
127
128 // Checksum variable
129 unsigned char check_sum = (rdat[nnn - 1] >> 24) & 0xFF;
130 unsigned char data_sum = 0;
131 unsigned char kdat[4] = {0};
132 for (int j = nnn - 2; j > -1; j--) {
133 kdat[0] = rdat[j] & 0xff;
134 kdat[1] = (rdat[j] >> 8) & 0xff;
135 kdat[2] = (rdat[j] >> 16) & 0xff;
136 kdat[3] = (rdat[j] >> 24) & 0xff;
137 for (int k = 0; k < 4; k++) {
138 data_sum = data_sum + kdat[k];
139 }
140 }
141
142 int flag_checksum = 0;
143
144 if (check_sum == data_sum) {
145 flag_checksum = 0;
146 } else {
147 flag_checksum = 1;
148 }
149
150 // Information
151 int l1_revo = rdat[0] & 0x7ff;
152 int i = 0;
153 int window_num = 0;
154
155 // Summary
156 /* cppcheck-suppress variableScope */
157 int summary_data = 0;
158 int summary_recon = 0;
159 /* cppcheck-suppress variableScope */
160 int summary_revo = 0;
161 /* cppcheck-suppress variableScope */
162 bool summary_trg = false;
163 /* cppcheck-suppress variableScope */
164 int data_win = 0;
165
166 // TC
167 /* cppcheck-suppress variableScope */
168 int ntc_win = 0;
169 bool tc_trg = false;
170 // TC info
171 int tc_id = 0;
172 int tc_t = 0;
173 int tc_e = 0;
174 int conv_tc_t = 0;
175 int win3_revo = -9999;
176
177 vector<unsigned> sum_data;
178 vector<vector<unsigned>> sum_info;
179
180 vector<int> tc_data;
181 vector<vector<int>> tc_info;
182 vector<vector<int>> tc_info_FE1;
183 vector<vector<int>> tc_info_FE2;
184 vector<vector<int>> tc_info_BE1;
185 vector<vector<int>> tc_info_BE2;
186
187 // Unpacking ---->
188 while (i < nnn - 2) {
189 summary_data = rdat[i + 1];
190 summary_trg = (summary_data >> 23) & 0x1;
191 summary_revo = (summary_data >> 16) & 0x7f;
192 ntc_win = summary_data & 0x3ff;
193 if (ntc_win == 0) {
194 tc_trg = false;
195 } else {
196 tc_trg = true;
197 }
198 data_win = window_num;
199 if (window_num == 3) win3_revo = summary_revo;
200
201 if (summary_trg == true) { // Summary on
202 sum_data.push_back(data_win);
203 sum_data.push_back(summary_revo);
204 for (int j = 0; j < 12; j++) {
205 summary_recon =
206 (((rdat[i + j + 2] >> 0) & 0xFF) << 24) +
207 (((rdat[i + j + 2] >> 8) & 0xFF) << 16) +
208 (((rdat[i + j + 2] >> 16) & 0xFF) << 8) +
209 (((rdat[i + j + 2] >> 24) & 0xFF) << 0);
210 sum_data.push_back(summary_recon);
211 }
212 sum_info.push_back(sum_data);
213 sum_data.clear();
214 i = i + 14;
215
216 if (tc_trg == true) { // summary on & TC on
217 for (int j = 0; j < ntc_win; j++) {
218 tc_id = (rdat[i + j] >> 20) & 0x3FF;
219 tc_t = (rdat[i + j] >> 12) & 0x7F;
220 tc_e = rdat[i + j] & 0xFFF;
221 conv_tc_t = (data_win - 3) * 128 + tc_t;
222
223 // TC vector
224 tc_data.push_back(tc_id);
225 tc_data.push_back(conv_tc_t);
226 tc_data.push_back(tc_e);
227 tc_data.push_back(data_win);
228 if (tc_id < 81) {
229 if (tc_id > 75) {
230 tc_info_FE1.push_back(tc_data);
231 } else {
232 tc_info_FE2.push_back(tc_data);
233 }
234 } else if (tc_id > 512) {
235 if (tc_id > 572) {
236 tc_info_BE1.push_back(tc_data);
237 } else {
238 tc_info_BE2.push_back(tc_data);
239 }
240 } else {
241 tc_info.push_back(tc_data);
242 }
243 tc_data.clear();
244 }
245 i = i + ntc_win - 1;
246 }
247
248 } else { // Summary off
249 if (tc_trg == true) { // summary off & TC on
250 for (int j = 0; j < ntc_win; j++) {
251 tc_id = (rdat[i + j + 2] >> 20) & 0x3FF;
252 tc_t = (rdat[i + j + 2] >> 12) & 0x7F;
253 conv_tc_t = (data_win - 3) * 128 + tc_t;
254 tc_e = rdat[i + j + 2] & 0xFFF;
255
256 // TC vector
257 tc_data.push_back(tc_id);
258 tc_data.push_back(conv_tc_t);
259 tc_data.push_back(tc_e);
260 tc_data.push_back(data_win);
261 if (tc_id < 81) {
262 if (tc_id > 75) {
263 tc_info_FE1.push_back(tc_data);
264 } else {
265 tc_info_FE2.push_back(tc_data);
266 }
267 } else if (tc_id > 512) {
268 if (tc_id > 572) {
269 tc_info_BE1.push_back(tc_data);
270 } else {
271 tc_info_BE2.push_back(tc_data);
272 }
273 } else {
274 tc_info.push_back(tc_data);
275 }
276 tc_data.clear();
277 }
278 i = i + ntc_win + 1;
279 } else { // Summary off & TC off
280 i = i + 1;
281 }
282 }
283 window_num++;
284 }
285
286 // <---- Unpacking
287
288 // Summary
289 /* cppcheck-suppress variableScope */
290 int sum_num_ord = 0;
291 /* cppcheck-suppress variableScope */
292 int sum_num = 0;
293 /* cppcheck-suppress variableScope */
294 int sum_revo = 0;
295 int cl_theta[6] = {0};
296 int cl_phi[6] = {0};
297 int cl_time[6] = { -9999};
298 int cl_energy[6] = {0};
299 int cl_1gev[6] = {0};
300 int cl_2gev[6] = {0};
301 int cl_bha[6] = {0};
302 int ncl = 0;
303 int low_multi = 0;
304 int b2bhabha_v = 0;
305 int b2bhabha_s = 0;
306 int mumu = 0;
307 int prescale = 0;
308 int icn_over = 0;
309 int bg_veto = 0;
310 int icn = 0;
311 int etot_type = 0;
312 int etot = 0;
313 int b1_type = 0;
314 int b1bhabha = 0;
315 int physics = 0;
316 int time_type = 0;
317 int time = 0;
318 int ecl_bst = 0;
319
320 int m_sumNum = 0;
321
322 TrgEclDataBase database;
323 TrgEclMapping mapping;
324
325 vector<int> cl_1d;
326 vector<vector<int>> cl_2d;
327
328 vector<int> evt_1d_vector;
329 vector<vector<int>> evt_2d_vector;
330
331 // Store Summary
332 int sum_size = sum_info.size();
333 if (sum_size != 0) {
334 for (int j = 0; j < sum_size; j++) {
335 sum_num = sum_info[j][0];
336 sum_revo = sum_info[j][1];
337 // TRG
338 time = (sum_info[j][2]) & 0x7F;
339 time_type = (sum_info[j][2] >> 7) & 0x7;
340 physics = (sum_info[j][2] >> 10) & 0x1;
341 b1bhabha = (sum_info[j][2] >> 11) & 0x1;
342 b1_type = (sum_info[j][2] >> 12) & 0x3FFF;
343 etot = ((sum_info[j][3] & 0x7F) << 6) + ((sum_info[j][2] >> 26) & 0x3F);
344 etot_type = (sum_info[j][3] >> 7) & 0x7;
345 icn = (sum_info[j][3] >> 10) & 0x7F;
346 bg_veto = (sum_info[j][3] >> 17) & 0x7;
347 icn_over = (sum_info[j][3] >> 20) & 0x1;
348 b2bhabha_v = (sum_info[j][3] >> 21) & 0x1;
349 low_multi = (((sum_info[j][4] >> 6) & 0x3) << 12) + ((sum_info[j][4] & 0x3) << 10) + ((sum_info[j][3] >> 22) & 0x3FF);
350 b2bhabha_s = (sum_info[j][4] >> 2) & 0x1;
351 mumu = (sum_info[j][4] >> 3) & 0x1;
352 prescale = (sum_info[j][4] >> 4) & 0x1;
353 ecl_bst = (sum_info[j][4] >> 5) & 0x1;
354 // CL
355 cl_energy[0] = (sum_info[j][5]) & 0xFFF;
356 cl_time[0] = (sum_info[j][5] >> 12) & 0xFF;
357 cl_phi[0] = (sum_info[j][5] >> 20) & 0xFF;
358 cl_theta[0] = ((sum_info[j][6] & 0x7) << 4) + ((sum_info[j][5] >> 28) & 0xF);
359
360 cl_energy[1] = (sum_info[j][6] >> 3) & 0xFFF;
361 cl_time[1] = (sum_info[j][6] >> 15) & 0xFF;
362 cl_phi[1] = (sum_info[j][6] >> 23) & 0xFF;
363 cl_theta[1] = ((sum_info[j][7] & 0x3F) << 1) + ((sum_info[j][6] >> 31) & 0x1);
364
365 cl_energy[2] = (sum_info[j][7] >> 6) & 0xFFF;
366 cl_time[2] = (sum_info[j][7] >> 18) & 0xFF;
367 cl_phi[2] = ((sum_info[j][8] & 0x3) << 6) + ((sum_info[j][7] >> 26) & 0x3F);
368 cl_theta[2] = (sum_info[j][8] >> 2) & 0x7F;
369
370 cl_energy[3] = (sum_info[j][8] >> 9) & 0xFFF;
371 cl_time[3] = (sum_info[j][8] >> 21) & 0xFF;
372 cl_phi[3] = ((sum_info[j][9] & 0x1F) << 3) + ((sum_info[j][8] >> 29) & 0x7);
373 cl_theta[3] = (sum_info[j][9] >> 5) & 0x7F;
374
375 cl_energy[4] = (sum_info[j][ 9] >> 12) & 0xFFF;
376 cl_time[4] = (sum_info[j][ 9] >> 24) & 0xFF;
377 cl_phi[4] = (sum_info[j][10]) & 0xFF;
378 cl_theta[4] = (sum_info[j][10] >> 8) & 0x7F;
379
380 cl_energy[5] = (sum_info[j][10] >> 15) & 0xFFF;
381 cl_time[5] = ((sum_info[j][11] & 0x7) << 5) + ((sum_info[j][10] >> 27) & 0x1F);
382 cl_phi[5] = (sum_info[j][11] >> 3) & 0xFF;
383 cl_theta[5] = (sum_info[j][11] >> 11) & 0x7F;
384 // CL others
385 for (int k = 0; k < 6; k++) {
386 cl_1gev[k] = (sum_info[j][12] >> k) & 0x1;
387 cl_2gev[k] = (sum_info[j][12] >> (k + 6)) & 0x1;
388 cl_bha[k] = (sum_info[j][12] >> (k + 12)) & 0x1;
389 }
390 ncl = (sum_info[j][13]) & 0x7;
391
392 m_TRGECLSumArray.appendNew();
393 m_sumNum = m_TRGECLSumArray.getEntries() - 1;
394 m_TRGECLSumArray[m_sumNum]->setEventId(n_basf2evt);
395 m_TRGECLSumArray[m_sumNum]->setSumNum(sum_num);
396 m_TRGECLSumArray[m_sumNum]->setSumRevo(sum_revo);
397 m_TRGECLSumArray[m_sumNum]->setCLTheta(cl_theta);
398 m_TRGECLSumArray[m_sumNum]->setCLPhi(cl_phi);
399 m_TRGECLSumArray[m_sumNum]->setCLTime(cl_time);
400 m_TRGECLSumArray[m_sumNum]->setCLEnergy(cl_energy);
401 m_TRGECLSumArray[m_sumNum]->setCLF1GeV(cl_1gev);
402 m_TRGECLSumArray[m_sumNum]->setCLF2GeV(cl_2gev);
403 m_TRGECLSumArray[m_sumNum]->setCLFBha(cl_bha);
404 m_TRGECLSumArray[m_sumNum]->setNCL(ncl);
405 m_TRGECLSumArray[m_sumNum]->setICN(icn);
406 m_TRGECLSumArray[m_sumNum]->setICNOver(icn_over);
407 m_TRGECLSumArray[m_sumNum]->setLowMulti(low_multi);
408 m_TRGECLSumArray[m_sumNum]->set3DBhabhaV(b2bhabha_v);
409 m_TRGECLSumArray[m_sumNum]->set3DBhabhaS(b2bhabha_s);
410 m_TRGECLSumArray[m_sumNum]->setMumu(mumu);
411 m_TRGECLSumArray[m_sumNum]->setPrescale(prescale);
412 m_TRGECLSumArray[m_sumNum]->set2DBhabha(b1bhabha);
413 m_TRGECLSumArray[m_sumNum]->setBhabhaType(b1_type);
414 m_TRGECLSumArray[m_sumNum]->setPhysics(physics);
415 m_TRGECLSumArray[m_sumNum]->setBG(bg_veto);
416 m_TRGECLSumArray[m_sumNum]->setEtot(etot);
417 m_TRGECLSumArray[m_sumNum]->setEtotType(etot_type);
418 m_TRGECLSumArray[m_sumNum]->setECLBST(ecl_bst);
419 m_TRGECLSumArray[m_sumNum]->setTime(time);
420 m_TRGECLSumArray[m_sumNum]->setTimeType(time_type);
421
422 for (int k = 0; k < 6; k++) {
423 cl_1d.push_back(cl_theta[k]);
424 cl_1d.push_back(cl_phi[k]);
425 cl_1d.push_back(cl_time[k]);
426 cl_1d.push_back(cl_energy[k]);
427 cl_1d.push_back(cl_1gev[k]);
428 cl_1d.push_back(cl_2gev[k]);
429 cl_1d.push_back(cl_bha[k]);
430 cl_2d.push_back(cl_1d);
431 cl_1d.clear();
432 }
433 sort(cl_2d.begin(), cl_2d.end(),
434 [](const vector<int>& aa1, const vector<int>& aa2) {return aa1[3] > aa2[3];});
435
436 if (sum_num == -9999) {
437 sum_num_ord = -9999;
438 } else if (sum_num <= 3) {
439 sum_num_ord = 2 * abs(sum_num - 3);
440 } else {
441 sum_num_ord = (sum_num * 2) - 7;
442 }
443 evt_1d_vector.push_back(sum_num_ord);
444 evt_1d_vector.push_back(sum_revo);
445 evt_1d_vector.push_back(sum_num);
446 evt_1d_vector.push_back(time);
447 for (int k = 0; k < 6; k++) {
448 for (int l = 0; l < 7; l++) {
449 evt_1d_vector.push_back(cl_2d[k][l]);
450 }
451 }
452 evt_1d_vector.push_back(ncl);
453 evt_1d_vector.push_back(low_multi);
454 evt_1d_vector.push_back(b2bhabha_v);
455 evt_1d_vector.push_back(b2bhabha_s);
456 evt_1d_vector.push_back(mumu);
457 evt_1d_vector.push_back(prescale);
458 evt_1d_vector.push_back(icn);
459 evt_1d_vector.push_back(icn_over);
460 evt_1d_vector.push_back(etot_type);
461 evt_1d_vector.push_back(etot);
462 evt_1d_vector.push_back(ecl_bst);
463 evt_1d_vector.push_back(b1_type);
464 evt_1d_vector.push_back(b1bhabha);
465 evt_1d_vector.push_back(physics);
466 evt_1d_vector.push_back(time_type);
467 evt_2d_vector.push_back(evt_1d_vector);
468 evt_1d_vector.clear();
469 }
470 } else {
471 for (int k = 0; k < 6; k++) {
472 cl_theta[k] = 0;
473 cl_phi[k] = 0;
474 cl_time[k] = -9999;
475 cl_energy[k] = 0;
476 cl_1gev[k] = 0;
477 cl_2gev[k] = 0;
478 cl_bha[k] = 0;
479 }
480 ncl = 0;
481 low_multi = 0;
482 b2bhabha_v = 0;
483 b2bhabha_s = 0;
484 mumu = 0;
485 prescale = 0;
486 icn_over = 0;
487 bg_veto = 0;
488 icn = 0;
489 etot_type = 0;
490 etot = 0;
491 ecl_bst = 0;
492 b1_type = 0;
493 b1bhabha = 0;
494 physics = 0;
495 time_type = 0;
496 time = -9999;
497
498 m_TRGECLSumArray.appendNew();
499 m_sumNum = m_TRGECLSumArray.getEntries() - 1;
500 m_TRGECLSumArray[m_sumNum]->setEventId(n_basf2evt);
501 m_TRGECLSumArray[m_sumNum]->setSumNum(0);
502 m_TRGECLSumArray[m_sumNum]->setCLTheta(cl_theta);
503 m_TRGECLSumArray[m_sumNum]->setCLPhi(cl_phi);
504 m_TRGECLSumArray[m_sumNum]->setCLTime(cl_time);
505 m_TRGECLSumArray[m_sumNum]->setCLEnergy(cl_energy);
506 m_TRGECLSumArray[m_sumNum]->setCLF1GeV(cl_1gev);
507 m_TRGECLSumArray[m_sumNum]->setCLF2GeV(cl_2gev);
508 m_TRGECLSumArray[m_sumNum]->setCLFBha(cl_bha);
509 m_TRGECLSumArray[m_sumNum]->setNCL(ncl);
510 m_TRGECLSumArray[m_sumNum]->setICN(icn);
511 m_TRGECLSumArray[m_sumNum]->setICNOver(icn_over);
512 m_TRGECLSumArray[m_sumNum]->setLowMulti(low_multi);
513 m_TRGECLSumArray[m_sumNum]->set3DBhabhaV(b2bhabha_v);
514 m_TRGECLSumArray[m_sumNum]->set3DBhabhaS(b2bhabha_s);
515 m_TRGECLSumArray[m_sumNum]->setMumu(mumu);
516 m_TRGECLSumArray[m_sumNum]->setPrescale(prescale);
517 m_TRGECLSumArray[m_sumNum]->set2DBhabha(b1bhabha);
518 m_TRGECLSumArray[m_sumNum]->setBhabhaType(b1_type);
519 m_TRGECLSumArray[m_sumNum]->setPhysics(physics);
520 m_TRGECLSumArray[m_sumNum]->setBG(bg_veto);
521 m_TRGECLSumArray[m_sumNum]->setEtot(etot);
522 m_TRGECLSumArray[m_sumNum]->setEtotType(etot_type);
523 m_TRGECLSumArray[m_sumNum]->setECLBST(ecl_bst);
524 m_TRGECLSumArray[m_sumNum]->setTime(time);
525 m_TRGECLSumArray[m_sumNum]->setTimeType(time_type);
526 }
527
528 // TC & TRG
529 tc_info.insert(tc_info.end(), tc_info_FE1.begin(), tc_info_FE1.end());
530 tc_info.insert(tc_info.end(), tc_info_FE2.begin(), tc_info_FE2.end());
531 tc_info.insert(tc_info.end(), tc_info_BE1.begin(), tc_info_BE1.end());
532 tc_info.insert(tc_info.end(), tc_info_BE2.begin(), tc_info_BE2.end());
533
534 int m_evtNum = 0;
535
536 int m_tcNum = 0;
537 /* cppcheck-suppress variableScope */
538 int m_tcid = 0;
539 /* cppcheck-suppress variableScope */
540 int m_time = -9999;
541 /* cppcheck-suppress variableScope */
542 int m_energy = 0;
543 /* cppcheck-suppress variableScope */
544 int m_win = 0;
545 /* cppcheck-suppress variableScope */
546 int m_revo = 0;
547 /* cppcheck-suppress variableScope */
548 int m_caltime = -9999;
549
550 int tot_ntc = tc_info.size();
551 /* cppcheck-suppress variableScope */
552 int evt_ntc = 0;
553 /* cppcheck-suppress variableScope */
554 int evt_revo = -9999;
555 /* cppcheck-suppress variableScope */
556 int evt_win = 0;
557 /* cppcheck-suppress variableScope */
558 int evt_timing = -9999; // most energetic
559 int evt_cl_theta[6] = {0};
560 int evt_cl_phi[6] = {0};
561 int evt_cl_time[6] = { -9999};
562 int evt_cl_energy[6] = {0};
563 int evt_cl_1gev[6] = {0};
564 int evt_cl_2gev[6] = {0};
565 int evt_cl_bha[6] = {0};
566 /* cppcheck-suppress variableScope */
567 int evt_ncl = 0;
568 /* cppcheck-suppress variableScope */
569 int evt_low_multi = 0;
570 /* cppcheck-suppress variableScope */
571 int evt_b2bhabha_v = 0;
572 /* cppcheck-suppress variableScope */
573 int evt_b2bhabha_s = 0;
574 /* cppcheck-suppress variableScope */
575 int evt_mumu = 0;
576 /* cppcheck-suppress variableScope */
577 int evt_prescale = 0;
578 /* cppcheck-suppress variableScope */
579 int evt_icn = 0;
580 /* cppcheck-suppress variableScope */
581 int evt_icn_over = 0;
582 /* cppcheck-suppress variableScope */
583 int evt_etot_type = 0;
584 /* cppcheck-suppress variableScope */
585 int evt_etot = 0;
586 /* cppcheck-suppress variableScope */
587 int evt_ecl_bst = 0;
588 /* cppcheck-suppress variableScope */
589 int evt_b1_type = 0;
590 /* cppcheck-suppress variableScope */
591 int evt_b1bhabha = 0;
592 /* cppcheck-suppress variableScope */
593 int evt_physics = 0;
594 /* cppcheck-suppress variableScope */
595 int evt_time_type = 0;
596 /* cppcheck-suppress variableScope */
597 int evt_etot_all = 0;
598 /* cppcheck-suppress variableScope */
599 int evt_time_min = 0;
600 /* cppcheck-suppress variableScope */
601 int evt_time_max = 0;
602 /* cppcheck-suppress variableScope */
603 int evt_time_win = 0;
604 /* cppcheck-suppress variableScope */
605 int etot_i = 0;
606 /* cppcheck-suppress variableScope */
607 int etot_c = 0;
608 /* cppcheck-suppress variableScope */
609 int etot_f = 0;
610 /* cppcheck-suppress variableScope */
611 int cl_tcid = 0;
612 /* cppcheck-suppress variableScope */
613 int cl_thetaid = 0;
614 /* cppcheck-suppress variableScope */
615 int cl_phiid = 0;
616 /* cppcheck-suppress variableScope */
617 int m_clNum = 0;
618
619
620 int evt_v_size = evt_2d_vector.size();
621 if (evt_v_size != 0) {
622 // Sort window : 3 => 4 => 2 => 5 => 1 => 6 => 7
623 sort(evt_2d_vector.begin(), evt_2d_vector.end(),
624 [](const vector<int>& aa1, const vector<int>& aa2) {return aa1[0] < aa2[0];});
625 }
626
627 if (tot_ntc != 0 && flag_checksum == 0 && nnn > 7) {
628 if (evt_v_size == 0) {
629 // Find most energetic TC timing
630 sort(tc_info.begin(), tc_info.end(),
631 [](const vector<int>& aa1, const vector<int>& aa2) {return aa1[2] > aa2[2];});
632 evt_revo = win3_revo;
633 evt_win = tc_info[0][3];
634 evt_timing = tc_info[0][1];
635 for (int k = 0; k < 6; k++) {
636 evt_cl_theta[k] = 0;
637 evt_cl_phi[k] = 0;
638 evt_cl_time[k] = 0;
639 evt_cl_energy[k] = 0;
640 evt_cl_1gev[k] = 0;
641 evt_cl_2gev[k] = 0;
642 evt_cl_bha[k] = 0;
643 }
644 evt_ncl = 0;
645 evt_low_multi = 0;
646 evt_b2bhabha_v = 0;
647 evt_b2bhabha_s = 0;
648 evt_mumu = 0;
649 evt_prescale = 0;
650 evt_icn = 0;
651 evt_icn_over = 0;
652 evt_etot_type = 0;
653 evt_etot = 0;
654 evt_ecl_bst = 0;
655 evt_b1_type = 0;
656 evt_b1bhabha = 0;
657 evt_physics = 0;
658 evt_time_type = 0;
659 } else {
660 evt_revo = evt_2d_vector[0][1];
661 evt_win = evt_2d_vector[0][2];
662 evt_timing = evt_2d_vector[0][3];
663 for (int k = 0; k < 6; k++) {
664 evt_cl_theta[k] = evt_2d_vector[0][ 4 + k * 7];
665 evt_cl_phi[k] = evt_2d_vector[0][ 5 + k * 7];
666 evt_cl_time[k] = evt_2d_vector[0][ 6 + k * 7];
667 evt_cl_energy[k] = evt_2d_vector[0][ 7 + k * 7];
668 evt_cl_1gev[k] = evt_2d_vector[0][ 8 + k * 7];
669 evt_cl_2gev[k] = evt_2d_vector[0][ 9 + k * 7];
670 evt_cl_bha[k] = evt_2d_vector[0][10 + k * 7];
671 }
672 evt_ncl = evt_2d_vector[0][46];
673 evt_low_multi = evt_2d_vector[0][47];
674 evt_b2bhabha_v = evt_2d_vector[0][48];
675 evt_b2bhabha_s = evt_2d_vector[0][49];
676 evt_mumu = evt_2d_vector[0][50];
677 evt_prescale = evt_2d_vector[0][51];
678 evt_icn = evt_2d_vector[0][52];
679 evt_icn_over = evt_2d_vector[0][53];
680 evt_etot_type = evt_2d_vector[0][54];
681 evt_etot = evt_2d_vector[0][55];
682 evt_ecl_bst = evt_2d_vector[0][56];
683 evt_b1_type = evt_2d_vector[0][57];
684 evt_b1bhabha = evt_2d_vector[0][58];
685 evt_physics = evt_2d_vector[0][59];
686 evt_time_type = evt_2d_vector[0][60];
687 }
688 // Sort by TC number
689 sort(tc_info.begin(), tc_info.end(),
690 [](const vector<int>& aa1, const vector<int>& aa2) {return aa1[0] < aa2[0];});
691
692 for (int ii = 0; ii < tot_ntc; ii++) {
693 m_tcid = tc_info[ii][0];
694 m_time = tc_info[ii][1];
695 m_energy = tc_info[ii][2];
696 m_win = tc_info[ii][3];
697 m_revo = win3_revo;
698 m_caltime = m_time - ((evt_win - 3) * 128 + evt_timing);
699 m_TRGECLTCArray.appendNew();
700 m_tcNum = m_TRGECLTCArray.getEntries() - 1;
701 m_TRGECLTCArray[m_tcNum]->setEventId(n_basf2evt);
702 m_TRGECLTCArray[m_tcNum]->setTCId(m_tcid);
703 m_TRGECLTCArray[m_tcNum]->setTCTime(m_time);
704 m_TRGECLTCArray[m_tcNum]->setTCCALTime(m_caltime);
705 m_TRGECLTCArray[m_tcNum]->setHitWin(m_win);
706 m_TRGECLTCArray[m_tcNum]->setRevoFAM(m_revo);
707 m_TRGECLTCArray[m_tcNum]->setTCEnergy(m_energy);
708 m_TRGECLTCArray[m_tcNum]->setChecksum(flag_checksum);
709
710 if (m_win == evt_win || m_win == evt_win + 1) evt_ntc++;
711 if (m_win == evt_win - 1) {
712 etot_i += m_energy;
713 }
714 if (m_win == evt_win) {
715 etot_c += m_energy;
716 }
717 if (m_win == evt_win + 1) {
718 etot_f += m_energy;
719 }
720 }
721
722 if (etot_i == 0 && etot_f == 0) {
723 evt_etot_all = etot_c;
724 evt_time_min = - evt_timing;
725 evt_time_max = 256 - evt_timing;
726 evt_time_win = 1;
727 } else if (etot_i >= etot_f) {
728 evt_etot_all = etot_c + etot_i;
729 evt_time_min = -128 - evt_timing;
730 evt_time_max = 128 - evt_timing;
731 evt_time_win = -1;
732 } else {
733 evt_etot_all = etot_c + etot_f;
734 evt_time_min = - evt_timing;
735 evt_time_max = 256 - evt_timing;
736 evt_time_win = 1;
737 }
738
739 for (int icluster = 0; icluster < 6; icluster++) {
740 if (evt_cl_energy[icluster] == 0 || evt_cl_theta[icluster] == 0 || evt_cl_phi[icluster] == 0) {continue;}
741 cl_tcid = mapping.getTCIdFromPosition(evt_cl_theta[icluster], evt_cl_phi[icluster]);
742 if (cl_tcid == 0) {continue;}
743 cl_thetaid = mapping.getTCThetaIdFromTCId(cl_tcid);
744 cl_phiid = mapping.getTCPhiIdFromTCId(cl_tcid);
745
746 m_TRGECLClusterArray.appendNew();
747 m_clNum = m_TRGECLClusterArray.getEntries() - 1;
748 m_TRGECLClusterArray[m_clNum]->setEventId(n_basf2evt);
749 m_TRGECLClusterArray[m_clNum]->setClusterId(icluster);
750 m_TRGECLClusterArray[m_clNum]->setEventRevo(evt_revo);
751
752 m_TRGECLClusterArray[m_clNum]->setMaxTCId(cl_tcid); // center of Cluster
753 m_TRGECLClusterArray[m_clNum]->setMaxThetaId(cl_thetaid);
754 m_TRGECLClusterArray[m_clNum]->setMaxPhiId(cl_phiid);
755 m_TRGECLClusterArray[m_clNum]->setClusterId(icluster);
756 m_TRGECLClusterArray[m_clNum]->setEnergyDep((double)evt_cl_energy[icluster] * 5.25); // MeV
757 m_TRGECLClusterArray[m_clNum]->setTimeAve((double)evt_cl_time[icluster]);
758 m_TRGECLClusterArray[m_clNum]->setPositionX(mapping.getTCPosition(cl_tcid).X());
759 m_TRGECLClusterArray[m_clNum]->setPositionY(mapping.getTCPosition(cl_tcid).Y());
760 m_TRGECLClusterArray[m_clNum]->setPositionZ(mapping.getTCPosition(cl_tcid).Z());
761 }
762 m_TRGECLEvtArray.appendNew();
763 m_evtNum = m_TRGECLEvtArray.getEntries() - 1;
764 m_TRGECLEvtArray[m_evtNum]->setEventId(n_basf2evt);
765 m_TRGECLEvtArray[m_evtNum]->setETM(etm_version);
766 m_TRGECLEvtArray[m_evtNum]->setL1Revo(l1_revo);
767 m_TRGECLEvtArray[m_evtNum]->setEvtRevo(evt_revo);
768 m_TRGECLEvtArray[m_evtNum]->setEvtWin(evt_win);
769 m_TRGECLEvtArray[m_evtNum]->setEvtTime(evt_timing);
770 m_TRGECLEvtArray[m_evtNum]->setNTC(evt_ntc);
771 m_TRGECLEvtArray[m_evtNum]->setCLTheta(evt_cl_theta);
772 m_TRGECLEvtArray[m_evtNum]->setCLPhi(evt_cl_phi);
773 m_TRGECLEvtArray[m_evtNum]->setCLTime(evt_cl_time);
774 m_TRGECLEvtArray[m_evtNum]->setCLEnergy(evt_cl_energy);
775 m_TRGECLEvtArray[m_evtNum]->setCLF1GeV(evt_cl_1gev);
776 m_TRGECLEvtArray[m_evtNum]->setCLF2GeV(evt_cl_2gev);
777 m_TRGECLEvtArray[m_evtNum]->setCLFBha(evt_cl_bha);
778 m_TRGECLEvtArray[m_evtNum]->setNCL(evt_ncl);
779 m_TRGECLEvtArray[m_evtNum]->setLowMulti(evt_low_multi);
780 m_TRGECLEvtArray[m_evtNum]->set3DBhabhaV(evt_b2bhabha_v);
781 m_TRGECLEvtArray[m_evtNum]->set3DBhabhaS(evt_b2bhabha_s);
782 m_TRGECLEvtArray[m_evtNum]->setMumu(evt_mumu);
783 m_TRGECLEvtArray[m_evtNum]->setPrescale(evt_prescale);
784 m_TRGECLEvtArray[m_evtNum]->setICN(evt_icn);
785 m_TRGECLEvtArray[m_evtNum]->setICNOver(evt_icn_over);
786 m_TRGECLEvtArray[m_evtNum]->setEtotType(evt_etot_type);
787 m_TRGECLEvtArray[m_evtNum]->setEtot(evt_etot);
788 m_TRGECLEvtArray[m_evtNum]->setECLBST(evt_ecl_bst);
789 m_TRGECLEvtArray[m_evtNum]->set2DBhabha(evt_b1bhabha);
790 m_TRGECLEvtArray[m_evtNum]->setBhabhaType(evt_b1_type);
791 m_TRGECLEvtArray[m_evtNum]->setPhysics(evt_physics);
792 m_TRGECLEvtArray[m_evtNum]->setTimeType(evt_time_type);
793 m_TRGECLEvtArray[m_evtNum]->setCheckSum(flag_checksum);
794 m_TRGECLEvtArray[m_evtNum]->setEvtExist(1);
795 m_TRGECLEvtArray[m_evtNum]->setTRGTYPE(trgtype);
796 m_TRGECLEvtArray[m_evtNum]->setEtotAll(evt_etot_all);
797 m_TRGECLEvtArray[m_evtNum]->setEvtTimeMin(evt_time_min);
798 m_TRGECLEvtArray[m_evtNum]->setEvtTimeMax(evt_time_max);
799 m_TRGECLEvtArray[m_evtNum]->setEvtTimeWin(evt_time_win);
800 } else {
801 m_TRGECLTCArray.appendNew();
802 m_tcNum = m_TRGECLTCArray.getEntries() - 1;
803 m_TRGECLTCArray[m_tcNum]->setEventId(n_basf2evt);
804 m_TRGECLTCArray[m_tcNum]->setTCId(0);
805 m_TRGECLTCArray[m_tcNum]->setTCTime(-9999);
806 m_TRGECLTCArray[m_tcNum]->setTCCALTime(-9999);
807 m_TRGECLTCArray[m_tcNum]->setHitWin(-9999);
808 m_TRGECLTCArray[m_tcNum]->setRevoFAM(-9999);
809 m_TRGECLTCArray[m_tcNum]->setTCEnergy(0);
810 m_TRGECLTCArray[m_tcNum]->setChecksum(flag_checksum);
811
812 m_TRGECLEvtArray.appendNew();
813 m_evtNum = m_TRGECLEvtArray.getEntries() - 1;
814 m_TRGECLEvtArray[m_evtNum]->setEventId(n_basf2evt);
815 m_TRGECLEvtArray[m_evtNum]->setETM(etm_version);
816 m_TRGECLEvtArray[m_evtNum]->setL1Revo(-9999);
817 m_TRGECLEvtArray[m_evtNum]->setEvtTime(-9999);
818 m_TRGECLEvtArray[m_evtNum]->setEvtRevo(-9999);
819 m_TRGECLEvtArray[m_evtNum]->setEvtWin(-9999);
820 m_TRGECLEvtArray[m_evtNum]->setNTC(0);
821 for (int k = 0; k < 6; k++) {
822 evt_cl_theta[k] = 0;
823 evt_cl_phi[k] = 0;
824 evt_cl_time[k] = -9999;
825 evt_cl_energy[k] = 0;
826 evt_cl_1gev[k] = 0;
827 evt_cl_2gev[k] = 0;
828 evt_cl_bha[k] = 0;
829 }
830 m_TRGECLEvtArray[m_evtNum]->setCLTheta(evt_cl_theta);
831 m_TRGECLEvtArray[m_evtNum]->setCLPhi(evt_cl_phi);
832 m_TRGECLEvtArray[m_evtNum]->setCLTime(evt_cl_time);
833 m_TRGECLEvtArray[m_evtNum]->setCLEnergy(evt_cl_energy);
834 m_TRGECLEvtArray[m_evtNum]->setCLF1GeV(evt_cl_1gev);
835 m_TRGECLEvtArray[m_evtNum]->setCLF2GeV(evt_cl_2gev);
836 m_TRGECLEvtArray[m_evtNum]->setCLFBha(evt_cl_bha);
837 m_TRGECLEvtArray[m_evtNum]->setNCL(0);
838 m_TRGECLEvtArray[m_evtNum]->setLowMulti(0);
839 m_TRGECLEvtArray[m_evtNum]->set3DBhabhaV(0);
840 m_TRGECLEvtArray[m_evtNum]->set3DBhabhaS(0);
841 m_TRGECLEvtArray[m_evtNum]->setMumu(0);
842 m_TRGECLEvtArray[m_evtNum]->setPrescale(0);
843 m_TRGECLEvtArray[m_evtNum]->setICN(0);
844 m_TRGECLEvtArray[m_evtNum]->setICNOver(0);
845 m_TRGECLEvtArray[m_evtNum]->setEtotType(0);
846 m_TRGECLEvtArray[m_evtNum]->setEtot(0);
847 m_TRGECLEvtArray[m_evtNum]->setECLBST(0);
848 m_TRGECLEvtArray[m_evtNum]->set2DBhabha(0);
849 m_TRGECLEvtArray[m_evtNum]->setBhabhaType(0);
850 m_TRGECLEvtArray[m_evtNum]->setPhysics(0);
851 m_TRGECLEvtArray[m_evtNum]->setTimeType(0);
852 m_TRGECLEvtArray[m_evtNum]->setCheckSum(flag_checksum);
853 m_TRGECLEvtArray[m_evtNum]->setEvtExist(0);
854 m_TRGECLEvtArray[m_evtNum]->setTRGTYPE(trgtype);
855 m_TRGECLEvtArray[m_evtNum]->setEtotAll(0);
856 m_TRGECLEvtArray[m_evtNum]->setEvtTimeMin(-9999);
857 m_TRGECLEvtArray[m_evtNum]->setEvtTimeMax(-9999);
858 m_TRGECLEvtArray[m_evtNum]->setEvtTimeWin(0);
859 }
860
861 return;
862}
StoreArray< TRGECLUnpackerSumStore > m_TRGECLSumArray
ECL Trigger Unpacker Summary output.
StoreArray< TRGECLUnpackerStore > m_TRGECLTCArray
ECL Trigger Unpacker TC output.
StoreArray< TRGECLUnpackerEvtStore > m_TRGECLEvtArray
ECL Trigger Unpacker Event output.
StoreArray< TRGECLCluster > m_TRGECLClusterArray
ECL Trigger Cluster output.
class TrgEclDataBase;
A class of TC Mapping.
Definition: TrgEclMapping.h:26
int getTCThetaIdFromTCId(int)
get [TC Theta ID] from [TC ID]
int getTCIdFromPosition(int, int)
get TCId from phi and theta position(LSB = 1.4)
ROOT::Math::XYZVector getTCPosition(int)
TC position (cm)
int getTCPhiIdFromTCId(int)
get [TC Phi ID] from [TC ID]

◆ checkBuffer_v136()

void checkBuffer_v136 ( int *  rdat,
int  nnn 
)
virtual

Unpacker main function for upto version 136.

Definition at line 864 of file trgeclUnpackerModule.cc.

865{
866
867 int version_check = (rdat[0] >> 12) & 0xf;
868 if (version_check != 15) return;
869
870 // Checksum variable
871 unsigned char check_sum = (rdat[nnn - 1] >> 24) & 0xFF;
872 unsigned char data_sum = 0;
873 unsigned char kdat[4] = {0};
874 for (int j = nnn - 2; j > -1; j--) {
875 kdat[0] = rdat[j] & 0xff;
876 kdat[1] = (rdat[j] >> 8) & 0xff;
877 kdat[2] = (rdat[j] >> 16) & 0xff;
878 kdat[3] = (rdat[j] >> 24) & 0xff;
879 for (int k = 0; k < 4; k++) {
880 data_sum = data_sum + kdat[k];
881 }
882 }
883
884 int flag_checksum = 0;
885
886 if (check_sum == data_sum) {
887 flag_checksum = 0;
888 } else {
889 flag_checksum = 1;
890 }
891
892 // Information
893 int l1_revo = rdat[0] & 0x7ff;
894 int i = 0;
895 int window_num = 0;
896
897 // Summary
898 /* cppcheck-suppress variableScope */
899 int summary_data = 0;
900 /* cppcheck-suppress variableScope */
901 int summary_revo = 0;
902 /* cppcheck-suppress variableScope */
903 bool summary_trg = false;
904 /* cppcheck-suppress variableScope */
905 int data_win = 0;
906
907 // TC
908 /* cppcheck-suppress variableScope */
909 int ntc_win = 0;
910 bool tc_trg = false;
911 // TC info
912 int tc_id = 0;
913 int tc_t = 0;
914 int tc_e = 0;
915 int conv_tc_t = 0;
916 int win3_revo = -9999;
917
918 vector<unsigned> sum_data;
919 vector<vector<unsigned>> sum_info; //TODO can these be unsigned? (required for bit shifts shifts)
920
921 vector<int> tc_data;
922 vector<vector<int>> tc_info;
923 vector<vector<int>> tc_info_FE1;
924 vector<vector<int>> tc_info_FE2;
925 vector<vector<int>> tc_info_BE1;
926 vector<vector<int>> tc_info_BE2;
927
928 // Unpacking ---->
929 while (i < nnn - 2) {
930 summary_data = rdat[i + 1];
931 summary_trg = (summary_data >> 23) & 0x1;
932 summary_revo = (summary_data >> 16) & 0x7f;
933 ntc_win = summary_data & 0x3ff;
934 if (ntc_win == 0) {
935 tc_trg = false;
936 } else {
937 tc_trg = true;
938 }
939 data_win = window_num;
940 if (window_num == 3) win3_revo = summary_revo;
941
942 if (summary_trg == true) { // Summary on
943 sum_data.push_back(data_win);
944 sum_data.push_back(summary_revo);
945 for (int j = 0; j < 9; j++) {
946 sum_data.push_back(rdat[i + j + 2]);
947 }
948 sum_info.push_back(sum_data);
949 sum_data.clear();
950 i = i + 11;
951
952 if (tc_trg == true) { // summary on & TC on
953 for (int j = 0; j < ntc_win; j++) {
954 tc_id = (rdat[i + j] >> 20) & 0x3FF;
955 tc_t = (rdat[i + j] >> 12) & 0x7F;
956 tc_e = rdat[i + j] & 0xFFF;
957 conv_tc_t = (data_win - 3) * 128 + tc_t;
958
959 // TC vector
960 tc_data.push_back(tc_id);
961 tc_data.push_back(conv_tc_t);
962 tc_data.push_back(tc_e);
963 tc_data.push_back(data_win);
964 if (tc_id < 81) {
965 if (tc_id > 75) {
966 tc_info_FE1.push_back(tc_data);
967 } else {
968 tc_info_FE2.push_back(tc_data);
969 }
970 } else if (tc_id > 512) {
971 if (tc_id > 572) {
972 tc_info_BE1.push_back(tc_data);
973 } else {
974 tc_info_BE2.push_back(tc_data);
975 }
976 } else {
977 tc_info.push_back(tc_data);
978 }
979 tc_data.clear();
980 }
981 i = i + ntc_win - 1;
982 }
983
984 } else { // Summary off
985 if (tc_trg == true) { // summary off & TC on
986 for (int j = 0; j < ntc_win; j++) {
987 tc_id = (rdat[i + j + 2] >> 20) & 0x3FF;
988 tc_t = (rdat[i + j + 2] >> 12) & 0x7F;
989 conv_tc_t = (data_win - 3) * 128 + tc_t;
990 tc_e = rdat[i + j + 2] & 0xFFF;
991
992 // TC vector
993 tc_data.push_back(tc_id);
994 tc_data.push_back(conv_tc_t);
995 tc_data.push_back(tc_e);
996 tc_data.push_back(data_win);
997 if (tc_id < 81) {
998 if (tc_id > 75) {
999 tc_info_FE1.push_back(tc_data);
1000 } else {
1001 tc_info_FE2.push_back(tc_data);
1002 }
1003 } else if (tc_id > 512) {
1004 if (tc_id > 572) {
1005 tc_info_BE1.push_back(tc_data);
1006 } else {
1007 tc_info_BE2.push_back(tc_data);
1008 }
1009 } else {
1010 tc_info.push_back(tc_data);
1011 }
1012 tc_data.clear();
1013 }
1014 i = i + ntc_win + 1;
1015 } else { // Summary off & TC off
1016 i = i + 1;
1017 }
1018 }
1019 window_num++;
1020 }
1021
1022 // <---- Unpacking
1023
1024 // Summary
1025 /* cppcheck-suppress variableScope */
1026 int sum_num = 0;
1027 /* cppcheck-suppress variableScope */
1028 int sum_revo = 0;
1029 int cl_theta[6] = {0};
1030 int cl_phi[6] = {0};
1031 int cl_time[6] = { -9999};
1032 int cl_energy[6] = {0};
1033 int ncl = 0;
1034 int low_multi = 0;
1035 int b2bhabha_v = 0;
1036 int b2bhabha_s = 0;
1037 int mumu = 0;
1038 int prescale = 0;
1039 int icn_over = 0;
1040 int bg_veto = 0;
1041 int icn = 0;
1042 int etot_type = 0;
1043 int etot = 0;
1044 int b1_type = 0;
1045 int b1bhabha = 0;
1046 int physics = 0;
1047 int time_type = 0;
1048 int time = 0;
1049 int ecl_bst = 0;
1050
1051 int m_sumNum = 0;
1052
1053 TrgEclDataBase database;
1054 TrgEclMapping mapping;
1055
1056 vector<int> cl_1d;
1057 vector<vector<int>> cl_2d;
1058
1059 vector<int> evt_1d_vector;
1060 vector<vector<int>> evt_2d_vector;
1061
1062 // Store Summary
1063 int sum_size = sum_info.size();
1064 if (sum_size != 0) {
1065 for (int j = 0; j < sum_size; j++) {
1066 sum_num = sum_info[j][0];
1067 sum_revo = sum_info[j][1];
1068 if (etm_version >= 128) {
1069 ecl_bst = (sum_info[j][2] >> 26) & 0x1;
1070 }
1071 if (etm_version > 119) {
1072 cl_theta[5] = (sum_info[j][2] >> 19) & 0x7f;
1073 cl_phi[5] = (sum_info[j][2] >> 11) & 0xff;
1074 cl_time[5] = (sum_info[j][2] >> 3) & 0xff;
1075 cl_energy[5] = ((sum_info[j][2] & 0x7) << 9) + ((sum_info[j][3] >> 23) & 0x1ff);
1076
1077 cl_theta[4] = (sum_info[j][3] >> 16) & 0x7f;
1078 cl_phi[4] = (sum_info[j][3] >> 8) & 0xff;
1079 cl_time[4] = (sum_info[j][3]) & 0xff;
1080 cl_energy[4] = (sum_info[j][4] >> 20) & 0xfff;
1081
1082 cl_theta[3] = (sum_info[j][4] >> 13) & 0x7f;
1083 cl_phi[3] = (sum_info[j][4] >> 5) & 0xff;
1084 cl_time[3] = ((sum_info[j][4] & 0x1f) << 3) + ((sum_info[j][5] >> 29) & 0x7);
1085 cl_energy[3] = (sum_info[j][5] >> 17) & 0xfff;
1086
1087 cl_theta[2] = (sum_info[j][5] >> 10) & 0x7f;
1088 cl_phi[2] = (sum_info[j][5] >> 2) & 0xff;
1089 cl_time[2] = ((sum_info[j][5] & 0x3) << 6) + ((sum_info[j][6] >> 26) & 0x3f);
1090 cl_energy[2] = (sum_info[j][6] >> 14) & 0xfff;
1091
1092 cl_theta[1] = (sum_info[j][6] >> 7) & 0x7f;
1093 cl_phi[1] = ((sum_info[j][6] & 0x7f) << 1) + ((sum_info[j][7] >> 31) & 0x1);
1094 cl_time[1] = (sum_info[j][7] >> 23) & 0xff;
1095 cl_energy[1] = (sum_info[j][7] >> 11) & 0xfff;
1096
1097 cl_theta[0] = (sum_info[j][7] >> 4) & 0x7f;
1098 cl_phi[0] = ((sum_info[j][7] & 0xf) << 4) + ((sum_info[j][8] >> 28) & 0xf);
1099 cl_time[0] = (sum_info[j][8] >> 20) & 0xff;
1100 cl_energy[0] = (sum_info[j][8] >> 8) & 0xfff;
1101
1102 ncl = (sum_info[j][8] >> 5) & 0x7;
1103
1104 prescale = (sum_info[j][8] >> 4) & 0x1;
1105 mumu = (sum_info[j][8] >> 3) & 0x1;
1106 b2bhabha_s = (sum_info[j][8] >> 2) & 0x1;
1107 if (etm_version >= 135) {
1108 low_multi = (((sum_info[j][2] >> 27) & 0x3) << 12) + ((sum_info[j][8] & 0x3) << 10) + ((sum_info[j][9] >> 22) & 0x3ff);
1109 } else {
1110 low_multi = ((sum_info[j][8] & 0x3) << 10) + ((sum_info[j][9] >> 22) & 0x3ff);
1111 }
1112 b2bhabha_v = (sum_info[j][9] >> 21) & 0x1;
1113 icn_over = (sum_info[j][9] >> 20) & 0x1;
1114 bg_veto = (sum_info[j][9] >> 17) & 0x7;
1115 icn = (sum_info[j][9] >> 10) & 0x7f;
1116 etot_type = (sum_info[j][9] >> 7) & 0x7;
1117 etot = ((sum_info[j][9] & 0x7f) << 6) + ((sum_info[j][10] >> 26) & 0x3f);
1118 b1_type = (sum_info[j][10] >> 12) & 0x3fff;
1119 b1bhabha = (sum_info[j][10] >> 11) & 0x1;
1120 physics = (sum_info[j][10] >> 10) & 0x1;
1121 time_type = (sum_info[j][10] >> 7) & 0x7;
1122 time = (sum_info[j][10]) & 0x7f;
1123 } else {
1124 cl_theta[5] = (sum_info[j][2] >> 24) & 0x7f;
1125 cl_phi[5] = (sum_info[j][2] >> 16) & 0xff;
1126 cl_time[5] = (sum_info[j][2] >> 8) & 0xff;
1127 cl_energy[5] = ((sum_info[j][2] & 0xff) << 4) + ((sum_info[j][3] >> 28) & 0xf);
1128
1129 cl_theta[4] = (sum_info[j][3] >> 21) & 0x7f;
1130 cl_phi[4] = (sum_info[j][3] >> 13) & 0xff;
1131 cl_time[4] = (sum_info[j][3] >> 5) & 0xff;
1132 cl_energy[4] = ((sum_info[j][3] & 0x1f) << 7) + ((sum_info[j][4] >> 25) & 0x7f);
1133
1134 cl_theta[3] = (sum_info[j][4] >> 18) & 0x7f;
1135 cl_phi[3] = (sum_info[j][4] >> 10) & 0xff;
1136 cl_time[3] = (sum_info[j][4] >> 2) & 0xff;
1137 cl_energy[3] = ((sum_info[j][4] & 0x3) << 10) + ((sum_info[j][5] >> 22) & 0x3ff);
1138
1139 cl_theta[2] = (sum_info[j][5] >> 15) & 0x7f;
1140 cl_phi[2] = (sum_info[j][5] >> 7) & 0xff;
1141 cl_time[2] = ((sum_info[j][5] & 0x7f) << 1) + ((sum_info[j][6] >> 31) & 0x1);
1142 cl_energy[2] = (sum_info[j][6] >> 19) & 0xfff;
1143
1144 cl_theta[1] = (sum_info[j][6] >> 12) & 0x7f;
1145 cl_phi[1] = (sum_info[j][6] >> 4) & 0xff;
1146 cl_time[1] = ((sum_info[j][6] & 0xf) << 4) + ((sum_info[j][7] >> 28) & 0xf);
1147 cl_energy[1] = (sum_info[j][7] >> 16) & 0xfff;
1148
1149 cl_theta[0] = (sum_info[j][7] >> 9) & 0x7f;
1150 cl_phi[0] = (sum_info[j][7] >> 1) & 0xff;
1151 cl_time[0] = ((sum_info[j][7] & 0x1) << 7) + ((sum_info[j][8] >> 25) & 0x7f);
1152 cl_energy[0] = (sum_info[j][8] >> 13) & 0xfff;
1153
1154 ncl = (sum_info[j][8] >> 10) & 0x7;
1155
1156 low_multi = ((sum_info[j][8] & 0x3ff) << 2) + ((sum_info[j][9] >> 30) & 0x3);
1157 b2bhabha_v = (sum_info[j][9] >> 29) & 0x1;
1158 icn_over = (sum_info[j][9] >> 28) & 0x1;
1159 bg_veto = (sum_info[j][9] >> 25) & 0x7;
1160 icn = (sum_info[j][9] >> 18) & 0x7f;
1161 etot_type = (sum_info[j][9] >> 15) & 0x7;
1162 etot = (sum_info[j][9] >> 2) & 0x1fff;
1163
1164 b1_type = ((sum_info[j][9] & 0x3) << 12) + ((sum_info[j][10] >> 20) & 0xfff);
1165 b1bhabha = (sum_info[j][10] >> 19) & 0x1;
1166 physics = (sum_info[j][10] >> 18) & 0x1;
1167 time_type = (sum_info[j][10] >> 15) & 0x7;
1168 time = (sum_info[j][10] >> 8) & 0x7f;
1169
1170 b2bhabha_s = 0;
1171 mumu = 0;
1172 prescale = 0;
1173 }
1174
1175 m_TRGECLSumArray.appendNew();
1176 m_sumNum = m_TRGECLSumArray.getEntries() - 1;
1177 m_TRGECLSumArray[m_sumNum]->setEventId(n_basf2evt);
1178 m_TRGECLSumArray[m_sumNum]->setSumNum(sum_num);
1179 m_TRGECLSumArray[m_sumNum]->setSumRevo(sum_revo);
1180 m_TRGECLSumArray[m_sumNum]->setCLTheta(cl_theta);
1181 m_TRGECLSumArray[m_sumNum]->setCLPhi(cl_phi);
1182 m_TRGECLSumArray[m_sumNum]->setCLTime(cl_time);
1183 m_TRGECLSumArray[m_sumNum]->setCLEnergy(cl_energy);
1184 m_TRGECLSumArray[m_sumNum]->setNCL(ncl);
1185 m_TRGECLSumArray[m_sumNum]->setICN(icn);
1186 m_TRGECLSumArray[m_sumNum]->setICNOver(icn_over);
1187 m_TRGECLSumArray[m_sumNum]->setLowMulti(low_multi);
1188 m_TRGECLSumArray[m_sumNum]->set3DBhabhaV(b2bhabha_v);
1189 m_TRGECLSumArray[m_sumNum]->set3DBhabhaS(b2bhabha_s);
1190 m_TRGECLSumArray[m_sumNum]->setMumu(mumu);
1191 m_TRGECLSumArray[m_sumNum]->setPrescale(prescale);
1192 m_TRGECLSumArray[m_sumNum]->set2DBhabha(b1bhabha);
1193 m_TRGECLSumArray[m_sumNum]->setBhabhaType(b1_type);
1194 m_TRGECLSumArray[m_sumNum]->setPhysics(physics);
1195 m_TRGECLSumArray[m_sumNum]->setBG(bg_veto);
1196 m_TRGECLSumArray[m_sumNum]->setEtot(etot);
1197 m_TRGECLSumArray[m_sumNum]->setEtotType(etot_type);
1198 m_TRGECLSumArray[m_sumNum]->setECLBST(ecl_bst);
1199 m_TRGECLSumArray[m_sumNum]->setTime(time);
1200 m_TRGECLSumArray[m_sumNum]->setTimeType(time_type);
1201
1202 for (int k = 0; k < 6; k++) {
1203 cl_1d.push_back(cl_theta[k]);
1204 cl_1d.push_back(cl_phi[k]);
1205 cl_1d.push_back(cl_time[k]);
1206 cl_1d.push_back(cl_energy[k]);
1207 cl_2d.push_back(cl_1d);
1208 cl_1d.clear();
1209 }
1210 sort(cl_2d.begin(), cl_2d.end(),
1211 [](const vector<int>& aa1, const vector<int>& aa2) {return aa1[3] > aa2[3];});
1212
1213 evt_1d_vector.push_back(abs(sum_num - 3));
1214 evt_1d_vector.push_back(sum_revo);
1215 evt_1d_vector.push_back(sum_num);
1216 evt_1d_vector.push_back(time);
1217 for (int k = 0; k < 6; k++) {
1218 evt_1d_vector.push_back(cl_2d[k][0]);
1219 evt_1d_vector.push_back(cl_2d[k][1]);
1220 evt_1d_vector.push_back(cl_2d[k][2]);
1221 evt_1d_vector.push_back(cl_2d[k][3]);
1222 }
1223 evt_1d_vector.push_back(ncl);
1224 evt_1d_vector.push_back(low_multi);
1225 evt_1d_vector.push_back(b2bhabha_v);
1226 evt_1d_vector.push_back(b2bhabha_s);
1227 evt_1d_vector.push_back(mumu);
1228 evt_1d_vector.push_back(prescale);
1229 evt_1d_vector.push_back(icn);
1230 evt_1d_vector.push_back(icn_over);
1231 evt_1d_vector.push_back(etot_type);
1232 evt_1d_vector.push_back(etot);
1233 evt_1d_vector.push_back(ecl_bst);
1234 evt_1d_vector.push_back(b1_type);
1235 evt_1d_vector.push_back(b1bhabha);
1236 evt_1d_vector.push_back(physics);
1237 evt_1d_vector.push_back(time_type);
1238 evt_2d_vector.push_back(evt_1d_vector);
1239 evt_1d_vector.clear();
1240 }
1241 } else {
1242
1243 for (int k = 0; k < 6; k++) {
1244 cl_theta[k] = 0;
1245 cl_phi[k] = 0;
1246 cl_time[k] = -9999;
1247 cl_energy[k] = 0;
1248 }
1249 ncl = 0;
1250 low_multi = 0;
1251 b2bhabha_v = 0;
1252 b2bhabha_s = 0;
1253 mumu = 0;
1254 prescale = 0;
1255 icn_over = 0;
1256 bg_veto = 0;
1257 icn = 0;
1258 etot_type = 0;
1259 etot = 0;
1260 ecl_bst = 0;
1261 b1_type = 0;
1262 b1bhabha = 0;
1263 physics = 0;
1264 time_type = 0;
1265 time = -9999;
1266
1267 m_TRGECLSumArray.appendNew();
1268 m_sumNum = m_TRGECLSumArray.getEntries() - 1;
1269 m_TRGECLSumArray[m_sumNum]->setEventId(n_basf2evt);
1270 m_TRGECLSumArray[m_sumNum]->setSumNum(0);
1271 m_TRGECLSumArray[m_sumNum]->setCLTheta(cl_theta);
1272 m_TRGECLSumArray[m_sumNum]->setCLPhi(cl_phi);
1273 m_TRGECLSumArray[m_sumNum]->setCLTime(cl_time);
1274 m_TRGECLSumArray[m_sumNum]->setCLEnergy(cl_energy);
1275 m_TRGECLSumArray[m_sumNum]->setNCL(ncl);
1276 m_TRGECLSumArray[m_sumNum]->setICN(icn);
1277 m_TRGECLSumArray[m_sumNum]->setICNOver(icn_over);
1278 m_TRGECLSumArray[m_sumNum]->setLowMulti(low_multi);
1279 m_TRGECLSumArray[m_sumNum]->set3DBhabhaV(b2bhabha_v);
1280 m_TRGECLSumArray[m_sumNum]->set3DBhabhaS(b2bhabha_s);
1281 m_TRGECLSumArray[m_sumNum]->setMumu(mumu);
1282 m_TRGECLSumArray[m_sumNum]->setPrescale(prescale);
1283 m_TRGECLSumArray[m_sumNum]->set2DBhabha(b1bhabha);
1284 m_TRGECLSumArray[m_sumNum]->setBhabhaType(b1_type);
1285 m_TRGECLSumArray[m_sumNum]->setPhysics(physics);
1286 m_TRGECLSumArray[m_sumNum]->setBG(bg_veto);
1287 m_TRGECLSumArray[m_sumNum]->setEtot(etot);
1288 m_TRGECLSumArray[m_sumNum]->setEtotType(etot_type);
1289 m_TRGECLSumArray[m_sumNum]->setECLBST(ecl_bst);
1290 m_TRGECLSumArray[m_sumNum]->setTime(time);
1291 m_TRGECLSumArray[m_sumNum]->setTimeType(time_type);
1292 }
1293
1294 // TC & TRG
1295 tc_info.insert(tc_info.end(), tc_info_FE1.begin(), tc_info_FE1.end());
1296 tc_info.insert(tc_info.end(), tc_info_FE2.begin(), tc_info_FE2.end());
1297 tc_info.insert(tc_info.end(), tc_info_BE1.begin(), tc_info_BE1.end());
1298 tc_info.insert(tc_info.end(), tc_info_BE2.begin(), tc_info_BE2.end());
1299
1300 int m_evtNum = 0;
1301
1302 int m_tcNum = 0;
1303 /* cppcheck-suppress variableScope */
1304 int m_tcid = 0;
1305 /* cppcheck-suppress variableScope */
1306 int m_time = -9999;
1307 /* cppcheck-suppress variableScope */
1308 int m_energy = 0;
1309 /* cppcheck-suppress variableScope */
1310 int m_win = 0;
1311 /* cppcheck-suppress variableScope */
1312 int m_revo = 0;
1313 /* cppcheck-suppress variableScope */
1314 int m_caltime = -9999;
1315
1316 int tot_ntc = tc_info.size();
1317 /* cppcheck-suppress variableScope */
1318 int evt_ntc = 0;
1319 /* cppcheck-suppress variableScope */
1320 int evt_revo = -9999;
1321 /* cppcheck-suppress variableScope */
1322 int evt_win = 0;
1323 /* cppcheck-suppress variableScope */
1324 int evt_timing = -9999;
1325 int evt_cl_theta[6] = {0};
1326 int evt_cl_phi[6] = {0};
1327 int evt_cl_time[6] = { -9999};
1328 int evt_cl_energy[6] = {0};
1329 /* cppcheck-suppress variableScope */
1330 int evt_ncl = 0;
1331 /* cppcheck-suppress variableScope */
1332 int evt_low_multi = 0;
1333 /* cppcheck-suppress variableScope */
1334 int evt_b2bhabha_v = 0;
1335 /* cppcheck-suppress variableScope */
1336 int evt_b2bhabha_s = 0;
1337 /* cppcheck-suppress variableScope */
1338 int evt_mumu = 0;
1339 /* cppcheck-suppress variableScope */
1340 int evt_prescale = 0;
1341 /* cppcheck-suppress variableScope */
1342 int evt_icn = 0;
1343 /* cppcheck-suppress variableScope */
1344 int evt_icn_over = 0;
1345 /* cppcheck-suppress variableScope */
1346 int evt_etot_type = 0;
1347 /* cppcheck-suppress variableScope */
1348 int evt_etot = 0;
1349 /* cppcheck-suppress variableScope */
1350 int evt_ecl_bst = 0;
1351 /* cppcheck-suppress variableScope */
1352 int evt_b1_type = 0;
1353 /* cppcheck-suppress variableScope */
1354 int evt_b1bhabha = 0;
1355 /* cppcheck-suppress variableScope */
1356 int evt_physics = 0;
1357 /* cppcheck-suppress variableScope */
1358 int evt_time_type = 0;
1359 /* cppcheck-suppress variableScope */
1360 int evt_etot_all = 0;
1361 /* cppcheck-suppress variableScope */
1362 int evt_time_min = 0;
1363 /* cppcheck-suppress variableScope */
1364 int evt_time_max = 0;
1365 /* cppcheck-suppress variableScope */
1366 int evt_time_win = 0;
1367 /* cppcheck-suppress variableScope */
1368 int etot_i = 0;
1369 /* cppcheck-suppress variableScope */
1370 int etot_c = 0;
1371 /* cppcheck-suppress variableScope */
1372 int etot_f = 0;
1373 /* cppcheck-suppress variableScope */
1374 int cl_tcid = 0;
1375 /* cppcheck-suppress variableScope */
1376 int cl_thetaid = 0;
1377 /* cppcheck-suppress variableScope */
1378 int cl_phiid = 0;
1379 /* cppcheck-suppress variableScope */
1380 int m_clNum = 0;
1381
1382
1383 int evt_v_size = evt_2d_vector.size();
1384 if (evt_v_size != 0) {
1385 // Sort window : 3 => 4 => 2 => 5 => 1 => 6 => 7
1386 sort(evt_2d_vector.begin(), evt_2d_vector.end(),
1387 [](const vector<int>& aa1, const vector<int>& aa2) {return aa1[0] < aa2[0];});
1388 }
1389
1390 if (tot_ntc != 0 && flag_checksum == 0 && nnn > 7) {
1391 if (evt_v_size == 0) {
1392 // Find most energetic TC timing
1393 sort(tc_info.begin(), tc_info.end(),
1394 [](const vector<int>& aa1, const vector<int>& aa2) {return aa1[2] > aa2[2];});
1395 evt_revo = win3_revo;
1396 evt_win = tc_info[0][3];
1397 evt_timing = tc_info[0][1];
1398 for (int k = 0; k < 6; k++) {
1399 evt_cl_theta[k] = 0;
1400 evt_cl_phi[k] = 0;
1401 evt_cl_time[k] = 0;
1402 evt_cl_energy[k] = 0;
1403 }
1404 evt_ncl = 0;
1405 evt_low_multi = 0;
1406 evt_b2bhabha_v = 0;
1407 evt_b2bhabha_s = 0;
1408 evt_mumu = 0;
1409 evt_prescale = 0;
1410 evt_icn = 0;
1411 evt_icn_over = 0;
1412 evt_etot_type = 0;
1413 evt_etot = 0;
1414 evt_ecl_bst = 0;
1415 evt_b1_type = 0;
1416 evt_b1bhabha = 0;
1417 evt_physics = 0;
1418 evt_time_type = 0;
1419 } else {
1420 evt_revo = evt_2d_vector[0][1];
1421 evt_win = evt_2d_vector[0][2];
1422 evt_timing = evt_2d_vector[0][3];
1423 for (int k = 0; k < 6; k++) {
1424 evt_cl_theta[k] = evt_2d_vector[0][4 + k * 4];
1425 evt_cl_phi[k] = evt_2d_vector[0][5 + k * 4];
1426 evt_cl_time[k] = evt_2d_vector[0][6 + k * 4];
1427 evt_cl_energy[k] = evt_2d_vector[0][7 + k * 4];
1428 }
1429 evt_ncl = evt_2d_vector[0][28];
1430 evt_low_multi = evt_2d_vector[0][29];
1431 evt_b2bhabha_v = evt_2d_vector[0][30];
1432 evt_b2bhabha_s = evt_2d_vector[0][31];
1433 evt_mumu = evt_2d_vector[0][32];
1434 evt_prescale = evt_2d_vector[0][33];
1435 evt_icn = evt_2d_vector[0][34];
1436 evt_icn_over = evt_2d_vector[0][35];
1437 evt_etot_type = evt_2d_vector[0][36];
1438 evt_etot = evt_2d_vector[0][37];
1439 evt_ecl_bst = evt_2d_vector[0][38];
1440 evt_b1_type = evt_2d_vector[0][39];
1441 evt_b1bhabha = evt_2d_vector[0][40];
1442 evt_physics = evt_2d_vector[0][41];
1443 evt_time_type = evt_2d_vector[0][42];
1444 }
1445 // Sort by TC number
1446 sort(tc_info.begin(), tc_info.end(),
1447 [](const vector<int>& aa1, const vector<int>& aa2) {return aa1[0] < aa2[0];});
1448
1449 for (int ii = 0; ii < tot_ntc; ii++) {
1450 m_tcid = tc_info[ii][0];
1451 m_time = tc_info[ii][1];
1452 m_energy = tc_info[ii][2];
1453 m_win = tc_info[ii][3];
1454 m_revo = win3_revo;
1455 m_caltime = m_time - ((evt_win - 3) * 128 + evt_timing);
1456 m_TRGECLTCArray.appendNew();
1457 m_tcNum = m_TRGECLTCArray.getEntries() - 1;
1458 m_TRGECLTCArray[m_tcNum]->setEventId(n_basf2evt);
1459 m_TRGECLTCArray[m_tcNum]->setTCId(m_tcid);
1460 m_TRGECLTCArray[m_tcNum]->setTCTime(m_time);
1461 m_TRGECLTCArray[m_tcNum]->setTCCALTime(m_caltime);
1462 m_TRGECLTCArray[m_tcNum]->setHitWin(m_win);
1463 m_TRGECLTCArray[m_tcNum]->setRevoFAM(m_revo);
1464 m_TRGECLTCArray[m_tcNum]->setTCEnergy(m_energy);
1465 m_TRGECLTCArray[m_tcNum]->setChecksum(flag_checksum);
1466
1467 if (m_win == evt_win || m_win == evt_win + 1) evt_ntc++;
1468 if (m_win == evt_win - 1) {
1469 etot_i += m_energy;
1470 }
1471 if (m_win == evt_win) {
1472 etot_c += m_energy;
1473 }
1474 if (m_win == evt_win + 1) {
1475 etot_f += m_energy;
1476 }
1477 }
1478
1479 if (etot_i == 0 && etot_f == 0) {
1480 evt_etot_all = etot_c;
1481 evt_time_min = - evt_timing;
1482 evt_time_max = 256 - evt_timing;
1483 evt_time_win = 1;
1484 } else if (etot_i >= etot_f) {
1485 evt_etot_all = etot_c + etot_i;
1486 evt_time_min = -128 - evt_timing;
1487 evt_time_max = 128 - evt_timing;
1488 evt_time_win = -1;
1489 } else {
1490 evt_etot_all = etot_c + etot_f;
1491 evt_time_min = - evt_timing;
1492 evt_time_max = 256 - evt_timing;
1493 evt_time_win = 1;
1494 }
1495
1496 for (int icluster = 0; icluster < 6; icluster++) {
1497 if (evt_cl_energy[icluster] == 0 || evt_cl_theta[icluster] == 0 || evt_cl_phi[icluster] == 0) {continue;}
1498 cl_tcid = mapping.getTCIdFromPosition(evt_cl_theta[icluster], evt_cl_phi[icluster]);
1499 if (cl_tcid == 0) {continue;}
1500 cl_thetaid = mapping.getTCThetaIdFromTCId(cl_tcid);
1501 cl_phiid = mapping.getTCPhiIdFromTCId(cl_tcid);
1502
1503 m_TRGECLClusterArray.appendNew();
1504 m_clNum = m_TRGECLClusterArray.getEntries() - 1;
1505 m_TRGECLClusterArray[m_clNum]->setEventId(n_basf2evt);
1506 m_TRGECLClusterArray[m_clNum]->setClusterId(icluster);
1507 m_TRGECLClusterArray[m_clNum]->setEventRevo(evt_revo);
1508
1509 m_TRGECLClusterArray[m_clNum]->setMaxTCId(cl_tcid); // center of Cluster
1510 m_TRGECLClusterArray[m_clNum]->setMaxThetaId(cl_thetaid);
1511 m_TRGECLClusterArray[m_clNum]->setMaxPhiId(cl_phiid);
1512 m_TRGECLClusterArray[m_clNum]->setClusterId(icluster);
1513 m_TRGECLClusterArray[m_clNum]->setEnergyDep((double)evt_cl_energy[icluster] * 5.25); // MeV
1514 m_TRGECLClusterArray[m_clNum]->setTimeAve((double)evt_cl_time[icluster]);
1515 m_TRGECLClusterArray[m_clNum]->setPositionX(mapping.getTCPosition(cl_tcid).X());
1516 m_TRGECLClusterArray[m_clNum]->setPositionY(mapping.getTCPosition(cl_tcid).Y());
1517 m_TRGECLClusterArray[m_clNum]->setPositionZ(mapping.getTCPosition(cl_tcid).Z());
1518 }
1519 m_TRGECLEvtArray.appendNew();
1520 m_evtNum = m_TRGECLEvtArray.getEntries() - 1;
1521 m_TRGECLEvtArray[m_evtNum]->setEventId(n_basf2evt);
1522 m_TRGECLEvtArray[m_evtNum]->setETM(etm_version);
1523 m_TRGECLEvtArray[m_evtNum]->setL1Revo(l1_revo);
1524 m_TRGECLEvtArray[m_evtNum]->setEvtRevo(evt_revo);
1525 m_TRGECLEvtArray[m_evtNum]->setEvtWin(evt_win);
1526 m_TRGECLEvtArray[m_evtNum]->setEvtTime(evt_timing);
1527 m_TRGECLEvtArray[m_evtNum]->setNTC(evt_ntc);
1528 m_TRGECLEvtArray[m_evtNum]->setCLTheta(evt_cl_theta);
1529 m_TRGECLEvtArray[m_evtNum]->setCLPhi(evt_cl_phi);
1530 m_TRGECLEvtArray[m_evtNum]->setCLTime(evt_cl_time);
1531 m_TRGECLEvtArray[m_evtNum]->setCLEnergy(evt_cl_energy);
1532 m_TRGECLEvtArray[m_evtNum]->setNCL(evt_ncl);
1533 m_TRGECLEvtArray[m_evtNum]->setLowMulti(evt_low_multi);
1534 m_TRGECLEvtArray[m_evtNum]->set3DBhabhaV(evt_b2bhabha_v);
1535 m_TRGECLEvtArray[m_evtNum]->set3DBhabhaS(evt_b2bhabha_s);
1536 m_TRGECLEvtArray[m_evtNum]->setMumu(evt_mumu);
1537 m_TRGECLEvtArray[m_evtNum]->setPrescale(evt_prescale);
1538 m_TRGECLEvtArray[m_evtNum]->setICN(evt_icn);
1539 m_TRGECLEvtArray[m_evtNum]->setICNOver(evt_icn_over);
1540 m_TRGECLEvtArray[m_evtNum]->setEtotType(evt_etot_type);
1541 m_TRGECLEvtArray[m_evtNum]->setEtot(evt_etot);
1542 m_TRGECLEvtArray[m_evtNum]->setECLBST(evt_ecl_bst);
1543 m_TRGECLEvtArray[m_evtNum]->set2DBhabha(evt_b1bhabha);
1544 m_TRGECLEvtArray[m_evtNum]->setBhabhaType(evt_b1_type);
1545 m_TRGECLEvtArray[m_evtNum]->setPhysics(evt_physics);
1546 m_TRGECLEvtArray[m_evtNum]->setTimeType(evt_time_type);
1547 m_TRGECLEvtArray[m_evtNum]->setCheckSum(flag_checksum);
1548 m_TRGECLEvtArray[m_evtNum]->setEvtExist(1);
1549 m_TRGECLEvtArray[m_evtNum]->setTRGTYPE(trgtype);
1550 m_TRGECLEvtArray[m_evtNum]->setEtotAll(evt_etot_all);
1551 m_TRGECLEvtArray[m_evtNum]->setEvtTimeMin(evt_time_min);
1552 m_TRGECLEvtArray[m_evtNum]->setEvtTimeMax(evt_time_max);
1553 m_TRGECLEvtArray[m_evtNum]->setEvtTimeWin(evt_time_win);
1554 } else {
1555 m_TRGECLTCArray.appendNew();
1556 m_tcNum = m_TRGECLTCArray.getEntries() - 1;
1557 m_TRGECLTCArray[m_tcNum]->setEventId(n_basf2evt);
1558 m_TRGECLTCArray[m_tcNum]->setTCId(0);
1559 m_TRGECLTCArray[m_tcNum]->setTCTime(-9999);
1560 m_TRGECLTCArray[m_tcNum]->setTCCALTime(-9999);
1561 m_TRGECLTCArray[m_tcNum]->setHitWin(-9999);
1562 m_TRGECLTCArray[m_tcNum]->setRevoFAM(-9999);
1563 m_TRGECLTCArray[m_tcNum]->setTCEnergy(0);
1564 m_TRGECLTCArray[m_tcNum]->setChecksum(flag_checksum);
1565
1566 m_TRGECLEvtArray.appendNew();
1567 m_evtNum = m_TRGECLEvtArray.getEntries() - 1;
1568 m_TRGECLEvtArray[m_evtNum]->setEventId(n_basf2evt);
1569 m_TRGECLEvtArray[m_evtNum]->setETM(etm_version);
1570 m_TRGECLEvtArray[m_evtNum]->setL1Revo(-9999);
1571 m_TRGECLEvtArray[m_evtNum]->setEvtTime(-9999);
1572 m_TRGECLEvtArray[m_evtNum]->setEvtRevo(-9999);
1573 m_TRGECLEvtArray[m_evtNum]->setEvtWin(-9999);
1574 m_TRGECLEvtArray[m_evtNum]->setNTC(0);
1575 for (int k = 0; k < 6; k++) {
1576 evt_cl_theta[k] = 0;
1577 evt_cl_phi[k] = 0;
1578 evt_cl_time[k] = -9999;
1579 evt_cl_energy[k] = 0;
1580 }
1581 m_TRGECLEvtArray[m_evtNum]->setCLTheta(evt_cl_theta);
1582 m_TRGECLEvtArray[m_evtNum]->setCLPhi(evt_cl_phi);
1583 m_TRGECLEvtArray[m_evtNum]->setCLTime(evt_cl_time);
1584 m_TRGECLEvtArray[m_evtNum]->setCLEnergy(evt_cl_energy);
1585 m_TRGECLEvtArray[m_evtNum]->setNCL(0);
1586 m_TRGECLEvtArray[m_evtNum]->setLowMulti(0);
1587 m_TRGECLEvtArray[m_evtNum]->set3DBhabhaV(0);
1588 m_TRGECLEvtArray[m_evtNum]->set3DBhabhaS(0);
1589 m_TRGECLEvtArray[m_evtNum]->setMumu(0);
1590 m_TRGECLEvtArray[m_evtNum]->setPrescale(0);
1591 m_TRGECLEvtArray[m_evtNum]->setICN(0);
1592 m_TRGECLEvtArray[m_evtNum]->setICNOver(0);
1593 m_TRGECLEvtArray[m_evtNum]->setEtotType(0);
1594 m_TRGECLEvtArray[m_evtNum]->setEtot(0);
1595 m_TRGECLEvtArray[m_evtNum]->setECLBST(0);
1596 m_TRGECLEvtArray[m_evtNum]->set2DBhabha(0);
1597 m_TRGECLEvtArray[m_evtNum]->setBhabhaType(0);
1598 m_TRGECLEvtArray[m_evtNum]->setPhysics(0);
1599 m_TRGECLEvtArray[m_evtNum]->setTimeType(0);
1600 m_TRGECLEvtArray[m_evtNum]->setCheckSum(flag_checksum);
1601 m_TRGECLEvtArray[m_evtNum]->setEvtExist(0);
1602 m_TRGECLEvtArray[m_evtNum]->setTRGTYPE(trgtype);
1603 m_TRGECLEvtArray[m_evtNum]->setEtotAll(0);
1604 m_TRGECLEvtArray[m_evtNum]->setEvtTimeMin(-9999);
1605 m_TRGECLEvtArray[m_evtNum]->setEvtTimeMax(-9999);
1606 m_TRGECLEvtArray[m_evtNum]->setEvtTimeWin(0);
1607 }
1608
1609 return;
1610}

◆ clone()

std::shared_ptr< PathElement > clone ( ) const
overridevirtualinherited

Create an independent copy of this module.

Note that parameters are shared, so changing them on a cloned module will also affect the original module.

Implements PathElement.

Definition at line 179 of file Module.cc.

180{
182 newModule->m_moduleParamList.setParameters(getParamList());
183 newModule->setName(getName());
184 newModule->m_package = m_package;
185 newModule->m_propertyFlags = m_propertyFlags;
186 newModule->m_logConfig = m_logConfig;
187 newModule->m_conditions = m_conditions;
188
189 return newModule;
190}
std::shared_ptr< Module > registerModule(const std::string &moduleName, std::string sharedLibPath="") noexcept(false)
Creates an instance of a module and registers it to the ModuleManager.
static ModuleManager & Instance()
Exception is thrown if the requested module could not be created by the ModuleManager.
const ModuleParamList & getParamList() const
Return module param list.
Definition: Module.h:363
const std::string & getName() const
Returns the name of the module.
Definition: Module.h:187
const std::string & getType() const
Returns the type of the module (i.e.
Definition: Module.cc:41
unsigned int m_propertyFlags
The properties of the module as bitwise or (with |) of EModulePropFlags.
Definition: Module.h:512
LogConfig m_logConfig
The log system configuration of the module.
Definition: Module.h:514
std::vector< ModuleCondition > m_conditions
Module condition, only non-null if set.
Definition: Module.h:521
std::string m_package
Package this module is found in (may be empty).
Definition: Module.h:510
std::shared_ptr< Module > ModulePtr
Defines a pointer to a module object as a boost shared pointer.
Definition: Module.h:43

◆ def_beginRun()

virtual void def_beginRun ( )
inlineprotectedvirtualinherited

Wrapper method for the virtual function beginRun() that has the implementation to be used in a call from Python.

Reimplemented in PyModule.

Definition at line 426 of file Module.h.

426{ beginRun(); }
virtual void beginRun()
Called when entering a new run.
Definition: Module.h:147

◆ def_endRun()

virtual void def_endRun ( )
inlineprotectedvirtualinherited

This method can receive that the current run ends as a call from the Python side.

For regular C++-Modules that forwards the call to the regular endRun() method.

Reimplemented in PyModule.

Definition at line 439 of file Module.h.

439{ endRun(); }
virtual void endRun()
This method is called if the current run ends.
Definition: Module.h:166

◆ def_event()

virtual void def_event ( )
inlineprotectedvirtualinherited

Wrapper method for the virtual function event() that has the implementation to be used in a call from Python.

Reimplemented in PyModule.

Definition at line 432 of file Module.h.

432{ event(); }
virtual void event()
This method is the core of the module.
Definition: Module.h:157

◆ def_initialize()

virtual void def_initialize ( )
inlineprotectedvirtualinherited

Wrappers to make the methods without "def_" prefix callable from Python.

Overridden in PyModule. Wrapper method for the virtual function initialize() that has the implementation to be used in a call from Python.

Reimplemented in PyModule.

Definition at line 420 of file Module.h.

420{ initialize(); }
virtual void initialize()
Initialize the Module.
Definition: Module.h:109

◆ def_terminate()

virtual void def_terminate ( )
inlineprotectedvirtualinherited

Wrapper method for the virtual function terminate() that has the implementation to be used in a call from Python.

Reimplemented in PyModule.

Definition at line 445 of file Module.h.

445{ terminate(); }
virtual void terminate()
This method is called at the end of the event processing.
Definition: Module.h:176

◆ endRun()

void endRun ( void  )
overridevirtual

Called when run ended.

Reimplemented from Module.

Definition at line 55 of file trgeclUnpackerModule.cc.

55{}

◆ evalCondition()

bool evalCondition ( ) const
inherited

If at least one condition was set, it is evaluated and true returned if at least one condition returns true.

If no condition or result value was defined, the method returns false. Otherwise, the condition is evaluated and true returned, if at least one condition returns true. To speed up the evaluation, the condition strings were already parsed in the method if_value().

Returns
True if at least one condition and return value exists and at least one condition expression was evaluated to true.

Definition at line 96 of file Module.cc.

97{
98 if (m_conditions.empty()) return false;
99
100 //okay, a condition was set for this Module...
101 if (!m_hasReturnValue) {
102 B2FATAL("A condition was set for '" << getName() << "', but the module did not set a return value!");
103 }
104
105 for (const auto& condition : m_conditions) {
106 if (condition.evaluate(m_returnValue)) {
107 return true;
108 }
109 }
110 return false;
111}
int m_returnValue
The return value.
Definition: Module.h:519
bool m_hasReturnValue
True, if the return value is set.
Definition: Module.h:518

◆ event()

void event ( void  )
overridevirtual

Called event by event.

Reimplemented from Module.

Definition at line 56 of file trgeclUnpackerModule.cc.

57{
58
59 StoreArray<RawTRG> raw_trgarray;
60
61 for (int i = 0; i < raw_trgarray.getEntries(); i++) { // # of readout boards
62 iFiness = i;
63 for (int j = 0; j < raw_trgarray[i]->GetNumEntries(); j++) { // Basically 1 entry
64 nodeid = ((raw_trgarray[i]->GetNodeID(j)) >> 24) & 0x1F;
65 trgtype = raw_trgarray[i]->GetTRGType(j);
66 n_basf2evt = raw_trgarray[i]->GetEveNo(j);
67 if (nodeid == 0x13) {
68 for (int ch = 0; ch < raw_trgarray[i]->GetMaxNumOfCh(j); ch++) { // ch in a readout board
69 nwords = raw_trgarray[i]->GetDetectorNwords(j, ch);
70 if (nwords == 0) {
71 continue; // This channel might be masked.
72 } else if (nwords < 9) {
73 B2ERROR("Consistecy error in unpacker.");
74 B2ERROR("data length " << nwords << " nWord " << nwords);
75 B2ERROR("Node ID " << nodeid << ", Finness ID " << iFiness);
76 continue;
77 }
78 readCOPPEREvent(raw_trgarray[i], j, nwords, ch);
79 }
80 }
81 }
82 }
83
84 // Count number of trigger cells in each ECL region for EventLevelClusteringInfo
85 uint16_t nTCsPerRegion[3] = {};
86 const int firstBarrelId = 81; // First TCId in the barrel
87 const int lastBarrelId = 512; // Last TCId in the barrel
88 for (auto& trgeclhit : m_TRGECLTCArray) {
89 const int tcId = trgeclhit.getTCId();
90 if (tcId < firstBarrelId) {
91 nTCsPerRegion[0]++;
92 } else if (tcId > lastBarrelId) {
93 nTCsPerRegion[2]++;
94 } else {
95 nTCsPerRegion[1]++;
96 }
97 }
98
99 // Store
101 m_eventLevelClusteringInfo->setNECLTriggerCellsFWD(nTCsPerRegion[0]);
102 m_eventLevelClusteringInfo->setNECLTriggerCellsBarrel(nTCsPerRegion[1]);
103 m_eventLevelClusteringInfo->setNECLTriggerCellsBWD(nTCsPerRegion[2]);
104
105}
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
virtual void readCOPPEREvent(RawTRG *, int, int, int)
Read data from TRG copper.
StoreObjPtr< EventLevelClusteringInfo > m_eventLevelClusteringInfo
EventLevelClusteringInfo.

◆ exposePythonAPI()

void exposePythonAPI ( )
staticinherited

Exposes methods of the Module class to Python.

Definition at line 325 of file Module.cc.

326{
327 // to avoid confusion between std::arg and boost::python::arg we want a shorthand namespace as well
328 namespace bp = boost::python;
329
330 docstring_options options(true, true, false); //userdef, py sigs, c++ sigs
331
332 void (Module::*setReturnValueInt)(int) = &Module::setReturnValue;
333
334 enum_<Module::EAfterConditionPath>("AfterConditionPath",
335 R"(Determines execution behaviour after a conditional path has been executed:
336
337.. attribute:: END
338
339 End processing of this path after the conditional path. (this is the default for if_value() etc.)
340
341.. attribute:: CONTINUE
342
343 After the conditional path, resume execution after this module.)")
344 .value("END", Module::EAfterConditionPath::c_End)
345 .value("CONTINUE", Module::EAfterConditionPath::c_Continue)
346 ;
347
348 /* Do not change the names of >, <, ... we use them to serialize conditional pathes */
349 enum_<Belle2::ModuleCondition::EConditionOperators>("ConditionOperator")
356 ;
357
358 enum_<Module::EModulePropFlags>("ModulePropFlags",
359 R"(Flags to indicate certain low-level features of modules, see :func:`Module.set_property_flags()`, :func:`Module.has_properties()`. Most useful flags are:
360
361.. attribute:: PARALLELPROCESSINGCERTIFIED
362
363 This module can be run in parallel processing mode safely (All I/O must be done through the data store, in particular, the module must not write any files.)
364
365.. attribute:: HISTOGRAMMANAGER
366
367 This module is used to manage histograms accumulated by other modules
368
369.. attribute:: TERMINATEINALLPROCESSES
370
371 When using parallel processing, call this module's terminate() function in all processes. This will also ensure that there is exactly one process (single-core if no parallel modules found) or at least one input, one main and one output process.
372)")
373 .value("INPUT", Module::EModulePropFlags::c_Input)
374 .value("OUTPUT", Module::EModulePropFlags::c_Output)
375 .value("PARALLELPROCESSINGCERTIFIED", Module::EModulePropFlags::c_ParallelProcessingCertified)
376 .value("HISTOGRAMMANAGER", Module::EModulePropFlags::c_HistogramManager)
377 .value("INTERNALSERIALIZER", Module::EModulePropFlags::c_InternalSerializer)
378 .value("TERMINATEINALLPROCESSES", Module::EModulePropFlags::c_TerminateInAllProcesses)
379 ;
380
381 //Python class definition
382 class_<Module, PyModule> module("Module", R"(
383Base class for Modules.
384
385A module is the smallest building block of the framework.
386A typical event processing chain consists of a Path containing
387modules. By inheriting from this base class, various types of
388modules can be created. To use a module, please refer to
389:func:`Path.add_module()`. A list of modules is available by running
390``basf2 -m`` or ``basf2 -m package``, detailed information on parameters is
391given by e.g. ``basf2 -m RootInput``.
392
393The 'Module Development' section in the manual provides detailed information
394on how to create modules, setting parameters, or using return values/conditions:
395https://xwiki.desy.de/xwiki/rest/p/f4fa4/#HModuleDevelopment
396
397)");
398 module
399 .def("__str__", &Module::getPathString)
400 .def("name", &Module::getName, return_value_policy<copy_const_reference>(),
401 "Returns the name of the module. Can be changed via :func:`set_name() <Module.set_name()>`, use :func:`type() <Module.type()>` for identifying a particular module class.")
402 .def("type", &Module::getType, return_value_policy<copy_const_reference>(),
403 "Returns the type of the module (i.e. class name minus 'Module')")
404 .def("set_name", &Module::setName, args("name"), R"(
405Set custom name, e.g. to distinguish multiple modules of the same type.
406
407>>> path.add_module('EventInfoSetter')
408>>> ro = path.add_module('RootOutput', branchNames=['EventMetaData'])
409>>> ro.set_name('RootOutput_metadata_only')
410>>> print(path)
411[EventInfoSetter -> RootOutput_metadata_only]
412
413)")
414 .def("description", &Module::getDescription, return_value_policy<copy_const_reference>(),
415 "Returns the description of this module.")
416 .def("package", &Module::getPackage, return_value_policy<copy_const_reference>(),
417 "Returns the package this module belongs to.")
418 .def("available_params", &_getParamInfoListPython,
419 "Return list of all module parameters as `ModuleParamInfo` instances")
420 .def("has_properties", &Module::hasProperties, (bp::arg("properties")),
421 R"DOCSTRING(Allows to check if the module has the given properties out of `ModulePropFlags` set.
422
423>>> if module.has_properties(ModulePropFlags.PARALLELPROCESSINGCERTIFIED):
424>>> ...
425
426Parameters:
427 properties (int): bitmask of `ModulePropFlags` to check for.
428)DOCSTRING")
429 .def("set_property_flags", &Module::setPropertyFlags, args("property_mask"),
430 "Set module properties in the form of an OR combination of `ModulePropFlags`.");
431 {
432 // python signature is too crowded, make ourselves
433 docstring_options subOptions(true, false, false); //userdef, py sigs, c++ sigs
434 module
435 .def("if_value", &Module::if_value,
436 (bp::arg("expression"), bp::arg("condition_path"), bp::arg("after_condition_path")= Module::EAfterConditionPath::c_End),
437 R"DOCSTRING(if_value(expression, condition_path, after_condition_path=AfterConditionPath.END)
438
439Sets a conditional sub path which will be executed after this
440module if the return value set in the module passes the given ``expression``.
441
442Modules can define a return value (int or bool) using ``setReturnValue()``,
443which can be used in the steering file to split the Path based on this value, for example
444
445>>> module_with_condition.if_value("<1", another_path)
446
447In case the return value of the ``module_with_condition`` for a given event is
448less than 1, the execution will be diverted into ``another_path`` for this event.
449
450You could for example set a special return value if an error occurs, and divert
451the execution into a path containing :b2:mod:`RootOutput` if it is found;
452saving only the data producing/produced by the error.
453
454After a conditional path has executed, basf2 will by default stop processing
455the path for this event. This behaviour can be changed by setting the
456``after_condition_path`` argument.
457
458Parameters:
459 expression (str): Expression to determine if the conditional path should be executed.
460 This should be one of the comparison operators ``<``, ``>``, ``<=``,
461 ``>=``, ``==``, or ``!=`` followed by a numerical value for the return value
462 condition_path (Path): path to execute in case the expression is fulfilled
463 after_condition_path (AfterConditionPath): What to do once the ``condition_path`` has been executed.
464)DOCSTRING")
465 .def("if_false", &Module::if_false,
466 (bp::arg("condition_path"), bp::arg("after_condition_path")= Module::EAfterConditionPath::c_End),
467 R"DOC(if_false(condition_path, after_condition_path=AfterConditionPath.END)
468
469Sets a conditional sub path which will be executed after this module if
470the return value of the module evaluates to False. This is equivalent to
471calling `if_value` with ``expression=\"<1\"``)DOC")
472 .def("if_true", &Module::if_true,
473 (bp::arg("condition_path"), bp::arg("after_condition_path")= Module::EAfterConditionPath::c_End),
474 R"DOC(if_true(condition_path, after_condition_path=AfterConditionPath.END)
475
476Sets a conditional sub path which will be executed after this module if
477the return value of the module evaluates to True. It is equivalent to
478calling `if_value` with ``expression=\">=1\"``)DOC");
479 }
480 module
481 .def("has_condition", &Module::hasCondition,
482 "Return true if a conditional path has been set for this module "
483 "using `if_value`, `if_true` or `if_false`")
484 .def("get_all_condition_paths", &_getAllConditionPathsPython,
485 "Return a list of all conditional paths set for this module using "
486 "`if_value`, `if_true` or `if_false`")
487 .def("get_all_conditions", &_getAllConditionsPython,
488 "Return a list of all conditional path expressions set for this module using "
489 "`if_value`, `if_true` or `if_false`")
490 .add_property("logging", make_function(&Module::getLogConfig, return_value_policy<reference_existing_object>()),
@ c_GE
Greater or equal than: ">=".
@ c_SE
Smaller or equal than: "<=".
@ c_GT
Greater than: ">"
@ c_NE
Not equal: "!=".
@ c_EQ
Equal: "=" or "=="
@ c_ST
Smaller than: "<"
Base class for Modules.
Definition: Module.h:72
LogConfig & getLogConfig()
Returns the log system configuration.
Definition: Module.h:225
void if_value(const std::string &expression, const std::shared_ptr< Path > &path, EAfterConditionPath afterConditionPath=EAfterConditionPath::c_End)
Add a condition to the module.
Definition: Module.cc:79
void if_true(const std::shared_ptr< Path > &path, EAfterConditionPath afterConditionPath=EAfterConditionPath::c_End)
A simplified version to set the condition of the module.
Definition: Module.cc:90
void setReturnValue(int value)
Sets the return value for this module as integer.
Definition: Module.cc:220
void setLogConfig(const LogConfig &logConfig)
Set the log system configuration.
Definition: Module.h:230
const std::string & getDescription() const
Returns the description of the module.
Definition: Module.h:202
void if_false(const std::shared_ptr< Path > &path, EAfterConditionPath afterConditionPath=EAfterConditionPath::c_End)
A simplified version to add a condition to the module.
Definition: Module.cc:85
bool hasCondition() const
Returns true if at least one condition was set for the module.
Definition: Module.h:311
const std::string & getPackage() const
Returns the package this module is in.
Definition: Module.h:197
void setName(const std::string &name)
Set the name of the module.
Definition: Module.h:214
bool hasProperties(unsigned int propertyFlags) const
Returns true if all specified property flags are available in this module.
Definition: Module.cc:160
std::string getPathString() const override
return the module name.
Definition: Module.cc:192

◆ getAfterConditionPath()

Module::EAfterConditionPath getAfterConditionPath ( ) const
inherited

What to do after the conditional path is finished.

(defaults to c_End if no condition is set)

Definition at line 133 of file Module.cc.

134{
135 if (m_conditions.empty()) return EAfterConditionPath::c_End;
136
137 //okay, a condition was set for this Module...
138 if (!m_hasReturnValue) {
139 B2FATAL("A condition was set for '" << getName() << "', but the module did not set a return value!");
140 }
141
142 for (const auto& condition : m_conditions) {
143 if (condition.evaluate(m_returnValue)) {
144 return condition.getAfterConditionPath();
145 }
146 }
147
148 return EAfterConditionPath::c_End;
149}

◆ getAllConditionPaths()

std::vector< std::shared_ptr< Path > > getAllConditionPaths ( ) const
inherited

Return all condition paths currently set (no matter if the condition is true or not).

Definition at line 150 of file Module.cc.

151{
152 std::vector<std::shared_ptr<Path>> allConditionPaths;
153 for (const auto& condition : m_conditions) {
154 allConditionPaths.push_back(condition.getPath());
155 }
156
157 return allConditionPaths;
158}

◆ getAllConditions()

const std::vector< ModuleCondition > & getAllConditions ( ) const
inlineinherited

Return all set conditions for this module.

Definition at line 324 of file Module.h.

325 {
326 return m_conditions;
327 }

◆ getCondition()

const ModuleCondition * getCondition ( ) const
inlineinherited

Return a pointer to the first condition (or nullptr, if none was set)

Definition at line 314 of file Module.h.

315 {
316 if (m_conditions.empty()) {
317 return nullptr;
318 } else {
319 return &m_conditions.front();
320 }
321 }

◆ getConditionPath()

std::shared_ptr< Path > getConditionPath ( ) const
inherited

Returns the path of the last true condition (if there is at least one, else reaturn a null pointer).


Definition at line 113 of file Module.cc.

114{
115 PathPtr p;
116 if (m_conditions.empty()) return p;
117
118 //okay, a condition was set for this Module...
119 if (!m_hasReturnValue) {
120 B2FATAL("A condition was set for '" << getName() << "', but the module did not set a return value!");
121 }
122
123 for (const auto& condition : m_conditions) {
124 if (condition.evaluate(m_returnValue)) {
125 return condition.getPath();
126 }
127 }
128
129 // if none of the conditions were true, return a null pointer.
130 return p;
131}
std::shared_ptr< Path > PathPtr
Defines a pointer to a path object as a boost shared pointer.
Definition: Path.h:35

◆ getDescription()

const std::string & getDescription ( ) const
inlineinherited

Returns the description of the module.

Definition at line 202 of file Module.h.

202{return m_description;}
std::string m_description
The description of the module.
Definition: Module.h:511

◆ getFileNames()

virtual std::vector< std::string > getFileNames ( bool  outputFiles)
inlinevirtualinherited

Return a list of output filenames for this modules.

This will be called when basf2 is run with "--dry-run" if the module has set either the c_Input or c_Output properties.

If the parameter outputFiles is false (for modules with c_Input) the list of input filenames should be returned (if any). If outputFiles is true (for modules with c_Output) the list of output files should be returned (if any).

If a module has sat both properties this member is called twice, once for each property.

The module should return the actual list of requested input or produced output filenames (including handling of input/output overrides) so that the grid system can handle input/output files correctly.

This function should return the same value when called multiple times. This is especially important when taking the input/output overrides from Environment as they get consumed when obtained so the finalized list of output files should be stored for subsequent calls.

Reimplemented in RootInputModule, StorageRootOutputModule, and RootOutputModule.

Definition at line 134 of file Module.h.

135 {
136 return std::vector<std::string>();
137 }

◆ getLogConfig()

LogConfig & getLogConfig ( )
inlineinherited

Returns the log system configuration.

Definition at line 225 of file Module.h.

225{return m_logConfig;}

◆ getModules()

std::list< ModulePtr > getModules ( ) const
inlineoverrideprivatevirtualinherited

no submodules, return empty list

Implements PathElement.

Definition at line 506 of file Module.h.

506{ return std::list<ModulePtr>(); }

◆ getName()

const std::string & getName ( ) const
inlineinherited

Returns the name of the module.

This can be changed via e.g. set_name() in the steering file to give more useful names if there is more than one module of the same type.

For identifying the type of a module, using getType() (or type() in Python) is recommended.

Definition at line 187 of file Module.h.

187{return m_name;}
std::string m_name
The name of the module, saved as a string (user-modifiable)
Definition: Module.h:508

◆ getPackage()

const std::string & getPackage ( ) const
inlineinherited

Returns the package this module is in.

Definition at line 197 of file Module.h.

197{return m_package;}

◆ getParamInfoListPython()

std::shared_ptr< boost::python::list > getParamInfoListPython ( ) const
inherited

Returns a python list of all parameters.

Each item in the list consists of the name of the parameter, a string describing its type, a python list of all default values and the description of the parameter.

Returns
A python list containing the parameters of this parameter list.

Definition at line 279 of file Module.cc.

280{
282}
std::shared_ptr< boost::python::list > getParamInfoListPython() const
Returns a python list of all parameters.
ModuleParamList m_moduleParamList
List storing and managing all parameter of the module.
Definition: Module.h:516

◆ getParamList()

const ModuleParamList & getParamList ( ) const
inlineinherited

Return module param list.

Definition at line 363 of file Module.h.

363{ return m_moduleParamList; }

◆ getPathString()

std::string getPathString ( ) const
overrideprivatevirtualinherited

return the module name.

Implements PathElement.

Definition at line 192 of file Module.cc.

193{
194
195 std::string output = getName();
196
197 for (const auto& condition : m_conditions) {
198 output += condition.getString();
199 }
200
201 return output;
202}

◆ getReturnValue()

int getReturnValue ( ) const
inlineinherited

Return the return value set by this module.

This value is only meaningful if hasReturnValue() is true

Definition at line 381 of file Module.h.

381{ return m_returnValue; }

◆ getType()

const std::string & getType ( ) const
inherited

Returns the type of the module (i.e.

class name minus 'Module')

Definition at line 41 of file Module.cc.

42{
43 if (m_type.empty())
44 B2FATAL("Module type not set for " << getName());
45 return m_type;
46}
std::string m_type
The type of the module, saved as a string.
Definition: Module.h:509

◆ hasCondition()

bool hasCondition ( ) const
inlineinherited

Returns true if at least one condition was set for the module.

Definition at line 311 of file Module.h.

311{ return not m_conditions.empty(); };

◆ hasProperties()

bool hasProperties ( unsigned int  propertyFlags) const
inherited

Returns true if all specified property flags are available in this module.

Parameters
propertyFlagsOred EModulePropFlags which should be compared with the module flags.

Definition at line 160 of file Module.cc.

161{
162 return (propertyFlags & m_propertyFlags) == propertyFlags;
163}

◆ hasReturnValue()

bool hasReturnValue ( ) const
inlineinherited

Return true if this module has a valid return value set.

Definition at line 378 of file Module.h.

378{ return m_hasReturnValue; }

◆ hasUnsetForcedParams()

bool hasUnsetForcedParams ( ) const
inherited

Returns true and prints error message if the module has unset parameters which the user has to set in the steering file.

Definition at line 166 of file Module.cc.

167{
169 std::string allMissing = "";
170 for (const auto& s : missing)
171 allMissing += s + " ";
172 if (!missing.empty())
173 B2ERROR("The following required parameters of Module '" << getName() << "' were not specified: " << allMissing <<
174 "\nPlease add them to your steering file.");
175 return !missing.empty();
176}
std::vector< std::string > getUnsetForcedParams() const
Returns list of unset parameters (if they are required to have a value.

◆ if_false()

void if_false ( const std::shared_ptr< Path > &  path,
EAfterConditionPath  afterConditionPath = EAfterConditionPath::c_End 
)
inherited

A simplified version to add a condition to the module.

Please note that successive calls of this function will add more than one condition to the module. If more than one condition results in true, only the last of them will be used.

Please be careful: Avoid creating cyclic paths, e.g. by linking a condition to a path which is processed before the path where this module is located in.

It is equivalent to the if_value() method, using the expression "<1". This method is meant to be used together with the setReturnValue(bool value) method.

Parameters
pathShared pointer to the Path which will be executed if the return value is false.
afterConditionPathWhat to do after executing 'path'.

Definition at line 85 of file Module.cc.

86{
87 if_value("<1", path, afterConditionPath);
88}

◆ if_true()

void if_true ( const std::shared_ptr< Path > &  path,
EAfterConditionPath  afterConditionPath = EAfterConditionPath::c_End 
)
inherited

A simplified version to set the condition of the module.

Please note that successive calls of this function will add more than one condition to the module. If more than one condition results in true, only the last of them will be used.

Please be careful: Avoid creating cyclic paths, e.g. by linking a condition to a path which is processed before the path where this module is located in.

It is equivalent to the if_value() method, using the expression ">=1". This method is meant to be used together with the setReturnValue(bool value) method.

Parameters
pathShared pointer to the Path which will be executed if the return value is true.
afterConditionPathWhat to do after executing 'path'.

Definition at line 90 of file Module.cc.

91{
92 if_value(">=1", path, afterConditionPath);
93}

◆ if_value()

void if_value ( const std::string &  expression,
const std::shared_ptr< Path > &  path,
EAfterConditionPath  afterConditionPath = EAfterConditionPath::c_End 
)
inherited

Add a condition to the module.

Please note that successive calls of this function will add more than one condition to the module. If more than one condition results in true, only the last of them will be used.

See https://xwiki.desy.de/xwiki/rest/p/a94f2 or ModuleCondition for a description of the syntax.

Please be careful: Avoid creating cyclic paths, e.g. by linking a condition to a path which is processed before the path where this module is located in.

Parameters
expressionThe expression of the condition.
pathShared pointer to the Path which will be executed if the condition is evaluated to true.
afterConditionPathWhat to do after executing 'path'.

Definition at line 79 of file Module.cc.

80{
81 m_conditions.emplace_back(expression, path, afterConditionPath);
82}

◆ initialize()

void initialize ( void  )
overridevirtual

Initilizes TRGECLUnpackerModuel.

Reimplemented from Module.

Definition at line 42 of file trgeclUnpackerModule.cc.

43{
44 m_TRGECLSumArray.registerInDataStore();
45 m_TRGECLTCArray.registerInDataStore();
46 m_TRGECLEvtArray.registerInDataStore();
47 m_TRGECLClusterArray.registerInDataStore();
48 // This object is registered by few packages. Let's be agnostic about the
49 // execution order of the modules: the first package run registers the module
50 m_eventLevelClusteringInfo.isOptional() ? m_eventLevelClusteringInfo.isRequired() :
51 m_eventLevelClusteringInfo.registerInDataStore();
52}

◆ readCOPPEREvent()

void readCOPPEREvent ( RawTRG raw_copper,
int  i,
int  nnn,
int  ch 
)
virtual

Read data from TRG copper.

Definition at line 107 of file trgeclUnpackerModule.cc.

108{
109 /* cppcheck-suppress variableScope */
110 int* rdat;
111 if (raw_copper->GetDetectorNwords(i, ch) > 0) {
112 rdat = raw_copper->GetDetectorBuffer(i, ch);
113 etm_version = ((rdat[0] >> 16) & 0xffff);
114 if (etm_version > 136) {
115 checkBuffer(rdat, nnn);
116 } else {
117 checkBuffer_v136(rdat, nnn);
118 }
119 }
120}
virtual void checkBuffer_v136(int *, int)
Unpacker main function for upto version 136.
virtual void checkBuffer(int *, int)
Unpacker main function.
int GetDetectorNwords(int n, int finesse_num)
get Detector buffer length
Definition: RawCOPPER.h:657
int * GetDetectorBuffer(int n, int finesse_num)
get Detector buffer
Definition: RawCOPPER.h:681

◆ setAbortLevel()

void setAbortLevel ( int  abortLevel)
inherited

Configure the abort log level.

Definition at line 67 of file Module.cc.

68{
69 m_logConfig.setAbortLevel(static_cast<LogConfig::ELogLevel>(abortLevel));
70}
ELogLevel
Definition of the supported log levels.
Definition: LogConfig.h:26
void setAbortLevel(ELogLevel abortLevel)
Configure the abort level.
Definition: LogConfig.h:112

◆ setDebugLevel()

void setDebugLevel ( int  debugLevel)
inherited

Configure the debug messaging level.

Definition at line 61 of file Module.cc.

62{
63 m_logConfig.setDebugLevel(debugLevel);
64}
void setDebugLevel(int debugLevel)
Configure the debug messaging level.
Definition: LogConfig.h:98

◆ setDescription()

void setDescription ( const std::string &  description)
protectedinherited

Sets the description of the module.

Parameters
descriptionA description of the module.

Definition at line 214 of file Module.cc.

215{
216 m_description = description;
217}

◆ setLogConfig()

void setLogConfig ( const LogConfig logConfig)
inlineinherited

Set the log system configuration.

Definition at line 230 of file Module.h.

230{m_logConfig = logConfig;}

◆ setLogInfo()

void setLogInfo ( int  logLevel,
unsigned int  logInfo 
)
inherited

Configure the printed log information for the given level.

Parameters
logLevelThe log level (one of LogConfig::ELogLevel)
logInfoWhat kind of info should be printed? ORed combination of LogConfig::ELogInfo flags.

Definition at line 73 of file Module.cc.

74{
75 m_logConfig.setLogInfo(static_cast<LogConfig::ELogLevel>(logLevel), logInfo);
76}
void setLogInfo(ELogLevel logLevel, unsigned int logInfo)
Configure the printed log information for the given level.
Definition: LogConfig.h:127

◆ setLogLevel()

void setLogLevel ( int  logLevel)
inherited

Configure the log level.

Definition at line 55 of file Module.cc.

56{
57 m_logConfig.setLogLevel(static_cast<LogConfig::ELogLevel>(logLevel));
58}
void setLogLevel(ELogLevel logLevel)
Configure the log level.
Definition: LogConfig.cc:25

◆ setName()

void setName ( const std::string &  name)
inlineinherited

Set the name of the module.

Note
The module name is set when using the REG_MODULE macro, but the module can be renamed before calling process() using the set_name() function in your steering file.
Parameters
nameThe name of the module

Definition at line 214 of file Module.h.

214{ m_name = name; };

◆ setParamList()

void setParamList ( const ModuleParamList params)
inlineprotectedinherited

Replace existing parameter list.

Definition at line 501 of file Module.h.

501{ m_moduleParamList = params; }

◆ setParamPython()

void setParamPython ( const std::string &  name,
const boost::python::object &  pyObj 
)
privateinherited

Implements a method for setting boost::python objects.

The method supports the following types: list, dict, int, double, string, bool The conversion of the python object to the C++ type and the final storage of the parameter value is done in the ModuleParam class.

Parameters
nameThe unique name of the parameter.
pyObjThe object which should be converted and stored as the parameter value.

Definition at line 234 of file Module.cc.

235{
236 LogSystem& logSystem = LogSystem::Instance();
237 logSystem.updateModule(&(getLogConfig()), getName());
238 try {
240 } catch (std::runtime_error& e) {
241 throw std::runtime_error("Cannot set parameter '" + name + "' for module '"
242 + m_name + "': " + e.what());
243 }
244
245 logSystem.updateModule(nullptr);
246}
Class for logging debug, info and error messages.
Definition: LogSystem.h:46
void updateModule(const LogConfig *moduleLogConfig=nullptr, const std::string &moduleName="")
Sets the log configuration to the given module log configuration and sets the module name This method...
Definition: LogSystem.h:191
static LogSystem & Instance()
Static method to get a reference to the LogSystem instance.
Definition: LogSystem.cc:31
void setParamPython(const std::string &name, const PythonObject &pyObj)
Implements a method for setting boost::python objects.

◆ setParamPythonDict()

void setParamPythonDict ( const boost::python::dict &  dictionary)
privateinherited

Implements a method for reading the parameter values from a boost::python dictionary.

The key of the dictionary has to be the name of the parameter and the value has to be of one of the supported parameter types.

Parameters
dictionaryThe python dictionary from which the parameter values are read.

Definition at line 249 of file Module.cc.

250{
251
252 LogSystem& logSystem = LogSystem::Instance();
253 logSystem.updateModule(&(getLogConfig()), getName());
254
255 boost::python::list dictKeys = dictionary.keys();
256 int nKey = boost::python::len(dictKeys);
257
258 //Loop over all keys in the dictionary
259 for (int iKey = 0; iKey < nKey; ++iKey) {
260 boost::python::object currKey = dictKeys[iKey];
261 boost::python::extract<std::string> keyProxy(currKey);
262
263 if (keyProxy.check()) {
264 const boost::python::object& currValue = dictionary[currKey];
265 setParamPython(keyProxy, currValue);
266 } else {
267 B2ERROR("Setting the module parameters from a python dictionary: invalid key in dictionary!");
268 }
269 }
270
271 logSystem.updateModule(nullptr);
272}
void setParamPython(const std::string &name, const boost::python::object &pyObj)
Implements a method for setting boost::python objects.
Definition: Module.cc:234

◆ setPropertyFlags()

void setPropertyFlags ( unsigned int  propertyFlags)
inherited

Sets the flags for the module properties.

Parameters
propertyFlagsbitwise OR of EModulePropFlags

Definition at line 208 of file Module.cc.

209{
210 m_propertyFlags = propertyFlags;
211}

◆ setReturnValue() [1/2]

void setReturnValue ( bool  value)
protectedinherited

Sets the return value for this module as bool.

The bool value is saved as an integer with the convention 1 meaning true and 0 meaning false. The value can be used in the steering file to divide the analysis chain into several paths.

Parameters
valueThe value of the return value.

Definition at line 227 of file Module.cc.

228{
229 m_hasReturnValue = true;
230 m_returnValue = value;
231}

◆ setReturnValue() [2/2]

void setReturnValue ( int  value)
protectedinherited

Sets the return value for this module as integer.

The value can be used in the steering file to divide the analysis chain into several paths.

Parameters
valueThe value of the return value.

Definition at line 220 of file Module.cc.

221{
222 m_hasReturnValue = true;
223 m_returnValue = value;
224}

◆ setType()

void setType ( const std::string &  type)
protectedinherited

Set the module type.

Only for use by internal modules (which don't use the normal REG_MODULE mechanism).

Definition at line 48 of file Module.cc.

49{
50 if (!m_type.empty())
51 B2FATAL("Trying to change module type from " << m_type << " is not allowed, the value is assumed to be fixed.");
52 m_type = type;
53}

◆ terminate()

void terminate ( void  )
overridevirtual

Called when processing ended.

Reimplemented from Module.

Definition at line 38 of file trgeclUnpackerModule.cc.

39{
40}

◆ version()

string version ( ) const

returns version of TRGECLUnpackerModule.

Definition at line 20 of file trgeclUnpackerModule.cc.

21{
22 return string("4.02");
23}

Member Data Documentation

◆ etm_version

int etm_version = 0
protected

ETM Version.

Definition at line 74 of file trgeclUnpackerModule.h.

◆ iFiness

int iFiness = 0
protected

Finess.

Definition at line 80 of file trgeclUnpackerModule.h.

◆ m_conditions

std::vector<ModuleCondition> m_conditions
privateinherited

Module condition, only non-null if set.

Definition at line 521 of file Module.h.

◆ m_description

std::string m_description
privateinherited

The description of the module.

Definition at line 511 of file Module.h.

◆ m_eventLevelClusteringInfo

StoreObjPtr<EventLevelClusteringInfo> m_eventLevelClusteringInfo
private

EventLevelClusteringInfo.

Definition at line 95 of file trgeclUnpackerModule.h.

◆ m_hasReturnValue

bool m_hasReturnValue
privateinherited

True, if the return value is set.

Definition at line 518 of file Module.h.

◆ m_logConfig

LogConfig m_logConfig
privateinherited

The log system configuration of the module.

Definition at line 514 of file Module.h.

◆ m_moduleParamList

ModuleParamList m_moduleParamList
privateinherited

List storing and managing all parameter of the module.

Definition at line 516 of file Module.h.

◆ m_name

std::string m_name
privateinherited

The name of the module, saved as a string (user-modifiable)

Definition at line 508 of file Module.h.

◆ m_package

std::string m_package
privateinherited

Package this module is found in (may be empty).

Definition at line 510 of file Module.h.

◆ m_propertyFlags

unsigned int m_propertyFlags
privateinherited

The properties of the module as bitwise or (with |) of EModulePropFlags.

Definition at line 512 of file Module.h.

◆ m_returnValue

int m_returnValue
privateinherited

The return value.

Definition at line 519 of file Module.h.

◆ m_TRGECLClusterArray

StoreArray<TRGECLCluster> m_TRGECLClusterArray
private

ECL Trigger Cluster output.

Definition at line 92 of file trgeclUnpackerModule.h.

◆ m_TRGECLEvtArray

StoreArray<TRGECLUnpackerEvtStore> m_TRGECLEvtArray
private

ECL Trigger Unpacker Event output.

Definition at line 90 of file trgeclUnpackerModule.h.

◆ m_TRGECLSumArray

StoreArray<TRGECLUnpackerSumStore> m_TRGECLSumArray
private

ECL Trigger Unpacker Summary output.

Definition at line 88 of file trgeclUnpackerModule.h.

◆ m_TRGECLTCArray

StoreArray<TRGECLUnpackerStore> m_TRGECLTCArray
private

ECL Trigger Unpacker TC output.

Definition at line 86 of file trgeclUnpackerModule.h.

◆ m_type

std::string m_type
privateinherited

The type of the module, saved as a string.

Definition at line 509 of file Module.h.

◆ n_basf2evt

int n_basf2evt
protected

Event number.

Definition at line 72 of file trgeclUnpackerModule.h.

◆ nodeid

unsigned int nodeid = 0
protected

Node Id.

Definition at line 76 of file trgeclUnpackerModule.h.

◆ nwords

int nwords = 0
protected

N Word.

Definition at line 78 of file trgeclUnpackerModule.h.

◆ trgtype

int trgtype = 0
protected

Trigger Type.

Definition at line 82 of file trgeclUnpackerModule.h.


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