REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestEventProcess.cxx
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see http://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see http://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
52
53#include "TRestEventProcess.h"
54
55#include <TBuffer.h>
56#include <TClass.h>
57#include <TDataMember.h>
58
59#include "TRestManager.h"
60#include "TRestRun.h"
61
62using namespace std;
63
64ClassImp(TRestEventProcess);
65
70 fObservablesDefined.clear();
71 fObservablesUpdated.clear();
72 fSingleThreadOnly = false;
73 fCanvasSize = TVector2(800, 600);
74}
75
80
99vector<string> TRestEventProcess::ReadObservables() {
100 TiXmlElement* e = GetElement("observable");
101 vector<string> obsNames;
102 vector<string> obsTypes;
103 vector<string> obsDesc;
104
105 while (e != nullptr) {
106 const char* obschr = e->Attribute("name");
107 const char* _value = e->Attribute("value");
108 const char* _type = e->Attribute("type");
109 const char* _desc = e->Attribute("description");
110
111 string value;
112 if (_value == nullptr)
113 value = "ON";
114 else
115 value = _value;
116 string type;
117 if (_type == nullptr)
118 type = "double";
119 else
120 type = _type;
121 string description;
122 if (_desc == nullptr)
123 description = "";
124 else
125 description = _desc;
126
127 if (ToUpper(value) == "ON") {
128 if (obschr != nullptr) {
129 RESTDebug << this->ClassName() << " : setting observable \"" << obschr << "\"" << RESTendl;
130 // vector<string> tmp = Split(obsstring, ":");
131 obsNames.emplace_back(obschr);
132 obsTypes.push_back(type);
133 obsDesc.push_back(description);
134 }
135 }
136
137 e = e->NextSiblingElement("observable");
138
139 } // now we get a list of all observal names, we add them into fAnalysisTree and fObservableInfo
140
141 for (unsigned int i = 0; i < obsNames.size(); i++) {
142 string obsName = this->GetName() + (string) "_" + obsNames[i];
143 fAnalysisTree->AddObservable(obsName, obsTypes[i], obsDesc[i]);
144 int id = fAnalysisTree->GetObservableID(obsName);
145 if (id != -1) {
146 fObservablesDefined[(string)GetName() + "_" + obsNames[i]] = id;
147 }
148 }
149
150 return obsNames;
151}
152
154 fAnalysisTree = tree;
155 if (fAnalysisTree == nullptr) return;
156 ReadObservables();
157}
158
164 if (p == nullptr) return;
165 for (const auto& fFriendlyProcess : fFriendlyProcesses) {
166 if (fFriendlyProcess->GetName() == p->GetName()) {
167 return;
168 }
169 }
170 fFriendlyProcesses.push_back(p);
171}
172
175 if (p == nullptr) return;
176 for (const auto& fParallelProcess : fParallelProcesses) {
177 if (fParallelProcess == p || this == p) {
178 return;
179 }
180 }
181 fParallelProcesses.push_back(p);
182}
183
186Bool_t TRestEventProcess::OpenInputFiles(const vector<string>& files) { return false; }
187
194
195 if (ToUpper(GetParameter("observable", "")) == "ALL" ||
196 ToUpper(GetParameter("observables", "")) == "ALL") {
197 fDynamicObs = true;
198 }
199
200 // load cuts
201 fCuts.clear();
202 if (ToUpper(GetParameter("cutsEnabled", "false")) == "TRUE") {
203 TiXmlElement* ele = fElement->FirstChildElement();
204 while (ele != nullptr) {
205 if (ele->Value() != nullptr && (string)ele->Value() == "cut") {
206 if (ele->Attribute("name") != nullptr && ele->Attribute("value") != nullptr) {
207 string name = ele->Attribute("name");
208 name = (string)this->GetName() + "_" + name;
209 TVector2 value = StringTo2DVector(ele->Attribute("value"));
210 if (value.X() != value.Y()) fCuts.push_back(pair<string, TVector2>(name, value));
211 }
212 }
213
214 else if (ele->Value() != nullptr && (string)ele->Value() == "parameter") {
215 if (ele->Attribute("name") != nullptr && ele->Attribute("value") != nullptr) {
216 string name = ele->Attribute("name");
217 if (name.find("Cut") == name.size() - 3 || name.find("CutRange") == name.size() - 8) {
218 name = name.substr(0, name.find("Cut") + 3);
219 TVector2 value = StringTo2DVector(ele->Attribute("value"));
220 if (value.X() != value.Y()) fCuts.push_back(pair<string, TVector2>(name, value));
221 }
222 }
223 }
224
225 ele = ele->NextSiblingElement();
226 }
227 }
228
229 return 0;
230}
231
238 TRestMetadata* m = nullptr;
239 if (fRunInfo != nullptr) {
240 m = fRunInfo->GetMetadata(name);
241 if (m == nullptr) m = fRunInfo->GetMetadataClass(name);
242 }
243 if (fHostmgr != nullptr) {
244 if (m == nullptr) m = fHostmgr->GetMetadata(name);
245 if (m == nullptr) m = fHostmgr->GetMetadataClass(name);
246 }
247 return m;
248}
249
262TRestEventProcess* TRestEventProcess::GetFriend(const string& nameOrType) {
263 TRestEventProcess* proc = GetFriendLive(nameOrType);
264 if (proc == nullptr) {
265 TRestMetadata* friendfromfile = GetMetadata(nameOrType);
266 if (friendfromfile != nullptr && friendfromfile->InheritsFrom("TRestEventProcess")) {
267 return (TRestEventProcess*)friendfromfile;
268 }
269 return nullptr;
270 } else {
271 return proc;
272 }
273}
274
279TRestEventProcess* TRestEventProcess::GetFriendLive(const string& nameOrType) {
280 for (const auto& fFriendlyProcess : fFriendlyProcesses) {
281 if ((string)fFriendlyProcess->GetName() == nameOrType ||
282 (string)fFriendlyProcess->ClassName() == nameOrType) {
283 return fFriendlyProcess;
284 }
285 }
286 return nullptr;
287}
288
298TRestEventProcess* TRestEventProcess::GetParallel(int i) {
299 if (i >= 0 && i < (int)fParallelProcesses.size()) {
300 return fParallelProcesses[i];
301 }
302 return nullptr;
303}
304
309bool TRestEventProcess::ApplyCut() {
310 for (auto cut : fCuts) {
311 string type = (string)fAnalysisTree->GetObservableType(cut.first);
312 if (fAnalysisTree != nullptr && type == "double") {
313 double val = fAnalysisTree->GetObservableValue<double>(cut.first);
314 if (val > cut.second.Y() || val < cut.second.X()) {
315 return true;
316 }
317 }
318 if (fAnalysisTree != nullptr && type == "int") {
319 int val = fAnalysisTree->GetObservableValue<int>(cut.first);
320 if (val > cut.second.Y() || val < cut.second.X()) {
321 return true;
322 }
323 }
324 }
325 return false;
326}
327
328/*
329
330void TRestEventProcess::InitProcess()
331{
332// virtual function to be executed once at the beginning of process
333// (before starting the process of the events)
334cout << GetName() << ": Process initialization..." << endl;
335}
336*/
337
347 RESTDebug << "Entering " << ClassName() << "::BeginOfEventProcess, Initializing output event..."
348 << RESTendl;
349 if (inEv != nullptr && GetOutputEvent().address != nullptr && (TRestEvent*)GetOutputEvent() != inEv) {
350 TRestEvent* outEv = GetOutputEvent();
351 outEv->Initialize();
352
353 outEv->SetID(inEv->GetID());
354 outEv->SetSubID(inEv->GetSubID());
355 outEv->SetSubEventTag(inEv->GetSubEventTag());
356
357 outEv->SetRunOrigin(inEv->GetRunOrigin());
358 outEv->SetSubRunOrigin(inEv->GetSubRunOrigin());
359
360 outEv->SetTime(inEv->GetTime());
361 }
362
363 fObservablesUpdated.clear();
364
365 // TODO if fIsExternal and we already have defined the fAnalysisTree run#,
366 // evId#, timestamp, etc at the analysisTree we could stamp the output event
367 // here.
368}
369
370/*
371
372void TRestEventProcess::ProcessEvent( TRestEvent *inputEvent )
373{
374// virtual function to be executed for every event to be processed
375}
376*/
377
382 RESTDebug << "Entering TRestEventProcess::EndOfEventProcess (" << ClassName() << ")" << RESTendl;
384 if (fObservablesDefined.size() != fObservablesUpdated.size()) {
385 for (auto x : fObservablesDefined) {
386 if (fObservablesUpdated.count(x.first) == 0) {
387 // the observable is added through <observable section but not set in the process
388 RESTWarning
389 << "The observable '" << x.first << "' is defined but not set by "
390 << this->ClassName()
391 << ", the event is empty? Makesure all observables are set through ProcessEvent()"
392 << RESTendl;
393 }
394 }
395 }
396 }
397}
398
399/*
400
401void TRestEventProcess::EndProcess()
402{
403// virtual function to be executed once at the end of the process
404// (after all events have been processed)
405cout << GetName() << ": Process ending..." << endl;
406}
407*/
408
416 RESTMetadata.setcolor(COLOR_BOLDGREEN);
417 RESTMetadata.setborder("||");
418 RESTMetadata.setlength(100);
419 // metadata << " " << endl;
420 cout << endl;
421 RESTMetadata << "=" << RESTendl;
422 RESTMetadata << "Process : " << ClassName() << RESTendl;
423 RESTMetadata << "Name: " << GetName() << " Title: " << GetTitle()
424 << " VerboseLevel: " << GetVerboseLevelString() << RESTendl;
425 RESTMetadata << " ----------------------------------------------- " << RESTendl;
426 RESTMetadata << " " << RESTendl;
427
428 if (!fObservablesDefined.empty()) {
429 RESTMetadata << " Analysis tree observables added by this process " << RESTendl;
430 RESTMetadata << " +++++++++++++++++++++++++++++++++++++++++++++++ " << RESTendl;
431 }
432
433 auto iter = fObservablesDefined.begin();
434 while (iter != fObservablesDefined.end()) {
435 RESTMetadata << " ++ " << iter->first << RESTendl;
436 iter++;
437 }
438
439 if (!fObservablesDefined.empty()) {
440 RESTMetadata << " +++++++++++++++++++++++++++++++++++++++++++++++ " << RESTendl;
441 RESTMetadata << " " << RESTendl;
442 }
443}
444
450void TRestEventProcess::EndPrintProcess() {
451 if (fCuts.size() > 0) {
452 RESTMetadata << "Cuts enabled" << RESTendl;
453 RESTMetadata << "------------" << RESTendl;
454
455 auto iter = fCuts.begin();
456 while (iter != fCuts.end()) {
457 if (iter->second.X() != iter->second.Y())
458 RESTMetadata << iter->first << ", range : ( " << iter->second.X() << " , " << iter->second.Y()
459 << " ) " << RESTendl;
460 iter++;
461 }
462 }
463
464 RESTMetadata << " " << RESTendl;
465 RESTMetadata << "=" << RESTendl;
466 RESTMetadata << RESTendl;
467 RESTMetadata.resetcolor();
468 RESTMetadata.setborder("");
469 RESTMetadata.setlength(10000);
470}
471
479// Double_t TRestEventProcess::GetDoubleParameterFromFriends(string className, string parName) {
480// for (size_t i = 0; i < fFriendlyProcesses.size(); i++)
481// if ((string)fFriendlyProcesses[i]->ClassName() == (string)className)
482// return StringToDouble(fFriendlyProcesses[i]->GetParameter((string)parName));
483//
484// return PARAMETER_NOT_FOUND_DBL;
485//}
486//
494// Double_t TRestEventProcess::GetDoubleParameterFromFriendsWithUnits(string className, string parName) {
495// for (size_t i = 0; i < fFriendlyProcesses.size(); i++)
496// if ((string)fFriendlyProcesses[i]->ClassName() == (string)className)
497// return fFriendlyProcesses[i]->GetDblParameterWithUnits((string)parName);
498//
499// return PARAMETER_NOT_FOUND_DBL;
500//}
501
504TRestAnalysisTree* TRestEventProcess::GetFullAnalysisTree() {
505 if (fHostmgr != nullptr && fHostmgr->GetProcessRunner() != nullptr)
506 return fHostmgr->GetProcessRunner()->GetOutputAnalysisTree();
507 return nullptr;
508}
509
512std::vector<string> TRestEventProcess::GetListOfAddedObservables() {
513 vector<string> list;
514 auto iter = fObservablesDefined.begin();
515 while (iter != fObservablesDefined.end()) {
516 list.push_back(iter->first);
517 iter++;
518 }
519 return list;
520}
REST core data-saving helper based on TTree.
T GetObservableValue(Int_t n)
Get observable in a given type, according to its id.
Int_t GetObservableID(const std::string &obsName)
Get the index of the specified observable.
TString GetObservableType(Int_t n)
Get the observable type according to id.
A base class for any REST event process.
virtual void EndOfEventProcess(TRestEvent *inputEvent=nullptr)
End of event process. Nothing to do. Called directly after ProcessEvent()
std::vector< TRestEventProcess * > fFriendlyProcesses
/// not used
TRestAnalysisTree * fAnalysisTree
TVector2 fCanvasSize
Canvas size.
void BeginPrintProcess()
[name, cut range]
void SetAnalysisTree(TRestAnalysisTree *tree)
Set analysis tree of this process, then add observables to it.
bool fValidateObservables
It defines if observable names should be added to the validation list.
virtual void BeginOfEventProcess(TRestEvent *inputEvent=nullptr)
Begin of event process, preparation work. Called right before ProcessEvent()
void SetParallelProcess(TRestEventProcess *p)
Add parallel process to this process.
std::vector< std::pair< std::string, TVector2 > > fCuts
Stores cut definitions. Any listed observables should be in the range.
TRestRun * fRunInfo
< Pointer to TRestRun object where to find metadata.
Int_t LoadSectionMetadata() override
This method does some preparation of xml section.
bool fDynamicObs
It defines whether to use added observables only or all the observables appear in the code.
std::vector< TRestEventProcess * > fParallelProcesses
Stores a list of parallel processes if multithreading is enabled.
virtual RESTValue GetOutputEvent() const =0
Get pointer to output event. Must be implemented in the derived class.
std::map< std::string, int > fObservablesUpdated
Stores the list of process observables updated when processing this event.
std::map< std::string, int > fObservablesDefined
Stores the list of all the appeared process observables in the code.
T * GetMetadata()
Get a metadata object from the host TRestRun.
void SetFriendProcess(TRestEventProcess *p)
Add friendly process to this process.
A base class for any REST event.
Definition: TRestEvent.h:38
void SetTime(Double_t time)
Definition: TRestEvent.cxx:85
virtual void Initialize()=0
Definition: TRestEvent.cxx:73
TRestMetadata * GetMetadataClass(std::string type)
Get the application metadata class, according to the type.
TRestMetadata * GetMetadata(std::string name)
Get the application metadata class, according to the name.
A base class for any REST metadata class.
Definition: TRestMetadata.h:74
endl_t RESTendl
Termination flag object for TRestStringOutput.
TiXmlElement * GetElement(std::string eleDeclare, TiXmlElement *e=nullptr)
Get an xml element from a given parent element, according to its declaration.
TString GetVerboseLevelString()
returns the verbose level in type of TString
TRestManager * fHostmgr
All metadata classes can be initialized and managed by TRestManager.
virtual Int_t LoadSectionMetadata()
This method does some preparation of xml section.
std::string GetParameter(std::string parName, TiXmlElement *e, TString defaultValue=PARAMETER_NOT_FOUND_STR)
Returns the value for the parameter named parName in the given section.
TiXmlElement * fElement
Saving the sectional element together with global element.
std::string ToUpper(std::string in)
Convert string to its upper case. Alternative of TString::ToUpper.
TVector2 StringTo2DVector(std::string in)
Gets a 2D-vector from a string.