1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.5 2002/10/14 14:57:00 hristov
19 Merging the VirtualMC branch to the main development branch (HEAD)
21 Revision 1.3.4.1 2002/06/10 17:51:14 hristov
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).
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.
32 Revision 1.2 2002/03/15 17:26:40 nilsen
33 New SDigit version of ITS Digitizer.
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.
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
49 #include <Riostream.h>
50 #include <TObjArray.h>
56 #include <AliRunDigitizer.h>
58 #include "AliITSDigitizer.h"
59 #include "AliITSpList.h"
60 #include "AliITSmodule.h"
61 #include "AliITSsimulation.h"
62 #include "AliITSDetType.h"
63 #include "AliITSgeom.h"
65 ClassImp(AliITSDigitizer)
67 //______________________________________________________________________
68 AliITSDigitizer::AliITSDigitizer() : AliDigitizer(){
69 // Default constructor. Assign fITS since it is never written out from
72 // Option_t * opt Not used
76 // A blank AliITSDigitizer class.
84 //______________________________________________________________________
85 AliITSDigitizer::AliITSDigitizer(AliRunDigitizer *mngr) : AliDigitizer(mngr){
86 // Standard constructor. Assign fITS since it is never written out from
89 // Option_t * opt Not used
93 // An AliItSDigitizer class.
101 //______________________________________________________________________
102 AliITSDigitizer::~AliITSDigitizer(){
103 // Default destructor.
105 // Option_t * opt Not used
111 fITS = 0; // don't delete fITS. Done else where.
112 if(fActive) delete[] fActive;
114 //______________________________________________________________________
115 Bool_t AliITSDigitizer::Init(){
116 // Initialization. Set up region of interest, if switched on, and
117 // loads ITS and ITSgeom.
125 fInit = kTRUE; // Assume for now init will work.
132 Warning("Init","gAlice not found");
135 fITS = (AliITS *)(gAlice->GetDetector("ITS"));
141 Warning("Init","ITS not found");
143 } else if(fITS->GetITSgeom()){
144 //cout << "fRoif,fRoiifile="<<fRoif<<" "<<fRoiifile<<endl;
145 fActive = new Bool_t[fITS->GetITSgeom()->GetIndexMax()];
151 Warning("Init","ITS geometry not found");
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.
158 if(fRoif>=0 && fRoiifile>=0 && fRoiifile<GetManager()->GetNinputs()){
159 ts = GetManager()->GetInputTreeS(fRoiifile);
161 if(!gAlice) ts = gAlice->TreeS();
163 cout <<"The TTree TreeS needed to set by region not found."
164 " No region of interest cut will be applied."<< endl;
168 cout << "calling SetByReionOfInterest ts="<< ts <<endl;
169 SetByRegionOfInterest(ts);
174 //______________________________________________________________________
175 void AliITSDigitizer::Exec(Option_t* opt){
176 // Main digitization function.
178 // Option_t * opt list of sub detector to digitize. =0 all.
184 char name[20] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
186 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
188 if( !det[0] && !det[1] && !det[2] ) all = "All";
190 AliITSsimulation *sim = 0;
191 AliITSDetType *iDetType = 0;
192 static Bool_t setDef = kTRUE;
195 Error("Exec","Init not succesfull, aborting.");
199 if( setDef ) fITS->SetDefaultSimulation();
201 sprintf(name,"%s",fITS->GetName());
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;
208 Int_t *fl = new Int_t[nfiles];
211 for(id=0;id<nfiles;id++) if(id!=fRoiifile){
212 // just in case fRoiifile!=0.
216 TClonesArray * sdig = new TClonesArray( "AliITSpListItem",1000 );
219 fITS->MakeBranchInTreeD(GetManager()->GetTreeD());
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();
228 Error( "Exec", "The simulation class was not instanciated!" );
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 );
242 brchSDigits->SetAddress( &sdig );
244 Error( "Exec", "branch ITS not found in TreeS, input file %d ",
247 } // end if brchSDigits
249 mask = GetManager()->GetMask(ifiles);
250 // add summable digits to module
251 brchSDigits->GetEvent( module );
252 lmod = sim->AddSDigitsToModule(sdig,mask);
254 fActive[module] = lmod;
256 //cout << " fActive["<<module<<"]=";
257 //if(fActive[module]) cout << "kTRUE";
258 //else cout << "kFALSE";
260 //cout << " end ifiles loop" << endl;
261 // Digitize current module sum(SDigits)->Digits
262 sim->FinishSDigitiseModule();
264 // fills all branches - wasted disk space
265 GetManager()->GetTreeD()->Fill();
268 //cout << "end modules loop"<<endl;
270 GetManager()->GetTreeD()->AutoSave();
275 for(Int_t i=0;i<fITS->GetITSgeom()->GetIndexMax();i++) fActive[i] = kTRUE;
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.
288 // TTree *ts The tree in which the existing SDigits will define the
289 // region of interest.
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;
303 brchSDigits->SetAddress( &sdig );
305 Error( "SetByRegionOfInterest","branch ITS not found in TreeS");
307 } // end if brchSDigits
309 nm = fITS->GetITSgeom()->GetIndexMax();
311 //cout << " fActive["<<m<<"]=";
312 fActive[m] = kFALSE; // Not active by default
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.
321 } // end if. end for i.
322 //cout << fActive[m];