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