ARTEMIS-CRIB
 
Loading...
Searching...
No Matches
TSegmentOutputProcessor.cc
Go to the documentation of this file.
1/**
2 * @file TSegmentOutputProcessor.cc
3 * @brief from seg conf, output raw data TTree object
4 * @author Kodai Okawa <okawa@cns.s.u-tokyo.ac.jp>
5 * @date 2024-08-21 17:36:19
6 * @note last modified: 2024-08-23 20:56:40
7 * @details
8 */
9
11
12#include <TArtemisUtil.h>
13#include <TModuleData.h>
14#include <TModuleType.h>
15#include <TROOT.h>
16#include <TRawDataObject.h>
17#include <TSegmentInfo.h>
18#include <TSegmentedData.h>
19#include <TTree.h>
20#include <constant.h>
21
23
25
27 : fFile(nullptr), fTree(nullptr), fSegmentList(nullptr), fModuleList(nullptr),
28 fSegmentedData(nullptr) {
29 RegisterProcessorParameter("FileName", "The name of output file", fFileName, TString("tmp.root"));
30 RegisterProcessorParameter("TreeName", "The name of output tree", fTreeName, TString("tree"));
31
32 StringVec_t defaultignore;
33 defaultignore.push_back("");
34 RegisterInputCollection("SegmentedDataName", "name of the segmented data",
35 fSegmentedDataName, TString("segdata"),
36 &fSegmentedData, "art::TSegmentedData");
37
38 RegisterInputInfo("SegmentList", "name of the segment list",
39 fSegmentListName, TString("seglist"), &fSegmentList, "TClonesArray", "art::TSegmentInfo");
40 RegisterInputInfo("ModuleList", "name of the module list",
41 fModuleListName, TString("modlist"), &fModuleList, "TClonesArray", "art::TModuleType");
42 RegisterProcessorParameter("Ignore", "ignore segment list", fIgnore, defaultignore);
43}
44
46 if (gROOT->GetListOfFiles()->FindObject(fFile)) {
47 if (fTree)
48 fTree->GetUserInfo()->Clear();
49 fFile->Close();
50 }
51
52 // input processors are deleted elsewhere
53 fSegmentedData = nullptr;
54 fSegmentList = nullptr;
55 fModuleList = nullptr;
56}
57
59
61 if (this != &rhs) {
62 }
63 return *this;
64}
65
66void TSegmentOutputProcessor::Init(TEventCollection *) {
67 Util::PrepareDirectoryFor(fFileName);
68 fFile = TFile::Open(fFileName, "RECREATE");
69 if (!fFile) {
70 SetStateError(TString::Format("Cannot create file: %s", fFileName.Data()));
71 return;
72 }
73 fTree = new TTree(fTreeName, fTreeName);
74 Info("init", "Created %s", fFileName.Data());
75
76 Int_t nSeg = (*fSegmentList)->GetEntriesFast();
77 for (Int_t iSeg = 0; iSeg != nSeg; iSeg++) {
78 TSegmentInfo *seg = (TSegmentInfo *)(*fSegmentList)->At(iSeg);
79 Int_t nMod = seg->GetNumModules();
80
81 // skip if no module registered
82 if (!nMod)
83 continue;
84 if (fIgnore.end() != std::find(fIgnore.begin(), fIgnore.end(), seg->GetName()))
85 continue;
86
87 std::pair<int, std::vector<TModuleData *>> segment;
88 segment.first = seg->GetSegID();
89 std::vector<TModuleData *> &modules = fSegments.insert(segment).first->second;
90
91 for (Int_t iMod = 0; iMod != nMod; iMod++) {
92 TModuleInfo *mod_info = seg->GetModule(iMod);
93 TModuleData *mod = new TModuleData(*mod_info);
94 Int_t id = mod->GetID(); // geo ID
95 if (modules.size() < (unsigned int)id + 1)
96 modules.resize(id + 1);
97 modules.at(id) = mod;
98 TModuleType *type = (TModuleType *)(*fModuleList)->FindObject(seg->GetModuleType());
99 Int_t nCh = type->GetNch();
100 mod->SetCh(nCh);
101
102 Int_t mod_id = type->GetDecoderID();
103 mod->SetMod(mod_id);
104
105 // prepare branch (id == 24 or 25 -> multihit TDC)
106 if (mod_id == 24 || mod_id == 25) {
107 Info("Init", "Set %s_%d branch", seg->GetName(), id);
108 fTree->Branch(Form("%s_%d", seg->GetName(), id), &(mod->fData2D));
109 } else {
110 Info("Init", "Set %s_%d branch", seg->GetName(), id);
111 fTree->Branch(Form("%s_%d", seg->GetName(), id), &(mod->fData1D));
112 }
113 }
114 }
115
116 for (std::vector<TString>::iterator it = fIgnore.begin(); it != fIgnore.end(); it++) {
117 if ((*it) == "") {
118 continue;
119 }
120 printf("%s will be ignored\n", (*it).Data());
121 }
122}
123
125 std::map<Int_t, std::vector<TModuleData *>>::iterator it;
126
127 for (it = fSegments.begin(); it != fSegments.end(); it++) {
128 std::vector<TModuleData *> &modules = it->second;
129 for (Size_t i = 0; i < modules.size(); i++) {
130 Int_t mod = modules[i]->GetMod();
131 if (mod < 0) {
132 continue;
133 } else if (mod == 24 || mod == 25) {
134 std::vector<Int_t> v;
135 v.emplace_back(kInvalidI);
136 modules[i]->fData2D.assign(modules[i]->GetNCh(), v);
137 } else {
138 modules[i]->fData1D.assign(modules[i]->GetNCh(), kInvalidI);
139 }
140 }
141 }
142
143 for (it = fSegments.begin(); it != fSegments.end(); it++) {
144 Int_t segid = it->first;
145 TObjArray *arr = (*fSegmentedData)->FindSegmentByID(segid);
146 if (!arr) {
147 Warning("Process", "No segment having segid = %d", segid);
148 Warning("Process", " Add this segid to Ignore if this semgment is not valid temporarily");
149 SetStopLoop();
150 return;
151 }
152 Int_t nHit = arr->GetEntriesFast();
153 std::vector<TModuleData *> &modules = it->second;
154 for (Int_t iHit = 0; iHit != nHit; iHit++) {
155 TRawDataObject *data = (TRawDataObject *)arr->UncheckedAt(iHit);
156 Int_t geo = data->GetGeo();
157 Int_t ch = data->GetCh();
158
159 // nVal is always 2 (see TRawTiming.h), so only take iVal=0
160#if 0
161 Int_t nVal = data->GetNumValues();
162 Int_t tmp = -100;
163 for (Int_t iVal = 0; iVal < nVal; iVal++) {
164 if (modules.size() > geo && modules[geo] != nullptr) {
165 if (iVal == 0) {
166 tmp = data->GetValue(iVal);
167 } else {
168 if (tmp != data->GetValue(iVal)) {
169 std::cout << tmp << std::endl;
170 std::cout << data->GetValue(iVal) << std::endl;
171 }
172 }
173 // std::cout << data->GetValue(iVal) << std::endl;
174 }
175 }
176#endif
177
178 if ((Int_t)modules.size() > geo && modules[geo] != nullptr) {
179 Int_t mod = modules[geo]->GetMod();
180 if (mod == 24 || mod == 25) {
181 if (modules[geo]->fData2D[ch][0] == kInvalidI) {
182 modules[geo]->fData2D[ch][0] = data->GetValue(0);
183 } else {
184 modules[geo]->fData2D[ch].emplace_back(data->GetValue(0));
185 }
186 } else {
187 modules[geo]->fData1D[ch] = data->GetValue(0);
188 }
189 }
190 }
191 }
192 fTree->Fill();
193}
194
197
199 if (!gDirectory)
200 gDirectory = gROOT;
201 TDirectory *saved = gDirectory;
202 fFile->cd();
203 fTree->Write(0, TFile::kOverwrite);
204 saved->cd();
205}
inherit from TModuleInfo
ClassImp(TSegmentOutputProcessor)
from seg conf, output raw data TTree object
std::vector< Int_t > fData1D
Definition TModuleData.h:37
void SetCh(Int_t Nch)
Definition TModuleData.h:28
void SetMod(Int_t mod)
Definition TModuleData.h:33
std::vector< std::vector< Int_t > > fData2D
Definition TModuleData.h:38
TSegmentOutputProcessor & operator=(const TSegmentOutputProcessor &rhs)
TString fSegmentListName
pure TTree object (not TArtTree)
void Init(TEventCollection *) override
std::map< Int_t, std::vector< TModuleData * > > fSegments
list of ignored segment
return to the guide