]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCReconstructor.cxx
Seting the default cotg range to +-1.05
[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"
32
33
34ClassImp(AliTPCReconstructor)
35
8018bb90 36Double_t AliTPCReconstructor::fgCtgRange = 1.05;
59697224 37
38//_____________________________________________________________________________
39void AliTPCReconstructor::Reconstruct(AliRunLoader* runLoader) const
40{
41// reconstruct clusters
42
43 AliLoader* loader = runLoader->GetLoader("TPCLoader");
44 if (!loader) {
45 Error("Reconstruct", "TPC loader not found");
46 return;
47 }
48 loader->LoadRecPoints("recreate");
49 loader->LoadDigits("read");
50
51 AliTPCParam* param = GetTPCParam(runLoader);
52 if (!param) return;
53 AliTPCclustererMI clusterer(param);
54 Int_t nEvents = runLoader->GetNumberOfEvents();
55
56 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
57 runLoader->GetEvent(iEvent);
58
59 TTree* treeClusters = loader->TreeR();
60 if (!treeClusters) {
61 loader->MakeTree("R");
62 treeClusters = loader->TreeR();
63 }
64 TTree* treeDigits = loader->TreeD();
65 if (!treeDigits) {
66 Error("Reconstruct", "Can't get digits tree !");
67 return;
68 }
69
70 clusterer.SetInput(treeDigits);
71 clusterer.SetOutput(treeClusters);
72 clusterer.Digits2Clusters();
73
74 loader->WriteRecPoints("OVERWRITE");
75 }
76
77 loader->UnloadRecPoints();
78 loader->UnloadDigits();
79}
80
38e6e547 81//_____________________________________________________________________________
82void AliTPCReconstructor::Reconstruct(AliRunLoader* runLoader,
83 AliRawReader* rawReader) const
84{
85// reconstruct clusters from raw data
86
87 AliLoader* loader = runLoader->GetLoader("TPCLoader");
88 if (!loader) {
89 Error("Reconstruct", "TPC loader not found");
90 return;
91 }
92 loader->LoadRecPoints("recreate");
93
94 AliTPCParam* param = GetTPCParam(runLoader);
95 if (!param) return;
96 AliTPCclustererMI clusterer(param);
97
98 Int_t iEvent = 0;
99 while (rawReader->NextEvent()) {
100 runLoader->GetEvent(iEvent++);
101
102 TTree* treeClusters = loader->TreeR();
103 if (!treeClusters) {
104 loader->MakeTree("R");
105 treeClusters = loader->TreeR();
106 }
107
108 clusterer.SetOutput(treeClusters);
109 clusterer.Digits2Clusters(rawReader);
110
111 loader->WriteRecPoints("OVERWRITE");
112 }
113
114 loader->UnloadRecPoints();
115}
116
59697224 117//_____________________________________________________________________________
118AliTracker* AliTPCReconstructor::CreateTracker(AliRunLoader* runLoader) const
119{
120// create a TPC tracker
121
122 AliTPCParam* param = GetTPCParam(runLoader);
123 if (!param) return NULL;
124 return new AliTPCtrackerMI(param);
125}
126
127//_____________________________________________________________________________
128void AliTPCReconstructor::FillESD(AliRunLoader* /*runLoader*/,
129 AliESD* esd) const
130{
131// make PID
132
133 Double_t parTPC[] = {47., 0.10, 10.};
134 AliTPCpidESD tpcPID(parTPC);
135 tpcPID.MakePID(esd);
136}
137
138
139//_____________________________________________________________________________
140AliTPCParam* AliTPCReconstructor::GetTPCParam(AliRunLoader* runLoader) const
141{
142// get the TPC parameters
143
6d75e4b6 144 TDirectory* saveDir = gDirectory;
145 runLoader->CdGAFile();
146
147 AliTPCParam* param = (AliTPCParam*) gDirectory->Get("75x40_100x60_150x60");
148 if (!param) Error("GetTPCParam", "no TPC parameters found");
149
150 saveDir->cd();
151 return param;
59697224 152}