BmnRoot
Loading...
Searching...
No Matches
BmnScWallReconstructor.cxx
Go to the documentation of this file.
1
9
10#include "TTreeFormula.h"
11
12BmnScWallReconstructor::BmnScWallReconstructor(bool isExp, bool isGlobalCoordinates)
13{
14 fIsExp = isExp;
15 fIsGlobal = isGlobalCoordinates;
16 LOG(detail) << Form("BmnScWallReconstructor: isExp %i; isGlobal %i", fIsExp, fIsGlobal);
17}
18
19void BmnScWallReconstructor::SetRecoCutsFile(TString reco_cuts_file)
20{
21 fRecoCutsFile = reco_cuts_file;
22 ParseCuts();
23 LOG(detail) << Form("BmnScWallReconstructor: Reco cuts %s", fRecoCutsFile.Data());
24}
25
27{
28 if (fBmnScWallEvent)
29 delete fBmnScWallEvent;
30}
31
33{
34 fworkTime = 0.;
35 fpFairRootMgr = FairRootManager::Instance();
36 (fIsExp) ? fArrayOfDigits = (TClonesArray*)fpFairRootMgr->GetObject("ScWallDigi")
37 : fArrayOfDigits = (TClonesArray*)fpFairRootMgr->GetObject("ScWallDigit");
38
39 if (!fArrayOfDigits) {
40 LOG(error) << "BmnScWallReconstructor::Init(): branch with Digits not found! Task will be deactivated";
41 SetActive(kFALSE);
42 return kERROR;
43 }
44
45 fBmnScWallEvent = new BmnScWallEvent();
46 fBmnScWallEvent->reset();
48 fpFairRootMgr->RegisterAny("ScWallEvent", fBmnScWallEvent, kTRUE);
49
50 LOG(debug) << "ScWall Reconstructor ready";
51 return kSUCCESS;
52}
53
55{
56 namespace po = boost::program_options;
57
58 TString dir = getenv("VMCWORKDIR");
59 TString path = dir + "/parameters/scwall/";
60
61 float version;
62 std::string comment;
63 std::string cutExp;
64 std::string cutSim;
65
66 // Setup options.
67 po::options_description desc("Options");
68 desc.add_options()("VERSION.id", po::value<float>(&version),
69 "version identificator")("COMMENT.str", po::value<std::string>(&comment), "comment")(
70 "EXPERIMENT.cut", po::value<std::string>(&cutExp),
71 "cut for Experimental data")("SIMULATION.cut", po::value<std::string>(&cutSim), "cut for Simulated data");
72
73 // Load config file.
74 po::variables_map vm;
75 std::ifstream config_file((path + fRecoCutsFile).Data(), std::ifstream::in);
76 if (!config_file.is_open()) {
77 LOG(error) << Form("BmnScWallReconstructor: Loading Cuts from file: %s - file open error!",
78 fRecoCutsFile.Data());
79 return;
80 }
81 LOG(debug) << Form("BmnScWallReconstructor: Loading Cuts from file: %s", fRecoCutsFile.Data());
82 po::store(po::parse_config_file(config_file, desc), vm);
83 config_file.close();
84 po::notify(vm);
85
86 if (fIsExp)
87 fSelectionString = cutExp;
88 else
89 fSelectionString = cutSim;
90}
91
93{
94 BmnScWallGeo* GeoHandler = new BmnScWallGeo();
95 GeoHandler->ReadGeometryFromGeoManager(fIsGlobal);
96 const auto& positionMap = GeoHandler->GetPositionMap();
97 for (auto it : positionMap) {
98 auto unique_address = it.first;
99 if (unique_address == 0)
100 continue;
101 auto cell_id = BmnScWallAddress::GetCellId(unique_address);
102 auto this_cell = fBmnScWallEvent->GetCell(cell_id);
103 this_cell->SetCellId(cell_id);
104 this_cell->SetSignal(0.);
105
106 this_cell->SetPosition(it.second.first);
107 this_cell->SetPositionError(it.second.second);
108 }
109 delete GeoHandler;
110}
111
113{
114
115 if (!IsActive())
116 return;
117
118 TStopwatch sw;
119 sw.Start();
120 fBmnScWallEvent->ResetCells();
121 auto digiVector = GetSelectedDigiVector(fSelectionString);
122
123 for (auto it : digiVector) {
124 BmnScWallDigit* ThisDigi = (BmnScWallDigit*)it;
125 auto cell_id = ThisDigi->GetCellId(); // 1 to 174
126 if (cell_id <= 0 || cell_id > BmnScWallEvent::fgkMaxCells) {
127 if (fVerbose)
128 LOG(warning) << Form("ScWall digi ignored. Cell %d ", cell_id);
129 continue;
130 }
131 fBmnScWallEvent->GetCell(cell_id)->SetSignal(ThisDigi->GetSignal());
132 }
133
134 fBmnScWallEvent->SummarizeEvent();
135 sw.Stop();
136 fworkTime += sw.RealTime();
137}
138
139std::vector<TObject*> BmnScWallReconstructor::GetSelectedDigiVector(TString formulaString)
140{
141 TTree* temp_tree = new TTree("temp_tree", "temp_tree");
142 TString temp_br_name = (fIsExp) ? "ScWallDigi" : "ScWallDigit";
143 fSelectedDigiVector.clear();
144
145 if (fIsExp) {
146 BmnScWallDigi* ptr = nullptr;
147 temp_tree->Branch(temp_br_name.Data(), &ptr);
148 for (int i = 0; i < fArrayOfDigits->GetEntriesFast(); i++) {
149 ptr = (BmnScWallDigi*)fArrayOfDigits->At(i);
150 temp_tree->Fill();
151 }
152 } else {
153 BmnScWallDigit* ptr = nullptr;
154 temp_tree->Branch(temp_br_name.Data(), &ptr);
155 for (int i = 0; i < fArrayOfDigits->GetEntriesFast(); i++) {
156 ptr = (BmnScWallDigit*)fArrayOfDigits->At(i);
157 temp_tree->Fill();
158 }
159 }
160
161 TTreeFormula* formula = new TTreeFormula("TTreeformula", formulaString.Data(), temp_tree);
162 if (!formula || formula->GetNdim() == 0)
163 LOG(error) << Form("BmnScWallReconstructor::GetSelectedDigiVector. Bad cuts. Check formula %s from file %s",
164 formulaString.Data(), fRecoCutsFile.Data());
165 // next line is needed. Reason unclear. without it will lose some digis.
166 TObject* ThisDigi = nullptr;
167 temp_tree->SetBranchAddress(temp_br_name.Data(), &ThisDigi);
168 for (int i = 0; i < temp_tree->GetEntries(); i++) {
169 temp_tree->GetEntry(i);
170 formula->GetNdata();
171 if (formula->EvalInstance() > 0) // passed selection
172 fSelectedDigiVector.push_back(fArrayOfDigits->At(i));
173 }
174 formula->Delete();
175 temp_tree->Delete();
176 return fSelectedDigiVector;
177}
178
179void BmnScWallReconstructor::OnlineWrite(const std::unique_ptr<TTree>& resultTree)
180{
181 if (!IsActive())
182 return;
183
184 resultTree->Branch("ScWallEvent", &fBmnScWallEvent);
185 resultTree->Fill();
186}
187
189{
190 printf("Work time of BmnScWallReconstructor: %4.2f sec.\n", fworkTime);
191}
int i
Definition P4_F32vec4.h:22
double GetSignal() const
const auto & GetPositionMap() const
BmnStatus ReadGeometryFromGeoManager(bool getGlobalPosition=true)
static uint32_t GetCellId(uint32_t address)
Return Cell id from address.
void SetSignal(float Signal)
void SetCellId(int CellId)
Position.
Class for experimental data at digi level.
uint32_t GetCellId() const
Class for Bmn ScWall data container in event.
BmnScWallCell * GetCell(uint8_t cell_id)
Cell info.
static const int fgkMaxCells
void reset()
Zero all fields.
BmnScWallReconstructor(bool isExp, bool isGlobalCoordinates)
void SetRecoCutsFile(TString reco_cuts_file)
virtual void Exec(Option_t *opt)
virtual void OnlineWrite(const std::unique_ptr< TTree > &resultTree)
Write task resul to tree.