]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDigitizer.cxx
New SDigit version of ITS Digitizer.
[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.1  2001/11/27 16:27:28  nilsen
19 Adding AliITSDigitizer class to do merging and digitization . Based on the
20 TTask method. AliITSDigitizer class added to the Makefile and ITSLinkDef.h
21 file. The following files required minor changes. AliITS, added functions
22 SetHitsAddressBranch, MakeBranchInTreeD and modified MakeBranchD.
23 AliITSsimulationSDD.cxx needed a Tree indepenent way of returning back to
24 the original Root Directory in function Compress1D. Now it uses gDirectory.
25
26 Revision 1.2  2002/03/01  E. Lopez
27 Diditization changed to start from SDigits instead of Hits.
28 The SDigits are reading as TClonesArray of AliITSpListItem
29 */
30
31 #include <stdlib.h>
32 #include <iostream.h>
33 #include <TObjArray.h>
34 #include <TTree.h>
35 #include <TBranch.h>
36
37 #include <AliRun.h>
38 #include <AliRunDigitizer.h>
39
40 #include "AliITSDigitizer.h"
41 #include "AliITShit.h"
42 #include "AliITSmodule.h"
43 #include "AliITSsimulation.h"
44 #include "AliITSDetType.h"
45 #include "AliITSgeom.h"
46
47 ClassImp(AliITSDigitizer)
48
49 //______________________________________________________________________
50 AliITSDigitizer::AliITSDigitizer() : AliDigitizer(){
51     // Default constructor. Assign fITS since it is never written out from
52     // here. 
53     // Inputs:
54     //      Option_t * opt   Not used
55     // Outputs:
56     //      none.
57     // Return:
58     //      A blank AliITSDigitizer class.
59
60     fITS    = 0;
61     fActive = 0;
62 }
63 //______________________________________________________________________
64 AliITSDigitizer::AliITSDigitizer(AliRunDigitizer *mngr) : AliDigitizer(mngr){
65     // Standard constructor. Assign fITS since it is never written out from
66     // here. 
67     // Inputs:
68     //      Option_t * opt   Not used
69     // Outputs:
70     //      none.
71     // Return:
72     //      An AliItSDigitizer class.
73
74     if(!gAlice) {
75         fITS = 0;
76         fActive = 0;
77         return;
78     } // end if
79     fITS = (AliITS *)(gAlice->GetDetector("ITS"));
80     if(!fITS){
81         fActive = 0;
82         return;
83     } else if(fITS->GetITSgeom()){
84         fActive = new Bool_t[fITS->GetITSgeom()->GetIndexMax()];
85     } else{
86         fActive = 0;
87         return;
88     } // end if
89     // fActive needs to be set to a default all kTRUE value
90     for(Int_t i=0;i<fITS->GetITSgeom()->GetIndexMax();i++) fActive[i] = kTRUE;
91 }
92 //______________________________________________________________________
93 AliITSDigitizer::~AliITSDigitizer(){
94     // Default destructor. 
95     // Inputs:
96     //      Option_t * opt   Not used
97     // Outputs:
98     //      none.
99     // Return:
100     //      none.
101
102     fITS = 0; // don't delete fITS. Done else where.
103     if(fActive) delete[] fActive;
104 }
105 //______________________________________________________________________
106 Bool_t AliITSDigitizer::Init(){
107     // Iniliztion 
108     // Inputs:
109     //      none.
110     // Outputs:
111     //      none.
112     // Return:
113     //      none.
114     return kTRUE;
115 }
116 //______________________________________________________________________
117 void AliITSDigitizer::Exec(Option_t* opt){
118     // Main digitizing function. 
119     // Inputs:
120     //      Option_t * opt   list of subdetector to digitize. =0 all.
121     // Outputs:
122     //      none.
123     // Return:
124     //      none.
125
126     char name[20] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
127     char *all;
128     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
129                           strstr(opt,"SSD")};
130     if( !det[0] && !det[1] && !det[2] ) all = "All";
131     else all = 0;
132     AliITSsimulation *sim      = 0;
133     AliITSDetType    *iDetType = 0;
134     static Bool_t    setDef    = kTRUE;
135
136     if( !fITS ) fITS = (AliITS*)(gAlice->GetDetector("ITS"));
137     if(!fITS){
138         Error("Exec","The ITS not found. aborting.");
139         return;
140     } // end if
141     if( !(fITS->GetITSgeom()) ) {
142         Warning( "Exec", "Need ITS geometry to be properly defined first." );
143         return; // need transformations to do digitization.
144     } // end if !GetITSgeom()
145
146     if( setDef ) fITS->SetDefaultSimulation();
147     setDef = kFALSE;
148     sprintf( name, "%s", fITS->GetName() );
149
150     Int_t nfiles = GetManager()->GetNinputs();
151     Int_t event  = GetManager()->GetOutputEventNr();
152     Int_t size   = fITS->GetITSgeom()->GetIndexMax();
153     TClonesArray * sdig = new TClonesArray( "AliITSpListItem",1000 );
154     
155     // Digitize 
156     fITS->MakeBranchInTreeD( GetManager()->GetTreeD() );
157
158     for( Int_t module=0; module<size; module++ ){
159         if(fActive) if(!fActive[module]) continue;
160         Int_t id = fITS->GetITSgeom()->GetModuleType( module );
161         if( !all && !det[id] ) continue;
162         iDetType = fITS->DetType( id );
163         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
164         if( !sim ) {
165             Error( "Exec", "The simulation class was not instanciated!" );
166             exit(1);
167         } // end if !sim
168         
169         // Fill the module with the sum of SDigits
170         sim->InitSimulationModule( module, event );
171         cout << "Module=" << module;
172         for( Int_t ifiles=0; ifiles<nfiles; ifiles++ ){
173             cout <<" ifiles=" << ifiles;
174             TTree *treeS = GetManager()->GetInputTreeS( ifiles );
175             if( !(treeS && fITS->GetSDigits()) ) continue;   
176             TBranch *brchSDigits = treeS->GetBranch( name );
177             if( brchSDigits ) {
178                 brchSDigits->SetAddress( &sdig ); 
179             } else {
180                 Error( "Exec", "branch ITS not found in TreeS, input file %d ",
181                        ifiles );
182             } // end if brchSDigits
183             sdig->Clear();
184             Int_t mask = GetManager()->GetMask( ifiles );
185             
186             // add summable digits to module
187             brchSDigits->GetEvent( module );
188             sim->AddSDigitsToModule( sdig, mask );    
189         } // end for ifiles
190         cout << " end ifiles loop" << endl;
191         // Digitise current module sum(SDigits)->Digits
192         sim->FinishSDigitiseModule();
193
194         // fills all branches - wasted disk space
195         GetManager()->GetTreeD()->Fill();
196         fITS->ResetDigits();
197     } // end for module
198     cout << "end modules loop"<<endl;
199  
200     GetManager()->GetTreeD()->GetEntries();
201     GetManager()->GetTreeD()->Write( 0, TObject::kOverwrite );
202     // reset tree
203     GetManager()->GetTreeD()->Reset();
204     
205     sdig->Clear();
206     delete sdig;
207 }