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