]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSreconstruction.cxx
Dependencies on MC truth in reconstruction limited to rec. of MC data for comparison...
[u/mrichter/AliRoot.git] / ITS / AliITSreconstruction.cxx
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 // Class for ITS RecPoint reconstruction                               //
20 //                                                                     //
21 ////////////////////////////////////////////////////////////////////////
22
23 #include <TString.h>
24 #include "AliRun.h"
25 #include "AliRunLoader.h"
26 #include "AliITSDetTypeRec.h"
27 #include "AliITSLoader.h"
28 #include "AliITSreconstruction.h"
29 #include "AliITSgeom.h"
30
31
32 ClassImp(AliITSreconstruction)
33
34 //______________________________________________________________________
35 AliITSreconstruction::AliITSreconstruction():
36  fInit(kFALSE),
37  fEnt(0),
38  fEnt0(0),
39  fDetTypeRec(0x0),
40  fDfArp(kFALSE),
41  fITSgeom(0x0),
42  fLoader(0x0),
43  fRunLoader(0x0)
44 {
45     // Default constructor.
46     // Inputs:
47     //  none.
48     // Outputs:
49     //   none.
50     // Return:
51     //    A zero-ed constructed AliITSreconstruction class.
52     fDet[0] = fDet[1] = fDet[2] = kTRUE;
53 }
54 //______________________________________________________________________
55
56 AliITSreconstruction::AliITSreconstruction(AliRunLoader *rl):
57  fInit(kFALSE),
58  fEnt(0),
59  fEnt0(0),
60  fDetTypeRec(0x0),
61  fDfArp(kFALSE),
62  fITSgeom(0x0),
63  fLoader(0x0),
64  fRunLoader(rl)
65 {
66   fDet[0] = fDet[1] = fDet[2] = kTRUE;
67 }
68 //______________________________________________________________________
69 AliITSreconstruction::AliITSreconstruction(const char* filename):
70  fInit(kFALSE),
71  fEnt(0),
72  fEnt0(0),
73  fDetTypeRec(0x0),
74  fDfArp(kFALSE),
75  fITSgeom(0x0),
76  fLoader(0x0),
77  fRunLoader(0x0)
78 {
79     // Standard constructor.
80     // Inputs:
81     //  const char* filename    filename containing the digits to be
82     //                          reconstructed. If filename = 0 (nil)
83     //                          then no file is opened but a file is
84     //                          assumed to already be opened. This 
85     //                          already opened file will be used.
86     // Outputs:
87     //   none.
88     // Return:
89     //    A standardly constructed AliITSreconstruction class.
90
91     fDet[0] = fDet[1] = fDet[2] = kTRUE;
92
93     fRunLoader = AliRunLoader::Open(filename);
94     if (fRunLoader == 0x0)
95      {
96        Error("AliITSreconstruction","Can not load the session",filename);
97        return;
98      }
99
100 }
101
102 //______________________________________________________________________
103 AliITSreconstruction::AliITSreconstruction(const AliITSreconstruction &rec):TTask(rec),
104 fInit(rec.fInit),
105 fEnt(rec.fEnt),
106 fEnt0(rec.fEnt0),
107 fDetTypeRec(rec.fDetTypeRec),
108 fDfArp(rec.fDfArp),
109 fITSgeom(rec.fITSgeom),
110 fLoader(rec.fLoader),
111 fRunLoader(rec.fRunLoader)
112 {
113     // Copy constructor. 
114
115   
116 }
117
118 //______________________________________________________________________
119 AliITSreconstruction& AliITSreconstruction::operator=(const AliITSreconstruction& source){
120     // Assignment operator. 
121   this->~AliITSreconstruction();
122   new(this) AliITSreconstruction(source);
123   return *this;
124
125 }
126
127 //______________________________________________________________________
128 AliITSreconstruction::~AliITSreconstruction(){
129     //    A destroyed AliITSreconstruction class.
130     
131     //fITS      = 0;
132     delete fRunLoader;
133     
134 }
135 //______________________________________________________________________
136 Bool_t AliITSreconstruction::Init(){
137     // Class Initilizer.
138     // Inputs:
139     //  none.
140     // Outputs:
141     //   none.
142     // Return:
143     //    kTRUE if no errors initilizing this class occurse else kFALSE
144     Info("Init","");
145     if (fRunLoader == 0x0)
146      {
147        Error("Init","Run Loader is NULL");
148        return kFALSE;
149      }
150     //  fRunLoader->LoadgAlice();
151     //   fRunLoader->LoadHeader();  
152
153     fLoader = (AliITSLoader*) fRunLoader->GetLoader("ITSLoader");
154     if(!fLoader) {
155       Error("Init","ITS loader not found");
156       fInit = kFALSE;
157     }
158
159     // Now ready to init.
160  
161     //fRunLoader->CdGAFile();
162     fITSgeom = fLoader->GetITSgeom();
163
164     fDetTypeRec = new AliITSDetTypeRec(fLoader);
165     fDetTypeRec->SetDefaults();
166     fDet[0] = fDet[1] = fDet[2] = kTRUE;
167     fEnt0 = 0;
168
169     //fEnt  = gAlice->GetEventsPerRun();
170     fEnt = Int_t(fRunLoader->GetNumberOfEvents());
171
172     fLoader->LoadDigits("read");
173     fLoader->LoadRecPoints("recreate");
174     fLoader->LoadRawClusters("recreate");
175     if (fLoader->TreeR() == 0x0) fLoader->MakeTree("R");
176     if (fLoader->TreeC() == 0x0) fLoader->MakeTree("C");
177  
178     fDetTypeRec->MakeBranchR(0);
179     fDetTypeRec->MakeBranchC();
180     fDetTypeRec->SetTreeAddress();
181     fDetTypeRec->SetTreeAddressR(fLoader->TreeR());
182
183     fInit = InitRec();
184
185     Info("Init","  Done\n\n\n");
186
187     return fInit;
188 }
189 //______________________________________________________________________
190 Bool_t AliITSreconstruction::InitRec(){
191     // Sets up Reconstruction part of AliITSDetType..
192     // Inputs:
193     //      none.
194     // Outputs:
195     //      none.
196     // Return:
197     //      none.
198   /*
199   //AliITSDetType *idt;
200   fDetTypeRec->SetLoader(fLoader);
201     // SPD
202   if(fDet[kSPD]){
203     Info("InitRec","SPD");
204     //idt = fITS->DetType(kSPD);
205     AliITSsegmentationSPD *segSPD = (AliITSsegmentationSPD*)fDetTypeRec->GetSegmentationModel(0);
206       TClonesArray *digSPD = fDetTypeRec->DigitsAddress(kSPD);
207       TClonesArray *recpSPD = fDetTypeRec->ClustersAddress(kSPD);
208       Info("InitRec","idt = %#x; digSPD = %#x; recpSPD = %#x",fDetTypeRec,digSPD,recpSPD);
209       AliITSClusterFinderSPD *recSPD = new AliITSClusterFinderSPD(segSPD,digSPD,recpSPD);
210       fDetTypeRec->SetReconstructionModel(kSPD,recSPD);
211     } // end if fDet[kSPD].
212   // SDD
213   if(fDet[kSDD]){
214     Info("InitRec","SDD");
215     //    idt = fITS->DetType(kSDD);
216     AliITSsegmentationSDD *segSDD = (AliITSsegmentationSDD*)
217       fDetTypeRec->GetSegmentationModel(1);
218     AliITSresponseSDD *resSDD = (AliITSresponseSDD*)
219       fDetTypeRec->GetCalibrationModel(fDetTypeRec->GetITSgeom()->GetStartSDD()); 
220     TClonesArray *digSDD = fDetTypeRec->DigitsAddress(kSDD);
221     TClonesArray *recpSDD = fDetTypeRec->ClustersAddress(kSDD);
222     AliITSClusterFinderSDD *recSDD =new AliITSClusterFinderSDD(segSDD,
223                                                                resSDD,
224                                                                digSDD,recpSDD);
225     fDetTypeRec->SetReconstructionModel(kSDD,recSDD);
226   } // end if fDet[kSDD]
227     // SSD
228   if(fDet[kSSD]){
229     Info("InitRec","SSD");
230     //idt = fITS->DetType(kSSD);
231     AliITSsegmentationSSD *segSSD = (AliITSsegmentationSSD*)
232                                        fDetTypeRec->GetSegmentationModel(2);
233       TClonesArray *digSSD = fDetTypeRec->DigitsAddress(kSSD);
234       AliITSClusterFinderSSD *recSSD =new AliITSClusterFinderSSD(segSSD,
235                                                                  digSSD);
236       recSSD->SetITSgeom(fDetTypeRec->GetITSgeom());
237       fDetTypeRec->SetReconstructionModel(kSSD,recSSD);
238     } // end if fDet[kSSD]
239   */
240   fDetTypeRec->SetDefaultClusterFinders();
241     Info("InitRec","    Done\n");
242     return kTRUE;
243 }
244 //______________________________________________________________________ 
245 void AliITSreconstruction::Exec(const Option_t *opt){
246     // Main reconstruction function.
247     // Inputs:
248     //      Option_t * opt   list of subdetector to digitize. =0 all.
249     // Outputs:
250     //      none.
251     // Return:
252     //      none.
253     Option_t *lopt;
254     Int_t evnt;
255
256     if(strstr(opt,"All")||strstr(opt,"ALL")||strstr(opt,"ITS")||opt==0){
257       fDet[0] = fDet[1] = fDet[2] = kTRUE;
258       lopt = "All";
259     }else{
260       fDet[0] = fDet[1] = fDet[2] = kFALSE;
261       if(strstr(opt,"SPD")) fDet[kSPD] = kTRUE;
262       if(strstr(opt,"SDD")) fDet[kSDD] = kTRUE;
263       if(strstr(opt,"SSD")) fDet[kSSD] = kTRUE;
264       if(fDet[kSPD] && fDet[kSDD] && fDet[kSSD]) lopt = "All";
265       else lopt = opt;
266     } // end if strstr(opt,...)
267
268     if(!fInit){
269       cout << "Initilization Failed, Can't run Exec." << endl;
270       return;
271     } // end if !fInit
272     for(evnt=0;evnt<fEnt;evnt++)
273      {
274       Info("Exec","");
275       Info("Exec","Processing Event %d",evnt);
276       Info("Exec","");
277
278       fRunLoader->GetEvent(evnt);
279       if (fLoader->TreeR() == 0x0) fLoader->MakeTree("R");
280       fDetTypeRec->MakeBranchR(0);
281       if (fLoader->TreeC() == 0x0){
282         fDetTypeRec->MakeTreeC();
283         fDetTypeRec->MakeBranchC();
284       }
285       fDetTypeRec->SetTreeAddressR(fLoader->TreeR());
286       fDetTypeRec->SetTreeAddressD(fLoader->TreeD());
287       fDetTypeRec->DigitsToRecPoints(evnt,0,lopt);
288     } // end for evnt
289 }
290 //______________________________________________________________________ 
291 void AliITSreconstruction::SetOutputFile(TString filename){
292   // Set a new file name for recpoints. 
293   // It must be called before Init()
294   if(!fLoader)fLoader = (AliITSLoader*) fRunLoader->GetLoader("ITSLoader");
295   if(fLoader){
296     Info("SetOutputFile","name for rec points is %s",filename.Data());
297     fLoader->SetRecPointsFileName(filename);
298   }
299   else {
300     Error("SetOutputFile",
301     "ITS loader not available. Not possible to set name: %s",filename.Data());
302   }
303 }