]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSreconstruction.cxx
Bug fix. SSD calibration objects were overwriting the SDD ones in case of partial...
[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();
165     fDetTypeRec->SetITSgeom(fITSgeom);
166     fDetTypeRec->SetDefaults();
167     fDet[0] = fDet[1] = fDet[2] = kTRUE;
168     fEnt0 = 0;
169
170     //fEnt  = gAlice->GetEventsPerRun();
171     fEnt = Int_t(fRunLoader->GetNumberOfEvents());
172
173     fLoader->LoadDigits("read");
174     fLoader->LoadRecPoints("recreate");
175     if (fLoader->TreeR() == 0x0) fLoader->MakeTree("R");
176  
177     fDetTypeRec->SetTreeAddressD(fLoader->TreeD());
178     fDetTypeRec->MakeBranchR(fLoader->TreeR());
179     fDetTypeRec->SetTreeAddressR(fLoader->TreeR());
180
181     fInit = InitRec();
182
183     Info("Init","  Done\n\n\n");
184
185     return fInit;
186 }
187 //______________________________________________________________________
188 Bool_t AliITSreconstruction::InitRec(){
189     // Sets up Reconstruction part of AliITSDetType..
190     // Inputs:
191     //      none.
192     // Outputs:
193     //      none.
194     // Return:
195     //      none.
196   /*
197   //AliITSDetType *idt;
198   fDetTypeRec->SetLoader(fLoader);
199     // SPD
200   if(fDet[kSPD]){
201     Info("InitRec","SPD");
202     //idt = fITS->DetType(kSPD);
203     AliITSsegmentationSPD *segSPD = (AliITSsegmentationSPD*)fDetTypeRec->GetSegmentationModel(0);
204       TClonesArray *digSPD = fDetTypeRec->DigitsAddress(kSPD);
205       TClonesArray *recpSPD = fDetTypeRec->ClustersAddress(kSPD);
206       Info("InitRec","idt = %#x; digSPD = %#x; recpSPD = %#x",fDetTypeRec,digSPD,recpSPD);
207       AliITSClusterFinderSPD *recSPD = new AliITSClusterFinderSPD(segSPD,digSPD,recpSPD);
208       fDetTypeRec->SetReconstructionModel(kSPD,recSPD);
209     } // end if fDet[kSPD].
210   // SDD
211   if(fDet[kSDD]){
212     Info("InitRec","SDD");
213     //    idt = fITS->DetType(kSDD);
214     AliITSsegmentationSDD *segSDD = (AliITSsegmentationSDD*)
215       fDetTypeRec->GetSegmentationModel(1);
216     AliITSresponseSDD *resSDD = (AliITSresponseSDD*)
217       fDetTypeRec->GetCalibrationModel(fDetTypeRec->GetITSgeom()->GetStartSDD()); 
218     TClonesArray *digSDD = fDetTypeRec->DigitsAddress(kSDD);
219     TClonesArray *recpSDD = fDetTypeRec->ClustersAddress(kSDD);
220     AliITSClusterFinderSDD *recSDD =new AliITSClusterFinderSDD(segSDD,
221                                                                resSDD,
222                                                                digSDD,recpSDD);
223     fDetTypeRec->SetReconstructionModel(kSDD,recSDD);
224   } // end if fDet[kSDD]
225     // SSD
226   if(fDet[kSSD]){
227     Info("InitRec","SSD");
228     //idt = fITS->DetType(kSSD);
229     AliITSsegmentationSSD *segSSD = (AliITSsegmentationSSD*)
230                                        fDetTypeRec->GetSegmentationModel(2);
231       TClonesArray *digSSD = fDetTypeRec->DigitsAddress(kSSD);
232       AliITSClusterFinderSSD *recSSD =new AliITSClusterFinderSSD(segSSD,
233                                                                  digSSD);
234       recSSD->SetITSgeom(fDetTypeRec->GetITSgeom());
235       fDetTypeRec->SetReconstructionModel(kSSD,recSSD);
236     } // end if fDet[kSSD]
237   */
238   fDetTypeRec->SetDefaultClusterFinders();
239     Info("InitRec","    Done\n");
240     return kTRUE;
241 }
242 //______________________________________________________________________ 
243 void AliITSreconstruction::Exec(const Option_t *opt){
244     // Main reconstruction function.
245     // Inputs:
246     //      Option_t * opt   list of subdetector to digitize. =0 all.
247     // Outputs:
248     //      none.
249     // Return:
250     //      none.
251     Option_t *lopt;
252     Int_t evnt;
253
254     if(strstr(opt,"All")||strstr(opt,"ALL")||strstr(opt,"ITS")||opt==0){
255       fDet[0] = fDet[1] = fDet[2] = kTRUE;
256       lopt = "All";
257     }else{
258       fDet[0] = fDet[1] = fDet[2] = kFALSE;
259       if(strstr(opt,"SPD")) fDet[kSPD] = kTRUE;
260       if(strstr(opt,"SDD")) fDet[kSDD] = kTRUE;
261       if(strstr(opt,"SSD")) fDet[kSSD] = kTRUE;
262       if(fDet[kSPD] && fDet[kSDD] && fDet[kSSD]) lopt = "All";
263       else lopt = opt;
264     } // end if strstr(opt,...)
265
266     if(!fInit){
267       cout << "Initilization Failed, Can't run Exec." << endl;
268       return;
269     } // end if !fInit
270     for(evnt=0;evnt<fEnt;evnt++)
271      {
272       Info("Exec","");
273       Info("Exec","Processing Event %d",evnt);
274       Info("Exec","");
275
276       fRunLoader->GetEvent(evnt);
277       if (fLoader->TreeR() == 0x0) fLoader->MakeTree("R");
278       fDetTypeRec->MakeBranchR(0);
279       fDetTypeRec->SetTreeAddressR(fLoader->TreeR());
280       fDetTypeRec->SetTreeAddressD(fLoader->TreeD());
281       fDetTypeRec->DigitsToRecPoints(fLoader->TreeD(),fLoader->TreeR(),0,lopt);
282     } // end for evnt
283 }
284 //______________________________________________________________________ 
285 void AliITSreconstruction::SetOutputFile(TString filename){
286   // Set a new file name for recpoints. 
287   // It must be called before Init()
288   if(!fLoader)fLoader = (AliITSLoader*) fRunLoader->GetLoader("ITSLoader");
289   if(fLoader){
290     Info("SetOutputFile","name for rec points is %s",filename.Data());
291     fLoader->SetRecPointsFileName(filename);
292   }
293   else {
294     Error("SetOutputFile",
295     "ITS loader not available. Not possible to set name: %s",filename.Data());
296   }
297 }