CbmRoot
Loading...
Searching...
No Matches
CbmConverterManager.cxx
Go to the documentation of this file.
1/* Copyright (C) 2020-2021 Physikalisches Institut, Eberhard Karls Universitaet Tuebingen, Tuebingen
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Viktor Klochkov [committer] */
4
6
7#include "AnalysisTree/DataHeader.hpp"
8#include "AnalysisTree/TaskManager.hpp"
9#include "CbmConverterTask.h"
10#include "CbmDefs.h"
11#include "CbmEvent.h"
12#include "TClonesArray.h"
13#include "TGeoBBox.h"
14#include "TGeoManager.h"
15#include "TString.h"
16
17#include <iostream>
18
20
22{
23 task_manager_->Init();
25 InitEvent();
26 return kSUCCESS;
27}
28
30{
31 tasks_.emplace_back(task);
32 task_manager_->AddTask(reinterpret_cast<AnalysisTree::Task*>(task));
33}
34
36{
37 index_map_.clear();
38
39 for (auto* task : tasks_) {
40 task->SetIndexesMap(&index_map_);
41 task->ProcessData(event);
42 index_map_.insert(std::make_pair(task->GetOutputBranchName(), task->GetOutIndexesMap()));
43 }
44 task_manager_->FillOutput();
45}
46
47void CbmConverterManager::Exec(Option_t* /*opt*/)
48{
49 if (events_ != nullptr) {
50 auto n_events = events_->GetEntriesFast();
51 for (int i_event = 0; i_event < n_events; ++i_event) {
52 auto* event = (CbmEvent*) events_->At(i_event);
53 ProcessData(event);
54 }
55 }
56 else {
57 LOG(info) << "Event based mode\n";
58 ProcessData(nullptr);
59 }
60}
61
62
64{
65 TDirectory* curr = gDirectory; // TODO check why this is needed
66 TFile* currentFile = gFile;
67
68 task_manager_->Finish();
69
70 gFile = currentFile;
71 gDirectory = curr;
72}
73
75{
76 // Force user to write data info //TODO is there a way to read it from a file automatically?
77 assert(!system_.empty() && beam_mom_);
78
79 auto* data_header = new AnalysisTree::DataHeader();
80
81 std::cout << "ReadDataHeader" << std::endl;
82 data_header->SetSystem(system_);
83 data_header->SetBeamMomentum(beam_mom_);
84 data_header->SetTimeSliceLength(ts_length_);
85
86 auto& specDet_mod_pos = data_header->AddDetector();
87 TString specDet_names[2] = {(TString) ToString(ECbmModuleId::kPsd), (TString) ToString(ECbmModuleId::kFsd)};
88 const char* module_name = "module"; // PSD and FSD modules always contains "module"
89
90 TVector3 frontFaceGlobal;
91 Int_t nSpecDetModules = 0;
92 TString specDetNameTag = "";
93
94 TGeoIterator geoIterator(gGeoManager->GetTopNode()->GetVolume());
95 TGeoNode* curNode;
96 geoIterator.Reset(); // safety to reset to "cave" before the loop starts
97 while ((curNode = geoIterator())) {
98 TString nodePath;
99 geoIterator.GetPath(nodePath);
100 if (!(nodePath.Contains(specDet_names[0], TString::kIgnoreCase)
101 || nodePath.Contains(specDet_names[1], TString::kIgnoreCase))) {
102 geoIterator.Skip(); // skip current branch when it is not FSD => should speed up
103 continue; // don't do anything for this branch
104 }
105
106 TString nodeName(curNode->GetName());
107
108 // spectator detector as a whole
109 if (nodeName.Contains(specDet_names[0], TString::kIgnoreCase)
110 || nodeName.Contains(specDet_names[1], TString::kIgnoreCase)) {
111 specDetNameTag = nodeName;
112
113 auto specDetGeoMatrix = curNode->GetMatrix();
114 auto specDetBox = (TGeoBBox*) curNode->GetVolume()->GetShape();
115 TVector3 frontFaceLocal(0, 0, -specDetBox->GetDZ());
116 specDetGeoMatrix->LocalToMaster(&frontFaceLocal[0], &frontFaceGlobal[0]);
117 }
118
119 // modules of spectator detector
120 if (nodeName.Contains(module_name)) {
121 nSpecDetModules++;
122
123 auto geoMatrix = curNode->GetMatrix();
124 TVector3 translation(geoMatrix->GetTranslation());
125 double x = translation.X();
126 double y = translation.Y();
127
128 auto* module = specDet_mod_pos.AddChannel();
129 module->SetPosition(x, y, frontFaceGlobal[2]);
130 LOG(info) << "Module " << nSpecDetModules << " : " << Form("(%.3f, %.3f)", x, y) << " id "
131 << curNode->GetNumber();
132 }
133 }
134
135 LOG(info) << "Detector " << specDetNameTag << " with " << nSpecDetModules << " modules";
136
137 task_manager_->SetOutputDataHeader(data_header);
138}
ClassImp(CbmConverterManager)
std::string ToString(ECbmModuleId modId)
Definition CbmDefs.cxx:70
@ kPsd
Projectile spectator detector.
@ kFsd
Forward spectator detector.
void Exec(Option_t *opt) override
void AddTask(CbmConverterTask *task)
std::map< std::string, std::map< int, int > > index_map_
map CbmRoot to AT of indexes for a given branch
AnalysisTree::TaskManager * task_manager_
InitStatus Init() override
~CbmConverterManager() override
std::vector< CbmConverterTask * > tasks_
void ProcessData(CbmEvent *event)
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34