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