]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/FindKrClusters.C
Implementation of copy constructor and assignement operators (Marian)
[u/mrichter/AliRoot.git] / TPC / FindKrClusters.C
CommitLineData
d19e1df9 1/****************************************************************************
2 * Origin: A.Matyja amatyja@cern.ch *
3 ****************************************************************************/
4
5/*
6
7 macro to create array of clusters from TPC digits
8 input files - galice.root
9 digits.root - file with digits - usualy use link to galice.root
10 - in splitted mode - neccesary to create link to proper file
11
12 output file - TPC.RecPoints.root
13
14
15// Warning - if cluster file AliTPCclusters.root already exist - macro exit and don't produce anything
16
17
18*/
19
20
21#ifndef __CINT__
22#include <iostream.h>
23#include "AliRun.h"
24#include "AliTPCv4.h"
25#include "AliTPCParam.h"
26#include "AliTPCclusterKr.h"
27#include "AliTPCclustererKr.h"
28#include "TFile.h"
29#include "TStopwatch.h"
30#include "TTree.h"
31#endif
32
33Int_t FindKrClusters(){
34
35 AliRunLoader* rl = AliRunLoader::Open("galice.root");
36 if (rl == 0x0) {
37 cerr<<"Can not open session"<<endl;
38 return 1;
39 }
40
41 AliTPCLoader *tpcl = (AliTPCLoader*)rl->GetLoader("TPCLoader");
42 if (tpcl == 0x0) {
43 cerr<<"Can not get TPC Loader"<<endl;
44 return 1;
45 }
46
47 if (tpcl->LoadDigits()) {
48 cerr<<"Error occured while loading digits"<<endl;
49 return 1;
50 }
51
52 if (rl->LoadgAlice()) {
53 cerr<<"Error occured while LoadgAlice"<<endl;
54 return 1;
55 }
56
57 gAlice=rl->GetAliRun();
58 if (!gAlice) {
59 cerr<<"Can't get gAlice !\n";
60 return 1;
61 }
62
63 TDirectory *cwd = gDirectory;
64
65 AliTPCv4 *tpc = (AliTPCv4*)gAlice->GetDetector("TPC");
66 Int_t ver = tpc->IsVersion();
67 cerr<<"TPC version "<<ver<<" has been found !\n";
68
69 rl->CdGAFile();
70
71 AliTPCParam *param=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60");
72 if (!param) {cerr<<"TPC parameters have not been found !\n"; return 4;}
73
74 AliTPCDigitsArray *digarr=new AliTPCDigitsArray;
75 digarr->Setup(param);
76
77 cerr<<"It has begun"<<endl;
78 TStopwatch timer;
79 timer.Start();
80 cwd->cd();
81
82 Int_t nevmax=rl->GetNumberOfEvents();//number of events in run
83 for(Int_t nev=0;nev<nevmax /*&& nev<1*/ ;nev++){
84 rl->GetEvent(nev);
85
86 TTree* input_tree= tpcl->TreeD();//tree with digits
87 if (input_tree == 0x0){
88 cerr << "Can not get TreeD for event " <<nev<<endl;
89 continue;
90 }
91
92 digarr->ConnectTree(input_tree);
93
94 TTree *output_tree =tpcl->TreeR();
95 if(output_tree==0x0){
96 tpcl->MakeTree("R");
97 output_tree = tpcl->TreeR();
98 if (output_tree == 0x0){
99 cerr << "Problems with output tree (TreeR) for event "<<nev<<endl;
100 continue;
101 }
102 }
103
104 cout<<"Processing event "<<nev<<endl;
105
106 //test
107 //cout<<"nentries"<<input_tree->GetEntries();
108
109 AliTPCclustererKr *clusters = new AliTPCclustererKr();
110 clusters->SetParam(param);
111 clusters->SetInput(input_tree);
112 clusters->SetOutput(output_tree);
113 clusters->SetDigArr(digarr);
114 clusters->finderIO();
115
116 tpcl->WriteRecPoints("OVERWRITE");
117 }
118
119
120 timer.Stop(); timer.Print();
121
122 delete rl;//cleans everything
123
124 return 0;
125}