]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSreconstruction.cxx
Changes by Massimo Masera to allow Recpoints and Clusters to be written
[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 /*
17 $Log$
18 Revision 1.3  2002/02/06 13:52:27  barbera
19 gAlice deletion corrected (from M. Masera)
20
21 Revision 1.2  2002/01/31 18:52:09  nilsen
22 Minor change to allow the use of files that are already open. grun.C macro
23 that also does ITS digitizationa and Reconstruction all in one go.
24
25 Revision 1.1  2002/01/30 22:20:22  nilsen
26 New TTask based method to do Digits To clusters. Works with files of multiple
27 events in the file. Macro added to show how to use. Also the changes made
28 in the nessesary complilation files.
29
30 */
31 #include <TROOT.h>
32 #include <TFile.h>
33 #include <TSeqCollection.h>
34 #include <TString.h>
35 #include <TClonesArray.h>
36
37 #include "AliRun.h"
38
39 #include "AliITS.h"
40 #include "AliITSDetType.h"
41 #include "AliITSreconstruction.h"
42 #include "AliITSsegmentationSPD.h"
43 #include "AliITSsegmentationSDD.h"
44 #include "AliITSsegmentationSSD.h"
45 #include "AliITSClusterFinderSPD.h"
46 #include "AliITSClusterFinderSDD.h"
47 #include "AliITSClusterFinderSSD.h"
48 #include "AliITSresponseSDD.h"
49 #include "AliITSgeom.h"
50
51 ClassImp(AliITSreconstruction)
52
53 //______________________________________________________________________
54 AliITSreconstruction::AliITSreconstruction(){
55     // Default constructor.
56     // Inputs:
57     //  none.
58     // Outputs:
59     //   none.
60     // Return:
61     //    A zero-ed constructed AliITSreconstruction class.
62
63     fFilename = "";
64     fFile     = 0;
65     fFile2    = 0;
66     fITS      = 0;
67     fDet[0] = fDet[1] = fDet[2] = kTRUE;
68     fInit     = kFALSE;
69 }
70 //______________________________________________________________________
71 AliITSreconstruction::AliITSreconstruction(const char* filename){
72     // Standard constructor.
73     // Inputs:
74     //  const char* filename    filename containing the digits to be
75     //                          reconstructed. If filename = 0 (nil)
76     //                          then no file is opened but a file is
77     //                          assumed to already be opened. This 
78     //                          already opened file will be used.
79     // Outputs:
80     //   none.
81     // Return:
82     //    A standardly constructed AliITSreconstruction class.
83
84     fFilename = filename;
85     fFile2 = 0;
86     fFile = 0;
87     if(filename){
88         fFile = (TFile*)gROOT->GetListOfFiles()->FindObject(fFilename.Data());
89         if(fFile) fFile->Close();
90         fFile = new TFile(fFilename.Data(),"UPDATE");
91     if(gAlice) {
92       delete gAlice;
93       gAlice = 0;
94     }
95         gAlice = (AliRun*)fFile->Get("gAlice");
96         if(!gAlice) {
97             cout << "gAlice not found on file. Aborting." << endl;
98             fInit = kFALSE;
99             return;
100         } // end if !gAlice
101     } // end if !filename.
102 }
103 //______________________________________________________________________
104 AliITSreconstruction::~AliITSreconstruction(){
105     // Default constructor.
106     // Inputs:
107     //  none.
108     // Outputs:
109     //   none.
110     // Return:
111     //    A destroyed AliITSreconstruction class.
112
113     if(fFile) fFile->Close();
114     fFile     = 0;
115     fITS      = 0;
116     
117 }
118 //______________________________________________________________________
119 Bool_t AliITSreconstruction::Init(){
120     // Class Initilizer.
121     // Inputs:
122     //  none.
123     // Outputs:
124     //   none.
125     // Return:
126     //    kTRUE if no errors initilizing this class occurse else kFALSE
127     Int_t nparticles;
128     fITS = (AliITS*) gAlice->GetDetector("ITS");
129     if(!fITS){
130         cout << "ITS not found aborting. fITS=" << fITS << endl;
131         fInit = kFALSE;
132         return fInit;
133     } // end if !fITS
134     if(!(fITS->GetITSgeom())){
135         cout << "ITSgeom not found aborting."<< endl;
136         fInit = kFALSE;
137         return fInit;
138     } // end if !GetITSgeom()
139     // Now ready to init.
140
141     fDet[0] = fDet[1] = fDet[2] = kTRUE;
142     fEnt0 = 0;
143     fEnt  = gAlice->GetEventsPerRun();
144     fITS->MakeTreeC();
145     nparticles = gAlice->GetEvent(fEnt0);
146     
147     // finished init.
148     fInit = InitRec();
149     return fInit;
150 }
151 //______________________________________________________________________
152 Bool_t AliITSreconstruction::InitRec(){
153     // Sets up Reconstruction part of AliITSDetType..
154     // Inputs:
155     //      none.
156     // Outputs:
157     //      none.
158     // Return:
159     //      none.
160     AliITSDetType *idt;
161
162     // SPD
163     if(fDet[kSPD]){
164         idt = fITS->DetType(kSPD);
165         AliITSsegmentationSPD *segSPD = (AliITSsegmentationSPD*)
166                                                idt->GetSegmentationModel();
167         TClonesArray *digSPD = fITS->DigitsAddress(kSPD);
168         TClonesArray *recpSPD = fITS->ClustersAddress(kSPD);
169         AliITSClusterFinderSPD *recSPD = new AliITSClusterFinderSPD(segSPD,
170                                                                     digSPD,
171                                                                     recpSPD);
172         fITS->SetReconstructionModel(kSPD,recSPD);
173     } // end if fDet[kSPD].
174     // SDD
175     if(fDet[kSDD]){
176         idt = fITS->DetType(kSDD);
177         AliITSsegmentationSDD *segSDD = (AliITSsegmentationSDD*)
178                                            idt->GetSegmentationModel();
179         AliITSresponseSDD *resSDD = (AliITSresponseSDD*)
180                                            idt->GetResponseModel();
181         TClonesArray *digSDD = fITS->DigitsAddress(kSDD);
182         TClonesArray *recpSDD = fITS->ClustersAddress(kSDD);
183         AliITSClusterFinderSDD *recSDD =new AliITSClusterFinderSDD(segSDD,
184                                                                    resSDD,
185                                                                digSDD,recpSDD);
186         fITS->SetReconstructionModel(kSDD,recSDD);
187     } // end if fDet[kSDD]
188     // SSD
189     if(fDet[kSSD]){
190         idt = fITS->DetType(kSSD);
191         AliITSsegmentationSSD *segSSD = (AliITSsegmentationSSD*)
192                                          idt->GetSegmentationModel();
193         TClonesArray *digSSD = fITS->DigitsAddress(kSSD);
194         AliITSClusterFinderSSD *recSSD =new AliITSClusterFinderSSD(segSSD,
195                                                                    digSSD);
196         fITS->SetReconstructionModel(kSSD,recSSD);
197     } // end if fDet[kSSD]
198
199     return kTRUE;
200 }
201 //______________________________________________________________________ 
202 void AliITSreconstruction::Exec(const Option_t *opt){
203     // Main reconstruction function.
204     // Inputs:
205     //      Option_t * opt   list of subdetector to digitize. =0 all.
206     // Outputs:
207     //      none.
208     // Return:
209     //      none.
210     Option_t *lopt;
211     Int_t nparticles,evnt;
212
213     if(strstr(opt,"All")||strstr(opt,"ALL")||strstr(opt,"ITS")||opt==0){
214         fDet[0] = fDet[1] = fDet[2] = kTRUE;
215         lopt = "All";
216     }else{
217         fDet[0] = fDet[1] = fDet[2] = kFALSE;
218         if(strstr(opt,"SPD")) fDet[kSPD] = kTRUE;
219         if(strstr(opt,"SDD")) fDet[kSDD] = kTRUE;
220         if(strstr(opt,"SSD")) fDet[kSSD] = kTRUE;
221         if(fDet[kSPD] && fDet[kSDD] && fDet[kSSD]) lopt = "All";
222         else lopt = opt;
223     } // end if strstr(opt,...)
224
225     if(!fInit){
226         cout << "Initilization Failed, Can't run Exec." << endl;
227         return;
228     } // end if !fInit
229     TDirectory *curr = gDirectory;
230     if(fFile2)fFile2->cd();
231     for(evnt=0;evnt<fEnt;evnt++){
232         nparticles = gAlice->GetEvent(evnt);
233         gAlice->SetEvent(evnt);
234         if(!gAlice->TreeR()){
235       if(fFile2){
236         gAlice->MakeTree("R",fFile2);
237       }
238       else {
239         gAlice->MakeTree("R");
240       }
241     }
242         fITS->MakeBranch("R");
243         fITS->DigitsToRecPoints(evnt,0,lopt);
244     cout << "end of recpoints finding *********************\n";
245     } // end for evnt
246     curr->cd();
247 }
248 //______________________________________________________________________ 
249 void AliITSreconstruction::SetOutputFile(TString filename){
250   // Set a file name for recpoints. Used only if this file is not the file
251   // containing digits. This obj is deleted by AliRun.
252   fFile2 = gAlice->InitTreeFile("R",filename);
253 }
254
255
256
257
258