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