]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCReconstructor.cxx
GetClusterFast function implemented (No getter before) (Marian)
[u/mrichter/AliRoot.git] / TPC / AliTPCReconstructor.cxx
CommitLineData
59697224 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// //
20// class for TPC reconstruction //
21// //
22///////////////////////////////////////////////////////////////////////////////
23
24
25#include "AliTPCReconstructor.h"
26#include "AliRunLoader.h"
27#include "AliRun.h"
38e6e547 28#include "AliRawReader.h"
59697224 29#include "AliTPCclustererMI.h"
30#include "AliTPCtrackerMI.h"
31#include "AliTPCpidESD.h"
4fb6310b 32#include "AliTPCParam.h"
90c7886e 33#include "AliTPCParamSR.h"
3d37b790 34#include "AliTPCcalibDB.h"
59697224 35
36ClassImp(AliTPCReconstructor)
37
8018bb90 38Double_t AliTPCReconstructor::fgCtgRange = 1.05;
3f82c4f2 39Double_t AliTPCReconstructor::fgMaxSnpTracker = 0.95; // max tangent in tracker - correspond to 3
40Double_t AliTPCReconstructor::fgMaxSnpTrack = 0.999; // tangent
34acb742 41Int_t AliTPCReconstructor::fgStreamLevel = 0; // stream (debug) level
59697224 42//_____________________________________________________________________________
43void AliTPCReconstructor::Reconstruct(AliRunLoader* runLoader) const
44{
45// reconstruct clusters
46
47 AliLoader* loader = runLoader->GetLoader("TPCLoader");
48 if (!loader) {
49 Error("Reconstruct", "TPC loader not found");
50 return;
51 }
52 loader->LoadRecPoints("recreate");
53 loader->LoadDigits("read");
54
55 AliTPCParam* param = GetTPCParam(runLoader);
56 if (!param) return;
57 AliTPCclustererMI clusterer(param);
58 Int_t nEvents = runLoader->GetNumberOfEvents();
59
60 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
61 runLoader->GetEvent(iEvent);
62
63 TTree* treeClusters = loader->TreeR();
64 if (!treeClusters) {
65 loader->MakeTree("R");
66 treeClusters = loader->TreeR();
67 }
68 TTree* treeDigits = loader->TreeD();
69 if (!treeDigits) {
70 Error("Reconstruct", "Can't get digits tree !");
71 return;
72 }
73
74 clusterer.SetInput(treeDigits);
75 clusterer.SetOutput(treeClusters);
76 clusterer.Digits2Clusters();
77
78 loader->WriteRecPoints("OVERWRITE");
79 }
80
81 loader->UnloadRecPoints();
82 loader->UnloadDigits();
83}
84
38e6e547 85//_____________________________________________________________________________
86void AliTPCReconstructor::Reconstruct(AliRunLoader* runLoader,
87 AliRawReader* rawReader) const
88{
89// reconstruct clusters from raw data
90
91 AliLoader* loader = runLoader->GetLoader("TPCLoader");
92 if (!loader) {
93 Error("Reconstruct", "TPC loader not found");
94 return;
95 }
96 loader->LoadRecPoints("recreate");
97
98 AliTPCParam* param = GetTPCParam(runLoader);
90c7886e 99 if (!param) {
100 AliWarning("Loading default TPC parameters !");
101 param = new AliTPCParamSR;
102 }
38e6e547 103 AliTPCclustererMI clusterer(param);
104
b6f060dc 105 TString option = GetOption();
106 if (option.Contains("PedestalSubtraction"))
107 clusterer.SetPedSubtraction(kTRUE);
c8c88e5d 108 if (option.Contains("OldRCUFormat"))
109 clusterer.SetOldRCUFormat(kTRUE);
b6f060dc 110
38e6e547 111 Int_t iEvent = 0;
3d37b790 112 while (rawReader->NextEvent() && iEvent<35) {
38e6e547 113 runLoader->GetEvent(iEvent++);
114
115 TTree* treeClusters = loader->TreeR();
116 if (!treeClusters) {
117 loader->MakeTree("R");
118 treeClusters = loader->TreeR();
119 }
120
121 clusterer.SetOutput(treeClusters);
122 clusterer.Digits2Clusters(rawReader);
123
124 loader->WriteRecPoints("OVERWRITE");
125 }
126
127 loader->UnloadRecPoints();
128}
129
59697224 130//_____________________________________________________________________________
131AliTracker* AliTPCReconstructor::CreateTracker(AliRunLoader* runLoader) const
132{
133// create a TPC tracker
134
135 AliTPCParam* param = GetTPCParam(runLoader);
90c7886e 136 if (!param) {
137 AliWarning("Loading default TPC parameters !");
138 param = new AliTPCParamSR;
139 }
4fb6310b 140 param->ReadGeoMatrices();
59697224 141 return new AliTPCtrackerMI(param);
142}
143
144//_____________________________________________________________________________
145void AliTPCReconstructor::FillESD(AliRunLoader* /*runLoader*/,
146 AliESD* esd) const
147{
148// make PID
149
150 Double_t parTPC[] = {47., 0.10, 10.};
151 AliTPCpidESD tpcPID(parTPC);
152 tpcPID.MakePID(esd);
153}
154
155
156//_____________________________________________________________________________
3d37b790 157AliTPCParam* AliTPCReconstructor::GetTPCParam(AliRunLoader* /*runLoader*/) const
59697224 158{
159// get the TPC parameters
160
3d37b790 161// TDirectory* saveDir = gDirectory;
162//runLoader->CdGAFile();
6d75e4b6 163
3d37b790 164 AliTPCParam* param = AliTPCcalibDB::Instance()->GetParameters();
6d75e4b6 165
3d37b790 166 // AliTPCParam* param = (AliTPCParam*) gDirectory->Get("75x40_100x60_150x60");
167// if (!param) Error("GetTPCParam", "no TPC parameters found");
168
169// saveDir->cd();
6d75e4b6 170 return param;
59697224 171}