]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDigitizer.cxx
New ITS SDigit merging with region of interest cut. Update for changes in
[u/mrichter/AliRoot.git] / ITS / AliITSDigitizer.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.2  2002/03/15 17:26:40  nilsen
19 New SDigit version of ITS Digitizer.
20
21 Revision 1.1  2001/11/27 16:27:28  nilsen
22 Adding AliITSDigitizer class to do merging and digitization . Based on the
23 TTask method. AliITSDigitizer class added to the Makefile and ITSLinkDef.h
24 file. The following files required minor changes. AliITS, added functions
25 SetHitsAddressBranch, MakeBranchInTreeD and modified MakeBranchD.
26 AliITSsimulationSDD.cxx needed a Tree independent way of returning back to
27 the original Root Directory in function Compress1D. Now it uses gDirectory.
28
29 Revision 1.2  2002/03/01  E. Lopez
30 Digitization changed to start from SDigits instead of Hits.
31 The SDigits are reading as TClonesArray of AliITSpListItem
32 */
33
34 #include <stdlib.h>
35 #include <iostream.h>
36 #include <TObjArray.h>
37 #include <TTree.h>
38 #include <TBranch.h>
39 #include <TFile.h>
40
41 #include <AliRun.h>
42 #include <AliRunDigitizer.h>
43
44 #include "AliITSDigitizer.h"
45 #include "AliITShit.h"
46 #include "AliITSmodule.h"
47 #include "AliITSsimulation.h"
48 #include "AliITSDetType.h"
49 #include "AliITSgeom.h"
50
51 ClassImp(AliITSDigitizer)
52
53 //______________________________________________________________________
54 AliITSDigitizer::AliITSDigitizer() : AliDigitizer(){
55     // Default constructor. Assign fITS since it is never written out from
56     // here. 
57     // Inputs:
58     //      Option_t * opt   Not used
59     // Outputs:
60     //      none.
61     // Return:
62     //      A blank AliITSDigitizer class.
63
64     fITS    = 0;
65     fActive = 0;
66     fRoif   = -1;
67     fRoiifile = 0;
68 }
69 //______________________________________________________________________
70 AliITSDigitizer::AliITSDigitizer(AliRunDigitizer *mngr) : AliDigitizer(mngr){
71     // Standard constructor. Assign fITS since it is never written out from
72     // here. 
73     // Inputs:
74     //      Option_t * opt   Not used
75     // Outputs:
76     //      none.
77     // Return:
78     //      An AliItSDigitizer class.
79
80     if(!gAlice) {
81         fITS = 0;
82         fActive = 0;
83         fRoif   = -1;
84         fRoiifile = 0;
85         return;
86     } // end if
87     fITS = (AliITS *)(gAlice->GetDetector("ITS"));
88     if(!fITS){
89         fActive = 0;
90         fRoif   = -1;
91         fRoiifile = 0;
92         return;
93     } else if(fITS->GetITSgeom()){
94         fActive = new Bool_t[fITS->GetITSgeom()->GetIndexMax()];
95         fRoif   = -1;
96         fRoiifile = 0;
97     } else{
98         fActive = 0;
99         fRoif   = -1;
100         fRoiifile = 0;
101         return;
102     } // end if
103     // fActive needs to be set to a default all kTRUE value
104     for(Int_t i=0;i<fITS->GetITSgeom()->GetIndexMax();i++) fActive[i] = kTRUE;
105 }
106 //______________________________________________________________________
107 AliITSDigitizer::~AliITSDigitizer(){
108     // Default destructor. 
109     // Inputs:
110     //      Option_t * opt   Not used
111     // Outputs:
112     //      none.
113     // Return:
114     //      none.
115
116     fITS = 0; // don't delete fITS. Done else where.
117     if(fActive) delete[] fActive;
118 }
119 //______________________________________________________________________
120 Bool_t AliITSDigitizer::Init(const char *filename){
121     // Initialization 
122     // Inputs:
123     //      none.
124     // Outputs:
125     //      none.
126     // Return:
127     //      none.
128      
129 //    TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
130 //    if (file) {file->Close(); delete file;}
131     TFile *file = new TFile(filename);
132     if (!file->IsOpen()) {
133         Error("Init","Can't open %s!",filename);
134         return kFALSE;
135     } // end if !file 
136     if (gAlice) delete gAlice;
137     gAlice = (AliRun*)file->Get("gAlice");
138     if (!gAlice) {
139         Error("Init","AliRun object not found on file");
140         return kFALSE;
141     } // end if !gAlice
142
143     this->Init();
144     delete gAlice;
145     gAlice=0;
146     return kTRUE;
147 }
148 //______________________________________________________________________
149 Bool_t AliITSDigitizer::Init(){
150     // Initialization 
151     // Inputs:
152     //      none.
153     // Outputs:
154     //      none.
155     // Return:
156
157     //      none.
158 //    TTree *ts;
159
160     if(!gAlice) {
161         fITS      = 0;
162         fActive   = 0;
163         fRoif     = -1;
164         fRoiifile = 0;
165 //      cout << "init !gAlice"<<endl;
166         return kFALSE;
167     } // end if
168     fITS = (AliITS *)(gAlice->GetDetector("ITS"));
169     if(!fITS){
170         fActive   = 0;
171         fRoif     = -1;
172         fRoiifile = 0;
173 //      cout << "init !fITS"<<endl;
174         return kFALSE;
175     } else if(fITS->GetITSgeom()){
176         fActive = new Bool_t[fITS->GetITSgeom()->GetIndexMax()];
177     } else{
178         fActive = 0;
179         fRoif   = -1;
180         fRoiifile = 0;
181 //      cout << "init !fITS-GetITSgeom()"<<endl;
182         return kFALSE;
183     } // end if
184     // fActive needs to be set to a default all kTRUE value
185     for(Int_t i=0;i<fITS->GetITSgeom()->GetIndexMax();i++) fActive[i] = kTRUE;
186 /*
187     if(fRoif>=0 && fRoiifile>=0 && fRoiifile<GetManager()->GetNinputs()){
188         ts = GetManager()->GetInputTreeS(fRoiifile);
189         cout << "calling SetByReionOfInterest"<<endl;
190         SetByRegionOfInterest(ts);
191     } // end if
192 */
193     return kTRUE;
194 }
195 //______________________________________________________________________
196 void AliITSDigitizer::Exec(Option_t* opt){
197     // Main digitization function. 
198     // Inputs:
199     //      Option_t * opt   list of sub detector to digitize. =0 all.
200     // Outputs:
201     //      none.
202     // Return:
203     //      none.
204
205     char name[20] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
206     char *all;
207     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
208                           strstr(opt,"SSD")};
209     if( !det[0] && !det[1] && !det[2] ) all = "All";
210     else all = 0;
211     AliITSsimulation *sim      = 0;
212     AliITSDetType    *iDetType = 0;
213     static Bool_t    setDef    = kTRUE;
214
215     if(!fITS) fITS = (AliITS*)(gAlice->GetDetector("ITS"));
216     if(!fITS){
217         Error("Exec","The ITS not found. aborting.");
218         return;
219     } // end if
220     if( !(fITS->GetITSgeom()) ) {
221         Warning( "Exec", "Need ITS geometry to be properly defined first." );
222         return; // need transformations to do digitization.
223     } // end if !GetITSgeom()
224
225     if( setDef ) fITS->SetDefaultSimulation();
226     setDef = kFALSE;
227     sprintf( name, "%s", fITS->GetName() );
228
229     Int_t nfiles = GetManager()->GetNinputs();
230     Int_t event  = GetManager()->GetOutputEventNr();
231     Int_t size   = fITS->GetITSgeom()->GetIndexMax();
232     TClonesArray * sdig = new TClonesArray( "AliITSpListItem",1000 );
233     
234     // Digitize 
235     fITS->MakeBranchInTreeD( GetManager()->GetTreeD() );
236     if(fRoif!=-1)SetByRegionOfInterest(GetManager()->GetInputTreeS(fRoiifile));
237
238     for( Int_t module=0; module<size; module++ ){
239         if(fActive && fRoif!=-1) if(!fActive[module]) continue;
240         Int_t id = fITS->GetITSgeom()->GetModuleType( module );
241         if( !all && !det[id] ) continue;
242         iDetType = fITS->DetType( id );
243         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
244         if(!sim) {
245             Error( "Exec", "The simulation class was not instanciated!" );
246             exit(1);
247         } // end if !sim
248
249         // Fill the module with the sum of SDigits
250         sim->InitSimulationModule( module, event );
251 //      cout << "Module=" << module;
252         for( Int_t ifiles=0; ifiles<nfiles; ifiles++ ){
253 //          cout <<" ifiles=" << ifiles;
254             TTree *treeS = GetManager()->GetInputTreeS(ifiles);
255             if( !(treeS && fITS->GetSDigits()) ) continue;   
256             TBranch *brchSDigits = treeS->GetBranch( name );
257             if( brchSDigits ) {
258                 brchSDigits->SetAddress( &sdig ); 
259             } else {
260                 Error( "Exec", "branch ITS not found in TreeS, input file %d ",
261                        ifiles );
262             } // end if brchSDigits
263             sdig->Clear();
264             Int_t mask = GetManager()->GetMask( ifiles );
265             
266             // add summable digits to module
267             brchSDigits->GetEvent( module );
268             sim->AddSDigitsToModule( sdig, mask );    
269         } // end for ifiles
270 //      cout << " end ifiles loop" << endl;
271         // Digitize current module sum(SDigits)->Digits
272         sim->FinishSDigitiseModule();
273
274         // fills all branches - wasted disk space
275         GetManager()->GetTreeD()->Fill();
276         fITS->ResetDigits();
277     } // end for module
278 //    cout << "end modules loop"<<endl;
279  
280     GetManager()->GetTreeD()->GetEntries();
281     GetManager()->GetTreeD()->Write( 0, TObject::kOverwrite );
282     // reset tree
283     GetManager()->GetTreeD()->Reset();
284     
285     sdig->Clear();
286     delete sdig;
287 }
288 //______________________________________________________________________
289 void AliITSDigitizer::SetByRegionOfInterest(TTree *ts){
290     // Scans through the ITS branch of the SDigits tree, ts, for modules
291     // which have SDigits in them. For these modules, a flag is set to
292     // digitize only these modules. The value of fRoif determines how many
293     // neighboring modules will also be turned on. fRoif=0 will turn on only
294     // those modules with SDigits in them. fRoif=1 will turn on, in addition,
295     // those modules that are +-1 module from the one with the SDigits. And
296     // So on. This last feature is not supported yet.
297     // Inputs:
298     //      TTree *ts  The tree in which the existing SDigits will define the
299     //                 region of interest.
300     // Outputs:
301     //      none.
302     // Return:
303     //      none.
304     Int_t m,nm;
305     if(gAlice==0 && ts==0) return;
306     if(ts==0) ts = gAlice->TreeS();
307     TBranch *brchSDigits = ts->GetBranch(fITS->GetName());
308     TClonesArray * sdig = new TClonesArray( "AliITSpListItem",1000 );
309 //    cout << "Region of Interest ts="<<ts<<" brchSDigits="<<brchSDigits<<" sdig="<<sdig<<endl;
310
311     if( brchSDigits ) {
312         brchSDigits->SetAddress( &sdig ); 
313     } else {
314         Error( "SetByRegionOfInterest","branch ITS not found in TreeS");
315         return;
316     } // end if brchSDigits
317
318     nm = fITS->GetITSgeom()->GetIndexMax();
319     for(m=0;m<nm;m++){
320 //      cout << " fActive["<<m<<"]=";
321         fActive[m] = kFALSE; // Not active by default
322         sdig->Clear();
323         brchSDigits->GetEvent(m);
324         if(sdig->GetLast()>=0){ // activate the necessary modules
325             if(fRoif>=0){
326                 fActive[m] = kTRUE;
327 //              cout << "kTRUE";
328             } // end if
329         } // end if.
330 //      cout << endl;
331     } // end for m
332     sdig->Clear();
333     delete sdig;
334 }