]>
Commit | Line | Data |
---|---|---|
fe8871b7 | 1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | /* $Id$ */ | |
17 | ||
18 | // | |
19 | // This class analyses TPC cosmics data from clusters | |
20 | // | |
21 | // Authors: Jan.Fiete.Grosse-Oetringhaus@cern.ch, Claus.Jorgensen@cern.ch | |
22 | // | |
23 | ||
24 | #include "AliROCClusterAnalysisSelector.h" | |
25 | ||
26 | #include <AliLog.h> | |
27 | #include <AliTPCclusterMI.h> | |
28 | #include <AliRunLoader.h> | |
29 | ||
30 | #include <AliTPCClustersRow.h> | |
31 | #include <AliESD.h> | |
32 | ||
33 | ||
34 | #include <TFile.h> | |
35 | #include <TMath.h> | |
36 | #include <TTree.h> | |
37 | #include <TCanvas.h> | |
38 | #include <TSystem.h> | |
39 | #include <TObjArray.h> | |
40 | #include <TTimeStamp.h> | |
41 | #include <TRandom.h> | |
42 | ||
43 | #include "TPC/AliTPCClusterHistograms.h" | |
44 | ||
45 | extern TSystem* gSystem; | |
46 | ||
47 | ClassImp(AliROCClusterAnalysisSelector) | |
48 | ||
49 | AliROCClusterAnalysisSelector::AliROCClusterAnalysisSelector() : | |
50 | AliSelectorRL(), | |
51 | fObjectsToSave(0) | |
52 | { | |
53 | // | |
54 | // Constructor. Initialization of pointers | |
55 | // | |
56 | ||
57 | fNMaxObjectsToSave = 50; | |
58 | fObjectsToSave = new TObjArray(); | |
59 | ||
60 | for (Int_t i=0; i<kTPCHists; i++) | |
61 | fClusterHistograms[i] = 0; | |
62 | } | |
63 | ||
64 | AliROCClusterAnalysisSelector::~AliROCClusterAnalysisSelector() | |
65 | { | |
66 | // | |
67 | // Destructor | |
68 | // | |
69 | } | |
70 | ||
71 | void AliROCClusterAnalysisSelector::SlaveBegin(TTree* tree) | |
72 | { | |
73 | // | |
74 | ||
75 | AliSelectorRL::SlaveBegin(tree); | |
76 | } | |
77 | ||
78 | void AliROCClusterAnalysisSelector::Init(TTree *tree) | |
79 | { | |
80 | // The Init() function is called when the selector needs to initialize | |
81 | // a new tree or chain. Typically here the branch addresses of the tree | |
82 | // will be set. It is normaly not necessary to make changes to the | |
83 | // generated code, but the routine can be extended by the user if needed. | |
84 | // Init() will be called many times when running with PROOF. | |
85 | ||
86 | AliSelectorRL::Init(tree); | |
87 | ||
88 | // Set branch address | |
89 | if (tree) { | |
90 | tree->SetBranchStatus("*", 0); | |
91 | tree->SetBranchStatus("fTimeStamp", 1); | |
92 | } | |
93 | } | |
94 | ||
95 | Bool_t AliROCClusterAnalysisSelector::Process(Long64_t entry) | |
96 | { | |
97 | // | |
98 | // Implement your analysis here. Do not forget to call the parent class Process by | |
99 | // if (AliSelectorRL::Process(entry) == kFALSE) | |
100 | // return kFALSE; | |
101 | // | |
102 | ||
103 | if (AliSelectorRL::Process(entry) == kFALSE) | |
104 | return kFALSE; | |
105 | ||
106 | ||
107 | // reset counters | |
108 | for (Int_t i=0; i<kTPCHists; i++) | |
109 | if (fClusterHistograms[i]) | |
110 | fClusterHistograms[i]->StartEvent(); | |
111 | ||
112 | // runLoader->Dump(); | |
113 | ||
114 | Int_t flag = ProcessEvent(entry, kFALSE); | |
115 | if (flag != 0) | |
116 | ProcessEvent(entry, kTRUE, ""); | |
117 | ||
118 | ||
119 | Int_t time = 0; | |
120 | if (fESD) | |
121 | if (fESD->GetTimeStamp()>1160000000) { | |
122 | time = fESD->GetTimeStamp(); | |
123 | } | |
124 | ||
125 | // finish event | |
126 | for (Int_t i=0; i<kTPCHists; i++) | |
127 | if (fClusterHistograms[i]) | |
128 | fClusterHistograms[i]->FinishEvent(time); | |
129 | ||
130 | ||
131 | // TODO This should not be needed, the TTree::GetEntry() should take care of this, maybe because it has a reference member, to be analyzed | |
132 | // if the ESDfriend is not deleted we get a major memory leak | |
133 | // here the esdfriend seems to be also deleted, very weird behaviour.... | |
134 | ||
135 | delete fESD; | |
136 | fESD = 0; | |
137 | ||
138 | return kTRUE; | |
139 | } | |
140 | ||
141 | Int_t AliROCClusterAnalysisSelector::ProcessEvent(Long64_t entry, Bool_t detailedHistogram, const Char_t* label) | |
142 | { | |
143 | // | |
144 | // Looping over clusters in event and filling histograms | |
145 | // | |
146 | // - if detailedHistogram = kTRUE special histograms are saved (in fObjectsToSave) | |
147 | // | |
148 | ||
149 | // save maximum N objects | |
150 | if (detailedHistogram) | |
151 | if (fObjectsToSave->GetEntries() > fNMaxObjectsToSave) | |
152 | return 0; | |
153 | ||
154 | // TODO: find a clever way to handle the time | |
155 | Int_t time = 0; | |
156 | if (fESD) | |
157 | if (fESD->GetTimeStamp()>1160000000) { | |
158 | time = fESD->GetTimeStamp(); | |
159 | } | |
160 | ||
161 | // for saving single events | |
162 | AliTPCClusterHistograms* clusterHistograms[kTPCSectors]; | |
163 | Bool_t keepEvent[kTPCSectors]; | |
164 | for (Int_t i=0; i<kTPCSectors; i++) { | |
165 | clusterHistograms[i] = 0; | |
166 | keepEvent[i] = kFALSE; | |
167 | ||
168 | if (fClusterHistograms[i]) { | |
169 | TString why; | |
170 | keepEvent[i] = fClusterHistograms[i]->KeepThisEvent(why); | |
171 | } | |
172 | } | |
173 | ||
174 | Bool_t intToReturn = 0; | |
175 | ||
176 | // -------------- | |
177 | ||
178 | AliRunLoader* runLoader = GetRunLoader(); | |
179 | runLoader->LoadRecPoints("TPC"); | |
180 | ||
181 | if (!runLoader) { | |
182 | AliDebug(AliLog::kError, " Run loader not found"); | |
183 | return kFALSE; | |
184 | } | |
185 | ||
186 | TTree* tree = runLoader->GetTreeR("TPC", kFALSE); | |
187 | ||
188 | // load clusters to the memory | |
189 | AliTPCClustersRow* clrow = new AliTPCClustersRow; | |
190 | clrow->SetClass("AliTPCclusterMI"); | |
191 | clrow->SetArray(0); | |
192 | clrow->GetArray()->ExpandCreateFast(10000); | |
193 | // | |
194 | ||
195 | if (!tree) { | |
196 | AliDebug(AliLog::kError, " TPC cluster tree not found"); | |
197 | return kFALSE; | |
198 | } | |
199 | ||
200 | TBranch* br = tree->GetBranch("Segment"); | |
201 | br->SetAddress(&clrow); | |
202 | // | |
203 | ||
204 | Int_t j = Int_t(tree->GetEntries()); | |
205 | for (Int_t i=0; i<j; i++) { | |
206 | br->GetEntry(i); | |
207 | // | |
208 | for (Int_t icl=0; icl<clrow->GetArray()->GetEntriesFast(); icl++){ | |
209 | ||
210 | AliTPCclusterMI* cluster = (AliTPCclusterMI*)clrow->GetArray()->At(icl); | |
211 | ||
212 | if (!cluster) | |
213 | continue; | |
214 | ||
215 | Int_t detector = cluster->GetDetector(); | |
216 | ||
217 | if (detector < 0 || detector >= kTPCSectors) | |
218 | { | |
219 | AliDebug(AliLog::kDebug, Form("We found a cluster from invalid sector %d", detector)); | |
220 | continue; | |
221 | } | |
222 | ||
223 | if (!detailedHistogram) { | |
224 | ||
225 | if (!fClusterHistograms[detector]) | |
226 | fClusterHistograms[detector] = new AliTPCClusterHistograms(detector,"",time,time+5*60*60); | |
227 | ||
228 | if (!fClusterHistograms[detector+kTPCSectors]) | |
229 | fClusterHistograms[detector+kTPCSectors] = new AliTPCClusterHistograms(detector,"",time,time+5*60*60, kTRUE); | |
230 | ||
231 | fClusterHistograms[detector]->FillCluster(cluster, time); | |
232 | fClusterHistograms[detector+kTPCSectors]->FillCluster(cluster, time); | |
233 | ||
234 | } // end of if !detailedHistograms | |
235 | else { | |
236 | // if we need the detailed histograms for this event | |
237 | ||
238 | if (keepEvent[detector]) { | |
239 | ||
240 | TString why(fClusterHistograms[detector]->WhyKeepEvent()); | |
241 | ||
242 | why.Append(Form("_entry_%d",entry)); | |
243 | why.Append(label); | |
244 | ||
245 | // if clusterHistograms for this event is not there, construct it | |
246 | if (!clusterHistograms[detector]) { | |
247 | clusterHistograms[detector] = new AliTPCClusterHistograms(detector, why.Data()); | |
248 | ||
249 | // adding file and time as comment | |
250 | TString comment = TString(Form("%s",fTree->GetCurrentFile()->GetName())); | |
251 | comment.Append(Form(" entry %d", entry)); | |
252 | if (time!=0) { | |
253 | TString timeStr(TTimeStamp(time).AsString()); | |
254 | timeStr.Remove(26); | |
255 | ||
256 | comment.Append(Form(" (%s)",timeStr.Data())); | |
257 | } | |
258 | clusterHistograms[detector]->SetCommentToHistograms(comment.Data()); | |
259 | ||
260 | } | |
261 | else | |
262 | clusterHistograms[detector]->FillCluster(cluster); | |
263 | } // end of (keep this event) | |
264 | } // | |
265 | } | |
266 | } | |
267 | ||
268 | if (!detailedHistogram) { | |
269 | for (Int_t i=0; i<kTPCSectors; i++) | |
270 | if (fClusterHistograms[i]) { | |
271 | TString why; | |
272 | if (fClusterHistograms[i]->KeepThisEvent(why)) | |
273 | intToReturn = 1; | |
274 | } | |
275 | } | |
276 | else { | |
277 | for (Int_t i=0; i< kTPCSectors; i++) { | |
278 | if (clusterHistograms[i]) { | |
279 | fObjectsToSave->Add(clusterHistograms[i]); | |
280 | } | |
281 | } | |
282 | } | |
283 | ||
284 | delete clrow; | |
285 | ||
286 | return intToReturn; | |
287 | } | |
288 | ||
289 | ||
290 | void AliROCClusterAnalysisSelector::SlaveTerminate() | |
291 | { | |
292 | // | |
293 | ||
294 | if (fOutput) | |
295 | { | |
296 | for (Int_t i=0; i<kTPCHists; i++) | |
297 | if (fClusterHistograms[i]) | |
298 | fOutput->Add(fClusterHistograms[i]); | |
299 | } | |
300 | } | |
301 | ||
302 | void AliROCClusterAnalysisSelector::Terminate() | |
303 | { | |
304 | // | |
305 | // read the objects from the output list and write them to a file | |
306 | // the filename is modified by the object comment passed in the tree info or input list | |
307 | // | |
308 | ||
309 | if (fOutput) | |
310 | { | |
311 | fOutput->Print(); | |
312 | ||
313 | for (Int_t i=0; i<kTPCSectors; i++) | |
314 | fClusterHistograms[i] = dynamic_cast<AliTPCClusterHistograms*> (fOutput->FindObject(AliTPCClusterHistograms::FormDetectorName(i, kFALSE))); | |
315 | for (Int_t i=0; i<kTPCSectors; i++) | |
316 | fClusterHistograms[kTPCSectors+i] = dynamic_cast<AliTPCClusterHistograms*> (fOutput->FindObject(AliTPCClusterHistograms::FormDetectorName(i, kTRUE))); | |
317 | } | |
318 | ||
319 | TNamed* comment = 0; | |
320 | if (fTree && fTree->GetUserInfo()) | |
321 | comment = dynamic_cast<TNamed*>(fTree->GetUserInfo()->FindObject("comment")); | |
322 | if (!comment && fInput) | |
323 | comment = dynamic_cast<TNamed*>(fInput->FindObject("comment")); | |
324 | ||
325 | if (comment) | |
326 | { | |
327 | AliDebug(AliLog::kInfo, Form("INFO: Found comment in input list: %s", comment->GetTitle())); | |
328 | } | |
329 | else | |
330 | return; | |
331 | ||
332 | TFile* file = TFile::Open(Form("rocCluster_%s.root",comment->GetTitle()), "RECREATE"); | |
333 | ||
334 | for (Int_t i=0; i<kTPCHists; i++) | |
335 | if (fClusterHistograms[i]) { | |
336 | fClusterHistograms[i]->SaveHistograms(); | |
337 | ||
338 | // TCanvas* c = fClusterHistograms[i]->DrawHistograms(); | |
339 | // TString dir; | |
340 | // dir.Form("WWW/%s/%s", comment->GetTitle(), c->GetName()); | |
341 | // gSystem->mkdir(dir, kTRUE); | |
342 | // c->SaveAs(Form("%s/plots_%s_%s.eps",dir.Data(),comment->GetTitle(),c->GetName())); | |
343 | // c->SaveAs(Form("%s/plots_%s_%s.gif",dir.Data(),comment->GetTitle(),c->GetName())); | |
344 | ||
345 | // c->Close(); | |
346 | // delete c; | |
347 | } | |
348 | ||
349 | gDirectory->mkdir("saved_objects"); | |
350 | gDirectory->cd("saved_objects"); | |
351 | ||
352 | for (Int_t i=0; i<fObjectsToSave->GetSize(); i++) { | |
353 | if (fObjectsToSave->At(i)) { | |
354 | AliTPCClusterHistograms* clusterHistograms = dynamic_cast<AliTPCClusterHistograms*> (fObjectsToSave->At(i)); | |
355 | if (clusterHistograms) | |
356 | clusterHistograms->SaveHistograms(); | |
357 | else | |
358 | fObjectsToSave->At(i)->Write(); | |
359 | } | |
360 | } | |
361 | ||
362 | gDirectory->cd("../"); | |
363 | ||
364 | ||
365 | file->Close(); | |
366 | } |