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