]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITS.cxx
Modifications associated with remerging the Ba/Sa and Dubna pixel simulations,
[u/mrichter/AliRoot.git] / ITS / AliITS.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 /* $Id$ */
17
18
19 ///////////////////////////////////////////////////////////////////////////////
20 //
21 //      An overview of the basic philosophy of the ITS code development
22 // and analysis is show in the figure below.
23 //Begin_Html
24 /*
25 <img src="picts/ITS/ITS_Analysis_schema.gif">
26 </pre>
27 <br clear=left>
28 <font size=+2 color=red>
29 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
30 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
31 </font>
32 <pre>
33 */
34 //End_Html
35 //
36 //  AliITS. Inner Traking System base class.
37 //  This class contains the base procedures for the Inner Tracking System
38 //
39 //Begin_Html
40 /*
41 <img src="picts/ITS/AliITS_Class_Diagram.gif">
42 </pre>
43 <br clear=left>
44 <font size=+2 color=red>
45 <p>This show the class diagram of the different elements that are part of
46 the AliITS class.
47 </font>
48 <pre>
49 */
50 //End_Html
51 //
52 // Version: 0
53 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
54 //
55 // Version: 1
56 // Modified and documented by Bjorn S. Nilsen
57 // July 11 1999
58 //
59 // Version: 2
60 // Modified and documented by A. Bologna
61 // October 18 1999
62 //
63 // AliITS is the general base class for the ITS. Also see AliDetector for
64 // futher information.
65 //
66 ///////////////////////////////////////////////////////////////////////////////
67
68 #include <Riostream.h>
69 #include <stdlib.h>
70
71 #include <TBranch.h>
72 #include <TClonesArray.h>
73 #include <TFile.h>
74 #include <TMath.h>
75 #include <TObjectTable.h>
76 #include <TROOT.h>
77 #include <TRandom.h>
78 #include <TString.h>
79 #include <TTree.h>
80 #include <TVector.h>
81 #include <TVirtualMC.h>
82
83 #include "AliConfig.h"
84 #include "AliHeader.h"
85 #include "AliITS.h"
86 #include "AliITSClusterFinderSDD.h"
87 #include "AliITSClusterFinderSPD.h"
88 #include "AliITSClusterFinderSSD.h"
89 #include "AliITSDetType.h"
90 #include "AliITSLoader.h"
91 #include "AliITSRawClusterSPD.h"
92 #include "AliITSRawClusterSDD.h"
93 #include "AliITSRawClusterSSD.h"
94 #include "AliITSRecPoint.h"
95 #include "AliITSdigitSPD.h"
96 #include "AliITSdigitSDD.h"
97 #include "AliITSdigitSSD.h"
98 #include "AliITSgeom.h"
99 #include "AliITShit.h"
100 #include "AliITSmodule.h"
101 #include "AliITSpList.h"
102 #include "AliITSresponseSDD.h"
103 #include "AliITSresponseSPD.h"
104 #include "AliITSresponseSSD.h"
105 #include "AliITSsegmentationSDD.h"
106 #include "AliITSsegmentationSPD.h"
107 #include "AliITSsegmentationSSD.h"
108 #include "AliITSsimulationSDD.h"
109 #include "AliITSsimulationSPD.h"
110 #include "AliITSsimulationSSD.h"
111 #include "AliMC.h"
112 #include "AliITSDigitizer.h"
113 #include "AliITSDDLRawData.h"
114 #include "AliRun.h"
115
116 ClassImp(AliITS)
117
118 //______________________________________________________________________
119 AliITS::AliITS() : AliDetector(),
120     fITSgeom(0),
121     fEuclidOut(0),
122     fITSmodules(0),
123     fOpt("All"),
124     fIdN(0),
125     fIdSens(0),
126     fIdName(0),
127     fNDetTypes(kNTYPES),
128     fDetTypes(0),
129     fSDigits(0),
130     fNSDigits(0),
131     fDtype(0),
132     fNdtype(0),
133     fCtype(0),
134     fNctype(0),
135     fRecPoints(0),
136     fNRecPoints(0),
137     fSelectedVertexer(0){
138     // Default initializer for ITS
139     //      The default constructor of the AliITS class. In addition to
140     // creating the AliITS class it zeros the variables fIshunt (a member
141     // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
142     // fITSpoints, fIdSens, and fIdName. The AliDetector default constructor
143     // is also called.
144     // Inputs:
145     //      none.
146     // Outputs:
147     //      none.
148     // Return:
149     //      Blank ITS class.
150
151     fIshunt     = 0;   // not zeroed in AliDetector.
152
153     // AliITS variables.
154 //    SetDetectors(); // default to fOpt="All". This variable not written out.
155     SetMarkerColor(kRed);
156     SelectVertexer(" ");
157 }
158 //______________________________________________________________________
159 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title),
160     fITSgeom(0),
161     fEuclidOut(0),
162     fITSmodules(0),
163     fOpt("All"),
164     fIdN(0),
165     fIdSens(0),
166     fIdName(0),
167     fNDetTypes(kNTYPES),
168     fDetTypes(0),
169     fSDigits(0),
170     fNSDigits(0),
171     fDtype(0),
172     fNdtype(0),
173     fCtype(0),
174     fNctype(0),
175     fRecPoints(0),
176     fNRecPoints(0),
177     fSelectedVertexer(0){
178     //     The standard Constructor for the ITS class. In addition to 
179     // creating the AliITS class, it allocates memory for the TClonesArrays 
180     // fHits, fSDigits, fDigits, fITSpoints, and the TObjArray of fCtype 
181     // (clusters). It also zeros the variables
182     // fIshunt (a member of AliDetector class), fEuclidOut, and fIdN, and zeros
183     // the pointers fIdSens and fIdName. To help in displaying hits via the
184     // ROOT macro display.C AliITS also sets the marker color to red. The
185     // variables passes with this constructor, const char *name and *title,
186     // are used by the constructor of AliDetector class. See AliDetector
187     // class for a description of these parameters and its constructor
188     // functions.
189     // Inputs:
190     //      const char *name      Detector name. Should always be "ITS"
191     //      const char *title     Detector title.
192     // Outputs:
193     //      none.
194     // Return:
195     //      An ITS class.
196
197     fIshunt     = 0;  // not zeroed in AliDetector
198     fHits       = new TClonesArray("AliITShit", 1560);//not done in AliDetector
199     // Not done in AliDetector.
200     if(gAlice->GetMCApp()) gAlice->GetMCApp()->AddHitList(fHits);
201
202     fEuclidOut  = 0;
203     fITSgeom    = 0;
204     fITSmodules = 0;
205     SetDetectors(); // default to fOpt="All". This variable not written out.
206
207     fIdN        = 0;
208     fIdName     = 0;
209     fIdSens     = 0;
210
211     fNDetTypes  = kNTYPES;
212     fDetTypes   = new TObjArray(fNDetTypes);
213
214     fSDigits    = new TClonesArray("AliITSpListItem",1000);
215     fNSDigits   = 0;
216
217     fNdtype     = new Int_t[fNDetTypes];
218     fDtype      = new TObjArray(fNDetTypes);
219
220     fCtype      = new TObjArray(fNDetTypes);
221     fNctype     = new Int_t[fNDetTypes];
222
223     fRecPoints  = new TClonesArray("AliITSRecPoint",1000);
224     fNRecPoints = 0;
225
226     Int_t i;
227     for(i=0;i<fNDetTypes;i++) {
228         fDetTypes->AddAt(new AliITSDetType(),i); 
229         fNdtype[i] = 0;
230         fNctype[i] = 0;
231     } // end for i
232
233     SetMarkerColor(kRed);
234     SelectVertexer(" ");
235 }
236 //______________________________________________________________________
237 AliITS::~AliITS(){
238     // Default destructor for ITS.
239     //     The default destructor of the AliITS class. In addition to deleting
240     // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
241     // fSDigits, fCtype, fITSmodules, fITSgeom, fRecPoints, fIdSens, fIdName, 
242     // fITSpoints, fDetType and it's contents.
243     // Inputs:
244     //      none.
245     // Outputs:
246     //      none.
247     // Return:
248     //      none.
249
250     if (fHits) {
251       fHits->Delete();
252       delete fHits;
253       fHits=0;
254     }
255     if (fSDigits) {
256       fSDigits->Delete();
257       delete fSDigits;
258       fSDigits=0;
259     }
260     if (fDigits) {
261       fDigits->Delete();
262       delete fDigits;
263       fDigits=0;
264     }
265     if (fRecPoints) {
266       fRecPoints->Delete();
267       delete fRecPoints;
268       fRecPoints=0;
269     }
270     delete[] fIdName;  // Array of TStrings
271     delete[] fIdSens;
272     if(fITSmodules) {
273         this->ClearModules();
274         delete fITSmodules;
275         fITSmodules = 0;
276     }// end if fITSmodules!=0
277
278     if(fDtype) {
279         fDtype->Delete();
280         delete fDtype;
281     } // end if fDtype
282     delete [] fNdtype;
283     if (fCtype) {
284         fCtype->Delete();
285         delete fCtype;
286         fCtype = 0;
287     } // end if fCtype
288     delete [] fNctype;
289
290     if (fDetTypes) {
291         fDetTypes->Delete();
292         delete fDetTypes;
293         fDetTypes = 0;
294     } // end if fDetTypes
295
296
297     if (fITSgeom) delete fITSgeom;
298 }
299 //______________________________________________________________________
300 AliITS::AliITS(const AliITS &source) : AliDetector(source){
301     // Copy constructor. This is a function which is not allowed to be
302     // done to the ITS. It exits with an error.
303     // Inputs:
304     //      AliITS &source  An AliITS class.
305     // Outputs:
306     //      none.
307     // Return:
308     //      none.
309
310     if(this==&source) return;
311     Error("Copy constructor",
312           "You are not allowed to make a copy of the AliITS");
313     exit(1);
314 }
315 //______________________________________________________________________
316 AliITS& AliITS::operator=(AliITS &source){
317     // Assignment operator. This is a function which is not allowed to be
318     // done to the ITS. It exits with an error.
319     // Inputs:
320     //      AliITS &source  An AliITS class.
321     // Outputs:
322     //      none.
323     // Return:
324     //      none.
325
326     if(this==&source) return *this;
327     Error("operator=","You are not allowed to make a copy of the AliITS");
328     exit(1);
329     return *this; //fake return
330 }
331 //______________________________________________________________________
332 AliDigitizer* AliITS::CreateDigitizer(AliRunDigitizer* manager)const{
333     // Creates the AliITSDigitizer in a standard way for use via AliModule.
334     // This function can not be included in the .h file because of problems
335     // with the order of inclusion (recursive).
336     // Inputs:
337     //    AliRunDigitizer *manager  The Manger class for Digitization
338     // Output:
339     //    none.
340     // Return:
341     //    A new AliITSRunDigitizer (cast as a AliDigitizer).
342
343      return new AliITSDigitizer(manager);
344 }
345 //______________________________________________________________________
346 void AliITS::Init(){
347     // Initializer ITS after it has been built
348     //     This routine initializes the AliITS class. It is intended to be
349     // called from the Init function in AliITSv?. Besides displaying a banner
350     // indicating that it has been called it initializes the array fIdSens
351     // and sets the default segmentation, response, digit and raw cluster
352     // classes therefore it should be called after a call to CreateGeometry.
353     // Inputs:
354     //      none.
355     // Outputs:
356     //      none.
357     // Return:
358     //      none.
359     Int_t i;
360
361     SetDefaults();
362     // Array of TStrings
363     if(gMC) for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
364 }
365 //______________________________________________________________________
366 void AliITS::SetDefaults(){
367     // sets the default segmentation, response, digit and raw cluster classes.
368     // Inputs:
369     //      none.
370     // Outputs:
371     //      none.
372     // Return:
373     //      none.
374     AliITSsegmentation *seg;
375     AliITSresponse    *resp;
376     AliITSDetType *iDetType;
377
378     if(fDebug) Info("SetDefauls","%s: SetDefaults",ClassName());
379
380     //SPD
381     iDetType = DetType(kSPD);
382     if(iDetType){
383         if (!iDetType->GetSegmentationModel()) {
384             seg = new AliITSsegmentationSPD(fITSgeom);
385             SetSegmentationModel(kSPD,seg); 
386         } // end if
387         if (!iDetType->GetResponseModel()) {
388             SetResponseModel(kSPD,new AliITSresponseSPD()); 
389         } // end if
390         // set digit and raw cluster classes to be used
391
392         const char *kData0=(iDetType->GetResponseModel())->DataType();
393         if (strstr(kData0,"real")) {
394             iDetType->ClassNames("AliITSdigit","AliITSRawClusterSPD");
395         } else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
396     } // end if iDetType
397
398     // SDD
399     iDetType = DetType(kSDD); 
400     if(iDetType){
401         if (!iDetType->GetResponseModel()) {
402             SetResponseModel(kSDD,new AliITSresponseSDD("simulated")); 
403         } // end if
404         resp = iDetType->GetResponseModel();
405         if (!iDetType->GetSegmentationModel()) {
406             seg = new AliITSsegmentationSDD(fITSgeom,resp);
407             SetSegmentationModel(kSDD,seg); 
408         } // end if
409         const char *kData1 = (iDetType->GetResponseModel())->DataType();
410         const char *kopt = iDetType->GetResponseModel()->ZeroSuppOption();
411         if((!strstr(kopt,"2D"))&&(!strstr(kopt,"1D")) || 
412                                    strstr(kData1,"real") ){
413             iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
414         } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
415     } // end if iDetType
416
417     // SSD
418     iDetType = DetType(kSSD); 
419     if(iDetType){
420         if (!iDetType->GetSegmentationModel()) {
421             seg = new AliITSsegmentationSSD(fITSgeom);
422             SetSegmentationModel(kSSD,seg); 
423         } // end if
424         if (!iDetType->GetResponseModel()) {
425             SetResponseModel(kSSD,new AliITSresponseSSD("simulated"));
426         } // end if
427         const char *kData2 = (iDetType->GetResponseModel())->DataType();
428         if (strstr(kData2,"real")) {
429             iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
430         } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
431     } // end if iDetType
432
433     if (kNTYPES>3) {
434         Warning("SetDefaults",
435                 "Only the three basic detector types are initialized!");
436     }  // end if
437 }
438 //______________________________________________________________________
439 void AliITS::SetDefaultSimulation(){
440     // sets the default simulation.
441     // Inputs:
442     //      none.
443     // Outputs:
444     //      none.
445     // Return:
446     //      none.
447     AliITSDetType *iDetType;
448     AliITSsimulation *sim;
449     AliITSsegmentation *seg;
450     AliITSresponse *res;
451
452     iDetType = DetType(kSPD);
453     if(iDetType){
454         sim = iDetType->GetSimulationModel();
455         if (!sim) {
456             seg = (AliITSsegmentation*)iDetType->GetSegmentationModel();
457             res = (AliITSresponse*)iDetType->GetResponseModel();
458             sim = new AliITSsimulationSPD(seg,res);
459             SetSimulationModel(kSPD,sim);
460         }else{ // simulation exists, make sure it is set up properly.
461             if(sim->GetResponseModel()==0) sim->SetResponseModel(
462                 (AliITSresponse*)iDetType->GetResponseModel());
463             if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
464                 (AliITSsegmentation*)iDetType->GetSegmentationModel());
465             sim->Init();
466         } // end if
467     }// end if iDetType
468     iDetType  = DetType(kSDD);
469     if(iDetType){
470         sim = iDetType->GetSimulationModel();
471         if (!sim) {
472             seg = (AliITSsegmentation*)iDetType->GetSegmentationModel();
473             res =  (AliITSresponse*)iDetType->GetResponseModel();
474             sim = new AliITSsimulationSDD(seg,res);
475             SetSimulationModel(kSDD,sim);
476         }else{ // simulation exists, make sure it is set up properly.
477             if(sim->GetResponseModel()==0) sim->SetResponseModel(
478                 (AliITSresponse*)iDetType->GetResponseModel());
479             if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
480                 (AliITSsegmentation*)iDetType->GetSegmentationModel());
481             sim->Init();
482         } //end if
483     }// end if iDetType
484     iDetType = DetType(kSSD);
485     if(iDetType){
486         sim = iDetType->GetSimulationModel();
487         if (!sim) {
488             seg =(AliITSsegmentation*)iDetType->GetSegmentationModel();
489             res =  (AliITSresponse*)iDetType->GetResponseModel();
490             sim = new AliITSsimulationSSD(seg,res);
491             SetSimulationModel(kSSD,sim);
492         }else{ // simulation exists, make sure it is set up properly.
493             if(sim->GetResponseModel()==0) sim->SetResponseModel(
494                 (AliITSresponse*)iDetType->GetResponseModel());
495             if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
496                 (AliITSsegmentation*)iDetType->GetSegmentationModel());
497             sim->Init();
498         } // end if
499     } // end if iDetType
500 }
501 //______________________________________________________________________
502 void AliITS::SetDefaultClusterFinders(){
503     // Sets the default cluster finders. Used in finding RecPoints.
504     // Inputs:
505     //      none.
506     // Outputs:
507     //      none.
508     // Return:
509     //      none.
510     AliITSDetType *iDetType;
511     AliITSsegmentation *seg;
512     AliITSClusterFinder *clf;
513     AliITSresponse *res;
514
515     MakeTreeC();
516
517     // SPD
518     iDetType=DetType(kSPD);
519     if(iDetType){
520         if (!iDetType->GetReconstructionModel()) {
521             seg =(AliITSsegmentation*)iDetType->GetSegmentationModel();
522             TClonesArray  *dig0 = DigitsAddress(0);
523             TClonesArray *recp0 = ClustersAddress(0);
524             clf = new AliITSClusterFinderSPD(seg,dig0,recp0);
525             SetReconstructionModel(kSPD,clf);
526         } // end if
527     } // end if iDetType
528
529     // SDD
530     iDetType=DetType(kSDD);
531     if(iDetType){
532         if (!iDetType->GetReconstructionModel()) {
533             seg = (AliITSsegmentation*)iDetType->GetSegmentationModel();
534             res = (AliITSresponse*)iDetType->GetResponseModel();
535             TClonesArray  *dig1 = DigitsAddress(1);
536             TClonesArray *recp1 = ClustersAddress(1);
537             clf = new AliITSClusterFinderSDD(seg,res,dig1,recp1);
538             SetReconstructionModel(kSDD,clf);
539         } // end if
540     } // end if iDetType
541
542     // SSD
543     iDetType=DetType(kSSD);
544     if(iDetType){
545         if (!iDetType->GetReconstructionModel()) {
546             seg = (AliITSsegmentation*)iDetType->GetSegmentationModel();
547             TClonesArray *dig2 = DigitsAddress(2);
548             clf = new AliITSClusterFinderSSD(seg,dig2);
549             SetReconstructionModel(kSSD,clf);
550         } // end if
551     } // end if iDetType
552 }
553 //______________________________________________________________________
554 void AliITS::MakeBranch(Option_t* option){
555     // Creates Tree branches for the ITS.
556     // Inputs:
557     //      Option_t *option    String of Tree types S,D, and/or R.
558     //      const char *file    String of the file name where these branches
559     //                          are to be stored. If blank then these branches
560     //                          are written to the same tree as the Hits were
561     //                          read from.
562     // Outputs:
563     //      none.
564     // Return:
565     //      none.
566     Bool_t cH = (strstr(option,"H")!=0);
567     Bool_t cS = (strstr(option,"S")!=0);
568     Bool_t cD = (strstr(option,"D")!=0);
569     Bool_t cR = (strstr(option,"R")!=0);
570     Bool_t cRF = (strstr(option,"RF")!=0);
571     
572     if(cRF)cR = kFALSE;
573     if(cH && (fHits == 0x0)) fHits  = new TClonesArray("AliITShit", 1560);
574     
575     AliDetector::MakeBranch(option);
576
577     if(cS) MakeBranchS(0);
578     if(cD) MakeBranchD(0);
579     if(cR) MakeBranchR(0);
580     if(cRF) MakeBranchRF(0);
581 }
582 //______________________________________________________________________
583 void AliITS::SetTreeAddress(){
584     // Set branch address for the Trees.
585     // Inputs:
586     //      none.
587     // Outputs:
588     //      none.
589     // Return:
590     //      none.
591     TTree *treeS = fLoader->TreeS();
592     TTree *treeD = fLoader->TreeD();
593     TTree *treeR = fLoader->TreeR();
594     if (fLoader->TreeH() && (fHits == 0x0)) 
595         fHits = new TClonesArray("AliITShit", 1560);
596
597     AliDetector::SetTreeAddress();
598
599     SetTreeAddressS(treeS);
600     SetTreeAddressD(treeD);
601     SetTreeAddressR(treeR);
602 }
603 //______________________________________________________________________
604 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
605     // Add an ITS hit
606     //     The function to add information to the AliITShit class. See the
607     // AliITShit class for a full description. This function allocates the
608     // necessary new space for the hit information and passes the variable
609     // track, and the pointers *vol and *hits to the AliITShit constructor
610     // function.
611     // Inputs:
612     //      Int_t   track   Track number which produced this hit.
613     //      Int_t   *vol    Array of Integer Hit information. See AliITShit.h
614     //      Float_t *hits   Array of Floating Hit information.  see AliITShit.h
615     // Outputs:
616     //      none.
617     // Return:
618     //      none.
619
620     TClonesArray &lhits = *fHits;
621     new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
622 }
623 //______________________________________________________________________
624 void AliITS::InitModules(Int_t size,Int_t &nmodules){
625     // Initialize the modules array.
626     // Inputs:
627     //      Int_t size  Size of array of the number of modules to be
628     //                  created. If size <=0 then the number of modules
629     //                  is gotten from AliITSgeom class kept in fITSgeom.
630     // Outputs:
631     //      Int_t &nmodules The number of modules existing.
632     // Return:
633     //      none.
634
635     if(fITSmodules){ 
636         fITSmodules->Delete();
637         delete fITSmodules;
638     } // end fir fITSmoudles
639
640     Int_t nl,indexMAX,index;
641
642     if(size<=0){ // default to using data stored in AliITSgeom
643         if(fITSgeom==0) {
644             Error("InitModules","fITSgeom not defined");
645             return;
646         } // end if fITSgeom==0
647         nl = fITSgeom->GetNlayers();
648         indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
649                                             fITSgeom->GetNdetectors(nl))+1;
650         nmodules = indexMAX;
651         fITSmodules = new TObjArray(indexMAX);
652         for(index=0;index<indexMAX;index++){
653             fITSmodules->AddAt( new AliITSmodule(index),index);
654         } // end for index
655     }else{
656         fITSmodules = new TObjArray(size);
657         for(index=0;index<size;index++) {
658             fITSmodules->AddAt( new AliITSmodule(index),index);
659         } // end for index
660
661         nmodules = size;
662     } // end i size<=0
663 }
664 //______________________________________________________________________
665 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,
666                          Option_t *option, const char *filename){
667     // fill the modules with the sorted by module hits; add hits from
668     // background if option=Add.
669     // Inputs:
670     //      Int_t evnt       Event to be processed.
671     //      Int_t bgrev      Background Hit tree number.
672     //      Int_t nmodules   Not used.
673     //      Option_t *option String indicating if merging hits or not. To
674     //                       merge hits set equal to "Add". Otherwise no
675     //                       background hits are considered.
676     //      Test_t *filename File name containing the background hits..
677     // Outputs:
678     //      none.
679     // Return:
680     //      none.
681     static TTree *trH1;                 //Tree with background hits
682     static Bool_t first=kTRUE;
683     static TFile *file;
684     const char *addBgr = strstr(option,"Add");
685
686     evnt = nmodules; // Dummy use of variables to remove warnings
687     if (addBgr ) {
688         if(first) {
689             file=new TFile(filename);
690         } // end if first
691         first=kFALSE;
692         file->cd();
693         file->ls();
694         // Get Hits Tree header from file
695         if(trH1) delete trH1;
696         trH1=0;
697
698         char treeName[20];
699         sprintf(treeName,"TreeH%d",bgrev);
700         trH1 = (TTree*)gDirectory->Get(treeName);
701         if (!trH1) {
702             Error("FillModules","cannot find Hits Tree for event:%d",bgrev);
703         } // end if !trH1
704         // Set branch addresses
705     } // end if addBgr
706
707     FillModules(fLoader->TreeH(),0); // fill from this file's tree.
708     
709     if (addBgr ) {
710         FillModules(trH1,10000000); // Default mask 10M.
711         TTree *fAli=fLoader->GetRunLoader()->TreeK();
712         TFile *fileAli=0;
713         if (fAli) fileAli =fAli->GetCurrentFile();
714         fileAli->cd();
715     } // end if add
716 }
717 //______________________________________________________________________
718 void AliITS::FillModules(TTree *treeH, Int_t mask) {
719     // fill the modules with the sorted by module hits; 
720     // can be called many times to do a merging
721     // Inputs:
722     //      TTree *treeH  The tree containing the hits to be copied into
723     //                    the modules.
724     //      Int_t mask    The track number mask to indecate which file
725     //                    this hits came from.
726     // Outputs:
727     //      none.
728     // Return:
729     //      none.
730
731     if (treeH == 0x0)
732      {
733        Error("FillModules","Tree is NULL");
734      }
735     Int_t lay,lad,det,index;
736     AliITShit *itsHit=0;
737     AliITSmodule *mod=0;
738     char branchname[20];
739     sprintf(branchname,"%s",GetName());
740     TBranch *branch = treeH->GetBranch(branchname);
741     if (!branch) {
742         Error("FillModules","%s branch in TreeH not found",branchname);
743         return;
744     } // end if !branch
745     branch->SetAddress(&fHits);
746     Int_t nTracks =(Int_t) treeH->GetEntries();
747     Int_t iPrimTrack,h;
748     for(iPrimTrack=0; iPrimTrack<nTracks; iPrimTrack++){
749         ResetHits();
750         Int_t nBytes = treeH->GetEvent(iPrimTrack);
751         if (nBytes <= 0) continue;
752         Int_t nHits = fHits->GetEntriesFast();
753         for(h=0; h<nHits; h++){
754             itsHit = (AliITShit *)fHits->UncheckedAt(h);
755             itsHit->GetDetectorID(lay,lad,det);
756             if (fITSgeom) {
757                 index = fITSgeom->GetModuleIndex(lay,lad,det);
758             } else {
759                 index=det-1; // This should not be used.
760             } // end if [You must have fITSgeom for this to work!]
761             mod = GetModule(index);
762             itsHit->SetTrack(itsHit->GetTrack()+mask); // Set track mask.
763             mod->AddHit(itsHit,iPrimTrack,h);
764         } // end loop over hits 
765     } // end loop over tracks
766 }
767 //______________________________________________________________________
768 void AliITS::MakeBranchS(const char *fl){
769     // Creates Tree Branch for the ITS summable digits.
770     // Inputs:
771     //      cont char *fl  File name where SDigits branch is to be written
772     //                     to. If blank it write the SDigits to the same
773     //                     file in which the Hits were found.
774     // Outputs:
775     //      none.
776     // Return:
777     //      none.
778     Int_t buffersize = 4000;
779     char branchname[30];
780
781     // only one branch for SDigits.
782     sprintf(branchname,"%s",GetName());
783     
784
785     if(fLoader->TreeS()){
786         if(fSDigits==0x0) fSDigits = new TClonesArray("AliITSpListItem",1000);
787         MakeBranchInTree(fLoader->TreeS(),branchname,&fSDigits,buffersize,fl);
788     } // end if
789 }
790 //______________________________________________________________________
791 void AliITS::SetTreeAddressS(TTree *treeS){
792     // Set branch address for the ITS summable digits Trees.
793     // Inputs:
794     //      TTree *treeS   Tree containing the SDigits.
795     // Outputs:
796     //      none.
797     // Return:
798     //      none.
799     char branchname[30];
800
801     if(!treeS) return;
802     if (fSDigits == 0x0)  fSDigits = new TClonesArray("AliITSpListItem",1000);
803     TBranch *branch;
804     sprintf(branchname,"%s",GetName());
805     branch = treeS->GetBranch(branchname);
806     if (branch) branch->SetAddress(&fSDigits);
807 }
808 //______________________________________________________________________
809 void AliITS::MakeBranchInTreeD(TTree *treeD,const char *file){
810     // Creates Tree branches for the ITS.
811     // Inputs:
812     //      TTree     *treeD Pointer to the Digits Tree.
813     //      cont char *file  File name where Digits branch is to be written
814     //                       to. If blank it write the SDigits to the same
815     //                       file in which the Hits were found.
816     // Outputs:
817     //      none.
818     // Return:
819     //      none.
820     // one branch for digits per type of detector
821     const Char_t *det[3] = {"SPD","SDD","SSD"};
822     TString digclass;
823     Int_t i;
824     Int_t buffersize = 4000;
825     Char_t branchname[30];
826
827     if (fDtype == 0x0) fDtype = new TObjArray(fNDetTypes);
828     for (i=0; i<kNTYPES ;i++) {
829         digclass = DetType(i)->GetDigitClassName();
830         // digits
831         if(!(fDtype->At(i))){
832             fDtype->AddAt(new TClonesArray(digclass.Data(),1000),i);
833         }else{
834             ResetDigits(i);
835         } // end if
836         if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
837         else  sprintf(branchname,"%sDigits%d",GetName(),i+1);      
838         if (fDtype && treeD) {
839             MakeBranchInTree(treeD,branchname,&((*fDtype)[i]),buffersize,file);
840         } // end if
841     } // end for i
842          /*
843     for (i=0; i<kNTYPES ;i++) {
844         if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
845         else  sprintf(branchname,"%sDigits%d",GetName(),i+1);      
846         if (fDtype && treeD) {
847             MakeBranchInTree(treeD,branchname,&((*fDtype)[i]),buffersize,file);
848         } // end if
849     } // end for i
850          */
851 }
852 //______________________________________________________________________
853 void AliITS::SetTreeAddressD(TTree *treeD){
854     // Set branch address for the Trees.
855     // Inputs:
856     //      TTree *treeD   Tree containing the Digits.
857     // Outputs:
858     //      none.
859     // Return:
860     //      none.
861     const char *det[3] = {"SPD","SDD","SSD"};
862     TBranch *branch;
863     TString digclass;
864     Int_t i;
865     char branchname[30];
866
867     if(!treeD) return;
868     if (fDtype == 0x0) fDtype = new TObjArray(fNDetTypes);
869     for (i=0; i<kNTYPES; i++) {
870         digclass = DetType(i)->GetDigitClassName();
871         // digits
872         if(!(fDtype->At(i))) {
873             fDtype->AddAt(new TClonesArray(digclass.Data(),1000),i);
874         }else{
875             ResetDigits(i);
876         } // end if
877         if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
878         else  sprintf(branchname,"%sDigits%d",GetName(),i+1);
879         if (fDtype) {
880             branch = treeD->GetBranch(branchname);
881             if (branch) branch->SetAddress(&((*fDtype)[i]));
882         } // end if fDtype
883     } // end for i
884 }
885 //______________________________________________________________________
886 void AliITS::Hits2SDigits(){
887     // Standard Hits to summable Digits function.
888     // Inputs:
889     //      none.
890     // Outputs:
891     //      none.
892
893     fLoader->LoadHits("read");
894     fLoader->LoadSDigits("recreate");
895     AliRunLoader* rl = fLoader->GetRunLoader(); 
896
897     for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
898         // Do the Hits to Digits operation. Use Standard input values.
899         // Event number from file, no background hit merging , use size from
900         // AliITSgeom class, option="All", input from this file only.
901         rl->GetEvent(iEvent);
902         if (!fLoader->TreeS()) fLoader->MakeTree("S");
903         MakeBranch("S");
904         SetTreeAddress();
905         HitsToPreDigits(iEvent,0,-1," ",fOpt," ");
906     } // end for iEvent
907     
908     fLoader->UnloadHits();
909 }
910 //______________________________________________________________________
911 void AliITS::SDigitsToDigits(Option_t *opt){
912     // Standard Summable digits to Digits function.
913     // Inputs:
914     //      none.
915     // Outputs:
916     //      none.
917     char name[20] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
918
919     if(!GetITSgeom()) return; // need transformations to do digitization.
920     AliITSgeom *geom = GetITSgeom();
921
922     const char *all = strstr(opt,"All");
923     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
924                           strstr(opt,"SSD")};
925     if( !det[0] && !det[1] && !det[2] ) all = "All";
926     else all = 0;
927     static Bool_t setDef = kTRUE;
928     if(setDef) SetDefaultSimulation();
929     setDef = kFALSE;
930
931     AliITSsimulation *sim      = 0;
932     AliITSDetType    *iDetType = 0;
933     TTree            *trees    = fLoader->TreeS();
934     if( !(trees && this->GetSDigits()) ){
935         Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
936         return;
937     } // end if
938     sprintf( name, "%s", this->GetName() );
939     TBranch *brchSDigits = trees->GetBranch( name );
940     
941     Int_t id,module;
942     for(module=0;module<geom->GetIndexMax();module++){
943         id       = geom->GetModuleType(module);
944         if (!all && !det[id]) continue;
945         iDetType = DetType(id);
946         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
947         if (!sim) {
948             Error("SDigit2Digits","The simulation class was not "
949                   "instanciated for module %d type %s!",module,
950                   geom->GetModuleTypeName(module));
951             exit(1);
952         } // end if !sim
953         sim->InitSimulationModule(module,gAlice->GetEvNumber());
954         //
955         // add summable digits to module
956         this->GetSDigits()->Clear();
957         brchSDigits->GetEvent(module);
958         sim->AddSDigitsToModule(GetSDigits(),0);
959         //
960         // Digitise current module sum(SDigits)->Digits
961         sim->FinishSDigitiseModule();
962
963         // fills all branches - wasted disk space
964         fLoader->TreeD()->Fill();
965         this->ResetDigits();
966     } // end for module
967
968     fLoader->TreeD()->GetEntries();
969
970     fLoader->TreeD()->AutoSave();
971     // reset tree
972     fLoader->TreeD()->Reset();
973     
974 }
975 //______________________________________________________________________
976 void AliITS::Hits2Digits(){
977     // Standard Hits to Digits function.
978     // Inputs:
979     //      none.
980     // Outputs:
981     //      none.
982
983     fLoader->LoadHits("read");
984     fLoader->LoadDigits("recreate");
985     AliRunLoader* rl = fLoader->GetRunLoader(); 
986
987     for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
988         // Do the Hits to Digits operation. Use Standard input values.
989         // Event number from file, no background hit merging , use size from
990         // AliITSgeom class, option="All", input from this file only.
991         rl->GetEvent(iEvent);
992         if (!fLoader->TreeD()) fLoader->MakeTree("D");
993         MakeBranch("D");
994         SetTreeAddress();   
995         HitsToDigits(iEvent,0,-1," ",fOpt," ");
996     } // end for iEvent
997
998     fLoader->UnloadHits();
999     fLoader->UnloadDigits();
1000 }
1001 //______________________________________________________________________
1002 void AliITS::HitsToPreDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1003                              Option_t *option,Option_t *opt,
1004                              const char *filename){
1005     //   Keep galice.root for signal and name differently the file for 
1006     // background when add! otherwise the track info for signal will be lost !
1007     // the condition below will disappear when the geom class will be
1008     // initialized for all versions - for the moment it is only for v5 !
1009     // 7 is the SDD beam test version.
1010     // Inputs:
1011     //      Int_t evnt       Event to be processed.
1012     //      Int_t bgrev      Background Hit tree number.
1013     //      Int_t nmodules   Not used.
1014     //      Option_t *option String indicating if merging hits or not. To
1015     //                       merge hits set equal to "Add". Otherwise no
1016     //                       background hits are considered.
1017     //      Test_t *filename File name containing the background hits..
1018     // Outputs:
1019     //      none.
1020     // Return:
1021     //      none.
1022
1023     if(!GetITSgeom()) return; // need transformations to do digitization.
1024     AliITSgeom *geom = GetITSgeom();
1025
1026     const char *all = strstr(opt,"All");
1027     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1028                           strstr(opt,"SSD")};
1029     static Bool_t setDef=kTRUE;
1030     if (setDef) SetDefaultSimulation();
1031     setDef=kFALSE;
1032
1033     Int_t nmodules;
1034     InitModules(size,nmodules);
1035     FillModules(evNumber,bgrev,nmodules,option,filename);
1036
1037     AliITSsimulation *sim      = 0;
1038     AliITSDetType    *iDetType = 0;
1039     AliITSmodule     *mod      = 0;
1040     Int_t id,module;
1041     for(module=0;module<geom->GetIndexMax();module++){
1042         id       = geom->GetModuleType(module);
1043         if (!all && !det[id]) continue;
1044         iDetType = DetType(id);
1045         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1046         if (!sim) {
1047             Error("HitsToSDigits","The simulation class was not "
1048                   "instanciated for module %d type %s!",module,
1049                   geom->GetModuleTypeName(module));
1050             exit(1);
1051         } // end if !sim
1052         mod      = (AliITSmodule *)fITSmodules->At(module);
1053         sim->SDigitiseModule(mod,module,evNumber);
1054         // fills all branches - wasted disk space
1055         fLoader->TreeS()->Fill(); 
1056         ResetSDigits();
1057     } // end for module
1058
1059     ClearModules();
1060
1061     fLoader->TreeS()->GetEntries();
1062     fLoader->TreeS()->AutoSave();
1063     fLoader->WriteSDigits("OVERWRITE");
1064     // reset tree
1065     fLoader->TreeS()->Reset();
1066 }
1067 //______________________________________________________________________
1068 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1069                           Option_t *option,Option_t *opt,
1070                           const char *filename){
1071     //   Keep galice.root for signal and name differently the file for 
1072     // background when add! otherwise the track info for signal will be lost !
1073     // the condition below will disappear when the geom class will be
1074     // initialized for all versions - for the moment it is only for v5 !
1075     // 7 is the SDD beam test version.
1076     // Inputs:
1077     //      Int_t evnt       Event to be processed.
1078     //      Int_t bgrev      Background Hit tree number.
1079     //      Int_t nmodules   Not used.
1080     //      Option_t *option String indicating if merging hits or not. To
1081     //                       merge hits set equal to "Add". Otherwise no
1082     //                       background hits are considered.
1083     //      Test_t *filename File name containing the background hits..
1084     // Outputs:
1085     //      none.
1086     // Return:
1087     //      none.
1088
1089     if(!GetITSgeom()) return; // need transformations to do digitization.
1090     AliITSgeom *geom = GetITSgeom();
1091
1092     const char *all = strstr(opt,"All");
1093     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1094                           strstr(opt,"SSD")};
1095     static Bool_t setDef=kTRUE;
1096     if (setDef) SetDefaultSimulation();
1097     setDef=kFALSE;
1098
1099     Int_t nmodules;
1100     InitModules(size,nmodules);
1101     FillModules(evNumber,bgrev,nmodules,option,filename);
1102
1103     AliITSsimulation *sim      = 0;
1104     AliITSDetType    *iDetType = 0;
1105     AliITSmodule     *mod      = 0;
1106     Int_t id,module;
1107     for(module=0;module<geom->GetIndexMax();module++){
1108         id       = geom->GetModuleType(module);
1109         if (!all && !det[id]) continue;
1110         iDetType = DetType(id);
1111         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1112         if (!sim) {
1113             Error("HitsToDigits","The simulation class was not "
1114                   "instanciated for module %d type %s!",module,
1115                   geom->GetModuleTypeName(module));
1116             exit(1);
1117         } // end if !sim
1118         mod      = (AliITSmodule *)fITSmodules->At(module);
1119         sim->DigitiseModule(mod,module,evNumber);
1120         // fills all branches - wasted disk space
1121         fLoader->TreeD()->Fill(); 
1122         ResetDigits();
1123     } // end for module
1124
1125     ClearModules();
1126
1127     fLoader->TreeD()->GetEntries();
1128     fLoader->TreeD()->AutoSave();
1129     // reset tree
1130     fLoader->TreeD()->Reset();
1131 }
1132 //______________________________________________________________________
1133 void AliITS::ResetDigits(){
1134     // Reset number of digits and the digits array for the ITS detector.
1135     // Inputs:
1136     //      none.
1137     // Outputs:
1138     //      none.
1139
1140     if (!fDtype) return;
1141
1142     Int_t i;
1143     for (i=0;i<kNTYPES;i++ ) {
1144         ResetDigits(i);
1145     } // end for i
1146 }
1147 //______________________________________________________________________
1148 void AliITS::ResetDigits(Int_t i){
1149     // Reset number of digits and the digits array for this branch.
1150     // Inputs:
1151     //      none.
1152     // Outputs:
1153     //      none.
1154
1155     if (fDtype->At(i)) ((TClonesArray*)fDtype->At(i))->Clear();
1156     if (fNdtype)       fNdtype[i]=0;
1157 }
1158 //______________________________________________________________________
1159 void AliITS::AddSumDigit(AliITSpListItem &sdig){
1160     // Adds the a module full of summable digits to the summable digits tree.
1161     // Inputs:
1162     //      AliITSpListItem &sdig   SDigit to be added to SDigits tree.
1163     // Outputs:
1164     //      none.
1165     // Return:
1166     //      none.
1167
1168     TClonesArray &lsdig = *fSDigits;
1169     new(lsdig[fNSDigits++]) AliITSpListItem(sdig);
1170 }
1171 //______________________________________________________________________
1172 void AliITS::AddRealDigit(Int_t id, Int_t *digits){
1173     //   Add a real digit - as coming from data.
1174     // Inputs:
1175     //      Int_t id        Detector type number.
1176     //      Int_t *digits   Integer array containing the digits info. See 
1177     //                      AliITSdigit.h
1178     // Outputs:
1179     //      none.
1180     // Return:
1181     //      none.
1182
1183     TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1184     new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
1185 }
1186 //______________________________________________________________________
1187 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d){
1188     //    Add a simulated digit.
1189     // Inputs:
1190     //      Int_t id        Detector type number.
1191     //      AliITSdigit *d  Digit to be added to the Digits Tree. See 
1192     //                      AliITSdigit.h
1193     // Outputs:
1194     //      none.
1195     // Return:
1196     //      none.
1197
1198     TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1199
1200     switch(id){
1201     case 0:
1202         new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
1203         break;
1204     case 1:
1205         new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
1206         break;
1207     case 2:
1208         new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
1209         break;
1210     } // end switch id
1211 }
1212 //______________________________________________________________________
1213 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,
1214                          Int_t *hits,Float_t *charges){
1215   //   Add a simulated digit to the list.
1216   // Inputs:
1217   //      Int_t id        Detector type number.
1218   //      Float_t phys    Physics indicator. See AliITSdigits.h
1219   //      Int_t *digits   Integer array containing the digits info. See 
1220   //                      AliITSdigit.h
1221   //      Int_t *tracks   Integer array [AliITSdigitS?D::GetNTracks()] 
1222   //                      containing the track numbers that contributed to
1223   //                      this digit.
1224   //      Int_t *hits     Integer array [AliITSdigitS?D::GetNTracks()]
1225   //                      containing the hit numbers, from AliITSmodule, that
1226   //                      contributed to this digit.
1227   //      Float_t *charge Floating point array of the signals contributed
1228   //                      to this digit by each track.
1229   // Outputs:
1230   //      none.
1231   // Return:
1232   //      none.
1233
1234   TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1235   AliITSresponseSDD *resp = 0;
1236   switch(id){
1237   case 0:
1238       new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
1239       break;
1240   case 1:
1241       resp = (AliITSresponseSDD*)DetType(1)->GetResponseModel();
1242       new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,
1243                                                  hits,charges,resp);
1244       break;
1245   case 2:
1246       new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
1247       break;
1248   } // end switch id
1249 }
1250 //______________________________________________________________________
1251 void AliITS::Digits2Raw(){
1252     // convert digits of the current event to raw data
1253
1254   fLoader->LoadDigits();
1255   TTree* digits = fLoader->TreeD();
1256   if (!digits) {
1257       Error("Digits2Raw", "no digits tree");
1258       return;
1259   }
1260   SetTreeAddressD(digits);
1261
1262   AliITSDDLRawData rawWriter;
1263   //Verbose level
1264   // 0: Silent
1265   // 1: cout messages
1266   // 2: txt files with digits 
1267   //BE CAREFUL, verbose level 2 MUST be used only for debugging and
1268   //it is highly suggested to use this mode only for debugging digits files
1269   //reasonably small, because otherwise the size of the txt files can reach
1270   //quickly several MB wasting time and disk space.
1271   rawWriter.SetVerbose(0);
1272     
1273   //SILICON PIXEL DETECTOR
1274   Info("Digits2Raw", "Formatting raw data for SPD");
1275   rawWriter.RawDataSPD(digits->GetBranch("ITSDigitsSPD"));
1276     
1277   //SILICON DRIFT DETECTOR
1278   Info("Digits2Raw", "Formatting raw data for SDD");
1279   rawWriter.RawDataSDD(digits->GetBranch("ITSDigitsSDD"));
1280     
1281   //SILICON STRIP DETECTOR
1282   Info("Digits2Raw", "Formatting raw data for SSD");
1283   rawWriter.RawDataSSD(digits->GetBranch("ITSDigitsSSD"));
1284
1285   fLoader->UnloadDigits();
1286 }
1287
1288 //______________________________________________________________________
1289 void AliITS::MakeTreeC(Option_t *option){
1290   //   Create a separate tree to store the clusters.
1291   // Inputs:
1292   //      Option_t *option  string which must contain "C" otherwise
1293   //                        no Cluster Tree is created.
1294   // Outputs:
1295   //      none.
1296   // Return:
1297   //      none.
1298
1299   AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;    
1300     
1301   if (pITSLoader == 0x0) {
1302       Error("MakeTreeC","fLoader == 0x0 option=%s",option);
1303       return;
1304   }
1305   if (pITSLoader->TreeC() == 0x0) pITSLoader->MakeTree("C");
1306   MakeBranchC();
1307 }
1308 //----------------------------------------------------------------------
1309 void AliITS::MakeBranchC(){
1310     //Makes barnches in treeC
1311     AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;    
1312     if (pITSLoader == 0x0) {
1313         Error("MakeTreeC","fLoader == 0x0");
1314         return;
1315     }
1316     TTree * lTC = pITSLoader->TreeC();
1317     if (lTC == 0x0){
1318         Error("MakeTreeC","Can not get TreeC from Loader");
1319         return;
1320     }
1321
1322     Int_t buffersize = 4000;
1323     char branchname[30];
1324     const char *det[3] = {"SPD","SDD","SSD"};
1325     char digclass[40];
1326     char clclass[40];
1327     
1328     // one branch for Clusters per type of detector
1329     Int_t i;   
1330     for (i=0; i<kNTYPES ;i++) {
1331         AliITSDetType *iDetType=DetType(i); 
1332         iDetType->GetClassNames(digclass,clclass);
1333         // clusters
1334         if (fCtype == 0x0) fCtype  = new TObjArray(fNDetTypes);
1335         if(!ClustersAddress(i)){
1336             fCtype->AddAt(new TClonesArray(clclass,1000),i);
1337         }
1338         if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
1339         else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
1340         if (fCtype  && lTC) {
1341             if (lTC->GetBranch(branchname)){
1342                 Warning("MakeBranchC","Branch %s alread exists in TreeC",
1343                         branchname);
1344             }else{
1345                 Info("MakeBranchC","Creating branch %s in TreeC",branchname);
1346                 lTC->Branch(branchname,&((*fCtype)[i]), buffersize);
1347             }
1348         } // end if fCtype && lTC
1349     } // end for i
1350 }
1351 //______________________________________________________________________
1352 void AliITS::GetTreeC(Int_t event){
1353     //    Get the clusters tree for this event and set the branch address.
1354     // Inputs:
1355     //      Int_t event    Event number for the cluster tree.
1356     // Outputs:
1357     //      none.
1358     // Return:
1359     //      none.
1360     char branchname[30];
1361     const char *det[3] = {"SPD","SDD","SSD"};
1362
1363     AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;
1364     TTree * lTC = pITSLoader->TreeC();
1365
1366     ResetClusters();
1367     if (lTC) {
1368         pITSLoader->CleanRawClusters();
1369     } // end if TreeC()
1370
1371     TBranch *branch;
1372
1373     if (lTC) {
1374         Int_t i;
1375         char digclass[40];
1376         char clclass[40];
1377         for (i=0; i<kNTYPES; i++) {
1378             AliITSDetType *iDetType=DetType(i); 
1379             iDetType->GetClassNames(digclass,clclass);
1380             // clusters
1381             if (fCtype == 0x0) fCtype  = new TObjArray(fNDetTypes);
1382             if(!fCtype->At(i)) fCtype->AddAt(new TClonesArray(clclass,1000),i);
1383             if(kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
1384             else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
1385             if (fCtype) {
1386                 branch = lTC->GetBranch(branchname);
1387                 if (branch) branch->SetAddress(&((*fCtype)[i]));
1388             } // end if fCtype
1389         } // end for i
1390     } else {
1391         Error("GetTreeC","cannot find Clusters Tree for event:%d",event);
1392     } // end if lTC
1393 }
1394 //______________________________________________________________________
1395 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c){
1396     //   Add a cluster to the list.
1397     // Inputs:
1398     //      Int_t id             Detector type number.
1399     //      AliITSRawCluster *c  Cluster class to be added to the tree of
1400     //                           clusters.
1401     // Outputs:
1402     //      none.
1403     // Return:
1404     //      none.
1405
1406     TClonesArray &lc = *((TClonesArray*)fCtype->At(id));
1407
1408     switch(id){
1409     case 0:
1410         new(lc[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
1411         break;
1412     case 1:
1413         new(lc[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
1414         break;
1415     case 2:
1416         new(lc[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
1417         break;
1418     } // end switch id
1419 }
1420 //______________________________________________________________________
1421 void AliITS::ResetClusters(Int_t i){
1422     //    Reset number of clusters and the clusters array for this branch.
1423     // Inputs:
1424     //      Int_t i        Detector type number.
1425     // Outputs:
1426     //      none.
1427     // Return:
1428     //      none.
1429
1430     if (fCtype->At(i))    ((TClonesArray*)fCtype->At(i))->Clear();
1431     if (fNctype)  fNctype[i]=0;
1432 }
1433 //______________________________________________________________________
1434 void AliITS::MakeBranchR(const char *file, Option_t *opt){
1435     // Creates Tree branches for the ITS Reconstructed points.
1436     // Inputs:
1437     //      cont char *file  File name where RecPoints branch is to be written
1438     //                       to. If blank it write the SDigits to the same
1439     //                       file in which the Hits were found.
1440     // Outputs:
1441     //      none.
1442     // Return:
1443     //      none.
1444     Int_t buffsz = 4000;
1445     char branchname[30];
1446
1447     // only one branch for rec points for all detector types
1448     Bool_t oFast= (strstr(opt,"Fast")!=0);
1449     if(oFast){
1450         sprintf(branchname,"%sRecPointsF",GetName());
1451     } else {
1452         sprintf(branchname,"%sRecPoints",GetName());
1453     }
1454
1455     if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1456     if (fLoader->TreeR()) {
1457         if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",
1458                                                           1000);
1459         MakeBranchInTree(fLoader->TreeR(),branchname,&fRecPoints,buffsz,file);
1460     } // end if
1461 }
1462 //______________________________________________________________________
1463 void AliITS::SetTreeAddressR(TTree *treeR){
1464     // Set branch address for the Reconstructed points Trees.
1465     // Inputs:
1466     //      TTree *treeR   Tree containing the RecPoints.
1467     // Outputs:
1468     //      none.
1469     // Return:
1470     //      none.
1471     char branchname[30];
1472
1473     if(!treeR) return;
1474     if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1475     TBranch *branch;
1476     sprintf(branchname,"%sRecPoints",GetName());
1477     branch = treeR->GetBranch(branchname);
1478     if (branch) {
1479         branch->SetAddress(&fRecPoints);
1480     }else {
1481         sprintf(branchname,"%sRecPointsF",GetName());
1482         branch = treeR->GetBranch(branchname);
1483         if (branch) {
1484             branch->SetAddress(&fRecPoints);
1485         }
1486     }
1487 }
1488 //______________________________________________________________________
1489 void AliITS::AddRecPoint(const AliITSRecPoint &r){
1490     // Add a reconstructed space point to the list
1491     // Inputs:
1492     //      const AliITSRecPoint &r RecPoint class to be added to the tree
1493     //                              of reconstructed points TreeR.
1494     // Outputs:
1495     //      none.
1496     // Return:
1497     //      none.
1498
1499     TClonesArray &lrecp = *fRecPoints;
1500     new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
1501 }
1502 //______________________________________________________________________
1503 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1504                                   Option_t *opt0,Option_t *opt1,
1505                                  const char *flnm){
1506     // keep galice.root for signal and name differently the file for 
1507     // background when add! otherwise the track info for signal will be lost !
1508     // the condition below will disappear when the geom class will be
1509     // initialized for all versions - for the moment it is only for v5 !
1510     // Inputs:
1511     //      Int_t evnt       Event to be processed.
1512     //      Int_t bgrev      Background Hit tree number.
1513     //      Int_t size       Size used by InitModules. See InitModules.
1514     //      Option_t *opt0   Option passed to FillModules. See FillModules.
1515     //      Option_t *opt1   String indicating if merging hits or not. To
1516     //                       merge hits set equal to "Add". Otherwise no
1517     //                       background hits are considered.
1518     //      Test_t *flnm     File name containing the background hits..
1519     // Outputs:
1520     //      none.
1521     // Return:
1522     //      none.
1523
1524     if(!GetITSgeom()) return;
1525     AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
1526     AliITSgeom *geom = GetITSgeom();
1527
1528     const char *all = strstr(opt1,"All");
1529     const char *det[3] ={strstr(opt1,"SPD"),strstr(opt1,"SDD"),
1530                          strstr(opt1,"SSD")};
1531     Int_t nmodules;
1532     InitModules(size,nmodules);
1533     FillModules(evNumber,bgrev,nmodules,opt0,flnm);
1534
1535     AliITSsimulation *sim      = 0;
1536     AliITSDetType    *iDetType = 0;
1537     AliITSmodule     *mod      = 0;
1538     Int_t id,module;
1539
1540     //m.b. : this change is nothing but a nice way to make sure
1541     //the CPU goes up !
1542     
1543     if(GetDebug()) cout<<"HitsToFastRecPoints: N mod = "<<
1544                        geom->GetIndexMax()<<endl;
1545     for(module=0;module<geom->GetIndexMax();module++){
1546         id       = geom->GetModuleType(module);
1547         if (!all && !det[id]) continue;
1548         iDetType = DetType(id);
1549         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1550         if (!sim) {
1551             Error("HitsToFastPoints","The simulation class was not "
1552                   "instanciated for module %d type %x!",module,
1553                   geom->GetModuleTypeName(module));
1554             exit(1);
1555         } // end if !sim
1556         mod      = (AliITSmodule *)fITSmodules->At(module);
1557         sim->CreateFastRecPoints(mod,module,gRandom);
1558         cout<<module<<"\r";fflush(0);
1559         //gAlice->TreeR()->Fill();
1560         TTree *lTR = pITSloader->TreeR();
1561         TBranch *br=lTR->GetBranch("ITSRecPointsF");
1562         br->Fill();
1563         ResetRecPoints();
1564     } // end for module
1565
1566     ClearModules();
1567     fLoader->WriteRecPoints("OVERWRITE");
1568 }
1569 //______________________________________________________________________
1570 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt){
1571   // cluster finding and reconstruction of space points
1572   // the condition below will disappear when the geom class will be
1573   // initialized for all versions - for the moment it is only for v5 !
1574   // 7 is the SDD beam test version
1575   // Inputs:
1576   //      Int_t evNumber   Event number to be processed.
1577   //      Int_t lastentry  Offset for module when not all of the modules
1578   //                       are processed.
1579   //      Option_t *opt    String indicating which ITS sub-detectors should
1580   //                       be processed. If ="All" then all of the ITS
1581   //                       sub detectors are processed.
1582   // Outputs:
1583   //      none.
1584   // Return:
1585   //      none.
1586
1587   if(!GetITSgeom()) return;
1588   AliITSgeom *geom = GetITSgeom();
1589     
1590   const char *all = strstr(opt,"All");
1591   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1592                         strstr(opt,"SSD")};
1593   static Bool_t setRec=kTRUE;
1594   if (setRec) SetDefaultClusterFinders();
1595   setRec=kFALSE;
1596
1597   AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
1598   TTree *treeC=pITSloader->TreeC();
1599   AliITSClusterFinder *rec     = 0;
1600   AliITSDetType      *iDetType = 0;
1601   Int_t id,module,first=0;
1602   for(module=0;module<geom->GetIndexMax();module++){
1603       id       = geom->GetModuleType(module);
1604       if (!all && !det[id]) continue;
1605       if(det[id]) first = geom->GetStartDet(id);
1606       iDetType = DetType(id);
1607       rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1608       TClonesArray *itsDigits  = this->DigitsAddress(id);
1609       if (!rec) {
1610           Error("DigitsToRecPoints",
1611                 "The reconstruction class was not instanciated! event=%d",
1612                 evNumber);
1613           exit(1);
1614       } // end if !rec
1615       this->ResetDigits();
1616       TTree *lTD = pITSloader->TreeD();
1617       if (all) {
1618           lTD->GetEvent(lastentry+module);
1619       }else {
1620           lTD->GetEvent(lastentry+(module-first));
1621       }
1622       Int_t ndigits = itsDigits->GetEntriesFast();
1623       if(ndigits>0){
1624           rec->SetDigits(DigitsAddress(id));
1625           rec->SetClusters(ClustersAddress(id));
1626           rec->FindRawClusters(module);
1627       } // end if
1628       pITSloader->TreeR()->Fill(); 
1629       ResetRecPoints();
1630       treeC->Fill();
1631       ResetClusters();
1632   } // end for module
1633
1634   pITSloader->WriteRecPoints("OVERWRITE");
1635   pITSloader->WriteRawClusters("OVERWRITE");
1636 }
1637 //______________________________________________________________________
1638 AliLoader* AliITS::MakeLoader(const char* topfoldername){ 
1639     //builds ITSgetter (AliLoader type)
1640     //if detector wants to use castomized getter, it must overload this method
1641
1642     Info("MakeLoader","Creating AliITSLoader. Top folder is %s.",
1643          topfoldername);
1644     fLoader = new AliITSLoader(GetName(),topfoldername);
1645     return fLoader;
1646 }