1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
19 // This class analyses TPC cosmics data from clusters
21 // Authors: Jan.Fiete.Grosse-Oetringhaus@cern.ch, Claus.Jorgensen@cern.ch
24 #include "AliROCClusterAnalysisSelector.h"
27 #include <AliTPCclusterMI.h>
28 #include <AliRunLoader.h>
30 #include <AliTPCClustersRow.h>
39 #include <TObjArray.h>
40 #include <TTimeStamp.h>
43 #include "TPC/AliTPCClusterHistograms.h"
45 extern TSystem* gSystem;
47 ClassImp(AliROCClusterAnalysisSelector)
49 AliROCClusterAnalysisSelector::AliROCClusterAnalysisSelector() :
54 // Constructor. Initialization of pointers
57 fNMaxObjectsToSave = 50;
58 fObjectsToSave = new TObjArray();
60 for (Int_t i=0; i<kTPCHists; i++)
61 fClusterHistograms[i] = 0;
64 AliROCClusterAnalysisSelector::~AliROCClusterAnalysisSelector()
71 void AliROCClusterAnalysisSelector::SlaveBegin(TTree* tree)
75 AliSelectorRL::SlaveBegin(tree);
78 void AliROCClusterAnalysisSelector::Init(TTree *tree)
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.
86 AliSelectorRL::Init(tree);
90 tree->SetBranchStatus("*", 0);
91 tree->SetBranchStatus("fTimeStamp", 1);
95 Bool_t AliROCClusterAnalysisSelector::Process(Long64_t entry)
98 // Implement your analysis here. Do not forget to call the parent class Process by
99 // if (AliSelectorRL::Process(entry) == kFALSE)
103 if (AliSelectorRL::Process(entry) == kFALSE)
108 for (Int_t i=0; i<kTPCHists; i++)
109 if (fClusterHistograms[i])
110 fClusterHistograms[i]->StartEvent();
112 // runLoader->Dump();
114 Int_t flag = ProcessEvent(entry, kFALSE);
116 ProcessEvent(entry, kTRUE, "");
121 if (fESD->GetTimeStamp()>1160000000) {
122 time = fESD->GetTimeStamp();
126 for (Int_t i=0; i<kTPCHists; i++)
127 if (fClusterHistograms[i])
128 fClusterHistograms[i]->FinishEvent(time);
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....
141 Int_t AliROCClusterAnalysisSelector::ProcessEvent(Long64_t entry, Bool_t detailedHistogram, const Char_t* label)
144 // Looping over clusters in event and filling histograms
146 // - if detailedHistogram = kTRUE special histograms are saved (in fObjectsToSave)
149 // save maximum N objects
150 if (detailedHistogram)
151 if (fObjectsToSave->GetEntries() > fNMaxObjectsToSave)
154 // TODO: find a clever way to handle the time
157 if (fESD->GetTimeStamp()>1160000000) {
158 time = fESD->GetTimeStamp();
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;
168 if (fClusterHistograms[i]) {
170 keepEvent[i] = fClusterHistograms[i]->KeepThisEvent(why);
174 Bool_t intToReturn = 0;
178 AliRunLoader* runLoader = GetRunLoader();
179 runLoader->LoadRecPoints("TPC");
182 AliDebug(AliLog::kError, " Run loader not found");
186 TTree* tree = runLoader->GetTreeR("TPC", kFALSE);
188 // load clusters to the memory
189 AliTPCClustersRow* clrow = new AliTPCClustersRow;
190 clrow->SetClass("AliTPCclusterMI");
192 clrow->GetArray()->ExpandCreateFast(10000);
196 AliDebug(AliLog::kError, " TPC cluster tree not found");
200 TBranch* br = tree->GetBranch("Segment");
201 br->SetAddress(&clrow);
204 Int_t j = Int_t(tree->GetEntries());
205 for (Int_t i=0; i<j; i++) {
208 for (Int_t icl=0; icl<clrow->GetArray()->GetEntriesFast(); icl++){
210 AliTPCclusterMI* cluster = (AliTPCclusterMI*)clrow->GetArray()->At(icl);
215 Int_t detector = cluster->GetDetector();
217 if (detector < 0 || detector >= kTPCSectors)
219 AliDebug(AliLog::kDebug, Form("We found a cluster from invalid sector %d", detector));
223 if (!detailedHistogram) {
225 if (!fClusterHistograms[detector])
226 fClusterHistograms[detector] = new AliTPCClusterHistograms(detector,"",time,time+5*60*60);
228 if (!fClusterHistograms[detector+kTPCSectors])
229 fClusterHistograms[detector+kTPCSectors] = new AliTPCClusterHistograms(detector,"",time,time+5*60*60, kTRUE);
231 fClusterHistograms[detector]->FillCluster(cluster, time);
232 fClusterHistograms[detector+kTPCSectors]->FillCluster(cluster, time);
234 } // end of if !detailedHistograms
236 // if we need the detailed histograms for this event
238 if (keepEvent[detector]) {
240 TString why(fClusterHistograms[detector]->WhyKeepEvent());
242 why.Append(Form("_entry_%d",entry));
245 // if clusterHistograms for this event is not there, construct it
246 if (!clusterHistograms[detector]) {
247 clusterHistograms[detector] = new AliTPCClusterHistograms(detector, why.Data());
249 // adding file and time as comment
250 TString comment = TString(Form("%s",fTree->GetCurrentFile()->GetName()));
251 comment.Append(Form(" entry %d", entry));
253 TString timeStr(TTimeStamp(time).AsString());
256 comment.Append(Form(" (%s)",timeStr.Data()));
258 clusterHistograms[detector]->SetCommentToHistograms(comment.Data());
262 clusterHistograms[detector]->FillCluster(cluster);
263 } // end of (keep this event)
268 if (!detailedHistogram) {
269 for (Int_t i=0; i<kTPCSectors; i++)
270 if (fClusterHistograms[i]) {
272 if (fClusterHistograms[i]->KeepThisEvent(why))
277 for (Int_t i=0; i< kTPCSectors; i++) {
278 if (clusterHistograms[i]) {
279 fObjectsToSave->Add(clusterHistograms[i]);
290 void AliROCClusterAnalysisSelector::SlaveTerminate()
296 for (Int_t i=0; i<kTPCHists; i++)
297 if (fClusterHistograms[i])
298 fOutput->Add(fClusterHistograms[i]);
302 void AliROCClusterAnalysisSelector::Terminate()
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
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)));
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"));
327 AliDebug(AliLog::kInfo, Form("INFO: Found comment in input list: %s", comment->GetTitle()));
332 TFile* file = TFile::Open(Form("rocCluster_%s.root",comment->GetTitle()), "RECREATE");
334 for (Int_t i=0; i<kTPCHists; i++)
335 if (fClusterHistograms[i]) {
336 fClusterHistograms[i]->SaveHistograms();
338 // TCanvas* c = fClusterHistograms[i]->DrawHistograms();
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()));
349 gDirectory->mkdir("saved_objects");
350 gDirectory->cd("saved_objects");
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();
358 fObjectsToSave->At(i)->Write();
362 gDirectory->cd("../");