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() :
51 fNMaxObjectsToSave(0),
55 // Constructor. Initialization of pointers
58 fNMaxObjectsToSave = 50;
59 fObjectsToSave = new TObjArray();
61 for (Int_t i=0; i<kTPCHists; i++)
62 fClusterHistograms[i] = 0;
65 AliROCClusterAnalysisSelector::~AliROCClusterAnalysisSelector()
72 void AliROCClusterAnalysisSelector::SlaveBegin(TTree* tree)
76 AliSelectorRL::SlaveBegin(tree);
79 void AliROCClusterAnalysisSelector::Init(TTree *tree)
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.
87 AliSelectorRL::Init(tree);
91 tree->SetBranchStatus("*", 0);
92 tree->SetBranchStatus("fTimeStamp", 1);
96 Bool_t AliROCClusterAnalysisSelector::Process(Long64_t entry)
99 // Implement your analysis here. Do not forget to call the parent class Process by
100 // if (AliSelectorRL::Process(entry) == kFALSE)
104 if (AliSelectorRL::Process(entry) == kFALSE)
109 for (Int_t i=0; i<kTPCHists; i++)
110 if (fClusterHistograms[i])
111 fClusterHistograms[i]->StartEvent();
113 // runLoader->Dump();
115 Int_t flag = ProcessEvent(entry, kFALSE);
117 ProcessEvent(entry, kTRUE, "");
122 if (fESD->GetTimeStamp()>1160000000) {
123 time = fESD->GetTimeStamp();
127 for (Int_t i=0; i<kTPCHists; i++)
128 if (fClusterHistograms[i])
129 fClusterHistograms[i]->FinishEvent(time);
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....
142 Int_t AliROCClusterAnalysisSelector::ProcessEvent(Long64_t entry, Bool_t detailedHistogram, const Char_t* label)
145 // Looping over clusters in event and filling histograms
147 // - if detailedHistogram = kTRUE special histograms are saved (in fObjectsToSave)
150 // save maximum N objects
151 if (detailedHistogram)
152 if (fObjectsToSave->GetEntries() > fNMaxObjectsToSave)
155 // TODO: find a clever way to handle the time
158 if (fESD->GetTimeStamp()>1160000000) {
159 time = fESD->GetTimeStamp();
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;
169 if (fClusterHistograms[i]) {
171 keepEvent[i] = fClusterHistograms[i]->KeepThisEvent(why);
175 Bool_t intToReturn = 0;
179 AliRunLoader* runLoader = GetRunLoader();
180 runLoader->LoadRecPoints("TPC");
183 AliDebug(AliLog::kError, " Run loader not found");
187 TTree* tree = runLoader->GetTreeR("TPC", kFALSE);
189 // load clusters to the memory
190 AliTPCClustersRow* clrow = new AliTPCClustersRow;
191 clrow->SetClass("AliTPCclusterMI");
193 clrow->GetArray()->ExpandCreateFast(10000);
197 AliDebug(AliLog::kError, " TPC cluster tree not found");
201 TBranch* br = tree->GetBranch("Segment");
202 br->SetAddress(&clrow);
205 Int_t j = Int_t(tree->GetEntries());
206 for (Int_t i=0; i<j; i++) {
209 for (Int_t icl=0; icl<clrow->GetArray()->GetEntriesFast(); icl++){
211 AliTPCclusterMI* cluster = (AliTPCclusterMI*)clrow->GetArray()->At(icl);
216 Int_t detector = cluster->GetDetector();
218 if (detector < 0 || detector >= kTPCSectors)
220 AliDebug(AliLog::kDebug, Form("We found a cluster from invalid sector %d", detector));
224 if (!detailedHistogram) {
226 if (!fClusterHistograms[detector])
227 fClusterHistograms[detector] = new AliTPCClusterHistograms(detector,"",time,time+5*60*60);
229 if (!fClusterHistograms[detector+kTPCSectors])
230 fClusterHistograms[detector+kTPCSectors] = new AliTPCClusterHistograms(detector,"",time,time+5*60*60, kTRUE);
232 fClusterHistograms[detector]->FillCluster(cluster, time);
233 fClusterHistograms[detector+kTPCSectors]->FillCluster(cluster, time);
235 } // end of if !detailedHistograms
237 // if we need the detailed histograms for this event
239 if (keepEvent[detector]) {
241 TString why(fClusterHistograms[detector]->WhyKeepEvent());
243 why.Append(Form("_entry_%d",entry));
246 // if clusterHistograms for this event is not there, construct it
247 if (!clusterHistograms[detector]) {
248 clusterHistograms[detector] = new AliTPCClusterHistograms(detector, why.Data());
250 // adding file and time as comment
251 TString comment = TString(Form("%s",fTree->GetCurrentFile()->GetName()));
252 comment.Append(Form(" entry %d", entry));
254 TString timeStr(TTimeStamp(time).AsString());
257 comment.Append(Form(" (%s)",timeStr.Data()));
259 clusterHistograms[detector]->SetCommentToHistograms(comment.Data());
263 clusterHistograms[detector]->FillCluster(cluster);
264 } // end of (keep this event)
269 if (!detailedHistogram) {
270 for (Int_t i=0; i<kTPCSectors; i++)
271 if (fClusterHistograms[i]) {
273 if (fClusterHistograms[i]->KeepThisEvent(why))
278 for (Int_t i=0; i< kTPCSectors; i++) {
279 if (clusterHistograms[i]) {
280 fObjectsToSave->Add(clusterHistograms[i]);
291 void AliROCClusterAnalysisSelector::SlaveTerminate()
297 for (Int_t i=0; i<kTPCHists; i++)
298 if (fClusterHistograms[i])
299 fOutput->Add(fClusterHistograms[i]);
303 void AliROCClusterAnalysisSelector::Terminate()
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
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)));
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"));
328 AliDebug(AliLog::kInfo, Form("INFO: Found comment in input list: %s", comment->GetTitle()));
333 TFile* file = TFile::Open(Form("rocCluster_%s.root",comment->GetTitle()), "RECREATE");
335 for (Int_t i=0; i<kTPCHists; i++)
336 if (fClusterHistograms[i]) {
337 fClusterHistograms[i]->SaveHistograms();
339 // TCanvas* c = fClusterHistograms[i]->DrawHistograms();
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()));
350 gDirectory->mkdir("saved_objects");
351 gDirectory->cd("saved_objects");
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();
359 fObjectsToSave->At(i)->Write();
363 gDirectory->cd("../");