]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITS.cxx
Bug fix (C.Cheshkov)
[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     fLoader->UnloadSDigits();
910 }
911 //______________________________________________________________________
912 void AliITS::SDigitsToDigits(Option_t *opt){
913     // Standard Summable digits to Digits function.
914     // Inputs:
915     //      none.
916     // Outputs:
917     //      none.
918     char name[20] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
919
920     if(!GetITSgeom()) return; // need transformations to do digitization.
921     AliITSgeom *geom = GetITSgeom();
922
923     const char *all = strstr(opt,"All");
924     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
925                           strstr(opt,"SSD")};
926     if( !det[0] && !det[1] && !det[2] ) all = "All";
927     else all = 0;
928     static Bool_t setDef = kTRUE;
929     if(setDef) SetDefaultSimulation();
930     setDef = kFALSE;
931
932     AliITSsimulation *sim      = 0;
933     AliITSDetType    *iDetType = 0;
934     TTree            *trees    = fLoader->TreeS();
935     if( !(trees && this->GetSDigits()) ){
936         Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
937         return;
938     } // end if
939     sprintf( name, "%s", this->GetName() );
940     TBranch *brchSDigits = trees->GetBranch( name );
941     
942     Int_t id,module;
943     for(module=0;module<geom->GetIndexMax();module++){
944         id       = geom->GetModuleType(module);
945         if (!all && !det[id]) continue;
946         iDetType = DetType(id);
947         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
948         if (!sim) {
949             Error("SDigit2Digits","The simulation class was not "
950                   "instanciated for module %d type %s!",module,
951                   geom->GetModuleTypeName(module));
952             exit(1);
953         } // end if !sim
954         sim->InitSimulationModule(module,gAlice->GetEvNumber());
955         //
956         // add summable digits to module
957         this->GetSDigits()->Clear();
958         brchSDigits->GetEvent(module);
959         sim->AddSDigitsToModule(GetSDigits(),0);
960         //
961         // Digitise current module sum(SDigits)->Digits
962         sim->FinishSDigitiseModule();
963
964         // fills all branches - wasted disk space
965         fLoader->TreeD()->Fill();
966         this->ResetDigits();
967     } // end for module
968
969     fLoader->TreeD()->GetEntries();
970
971     fLoader->TreeD()->AutoSave();
972     // reset tree
973     fLoader->TreeD()->Reset();
974     
975 }
976 //______________________________________________________________________
977 void AliITS::Hits2Digits(){
978     // Standard Hits to Digits function.
979     // Inputs:
980     //      none.
981     // Outputs:
982     //      none.
983
984     fLoader->LoadHits("read");
985     fLoader->LoadDigits("recreate");
986     AliRunLoader* rl = fLoader->GetRunLoader(); 
987
988     for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
989         // Do the Hits to Digits operation. Use Standard input values.
990         // Event number from file, no background hit merging , use size from
991         // AliITSgeom class, option="All", input from this file only.
992         rl->GetEvent(iEvent);
993         if (!fLoader->TreeD()) fLoader->MakeTree("D");
994         MakeBranch("D");
995         SetTreeAddress();   
996         HitsToDigits(iEvent,0,-1," ",fOpt," ");
997     } // end for iEvent
998
999     fLoader->UnloadHits();
1000     fLoader->UnloadDigits();
1001 }
1002 //______________________________________________________________________
1003 void AliITS::HitsToPreDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1004                              Option_t *option,Option_t *opt,
1005                              const char *filename){
1006     //   Keep galice.root for signal and name differently the file for 
1007     // background when add! otherwise the track info for signal will be lost !
1008     // the condition below will disappear when the geom class will be
1009     // initialized for all versions - for the moment it is only for v5 !
1010     // 7 is the SDD beam test version.
1011     // Inputs:
1012     //      Int_t evnt       Event to be processed.
1013     //      Int_t bgrev      Background Hit tree number.
1014     //      Int_t nmodules   Not used.
1015     //      Option_t *option String indicating if merging hits or not. To
1016     //                       merge hits set equal to "Add". Otherwise no
1017     //                       background hits are considered.
1018     //      Test_t *filename File name containing the background hits..
1019     // Outputs:
1020     //      none.
1021     // Return:
1022     //      none.
1023
1024     if(!GetITSgeom()) return; // need transformations to do digitization.
1025     AliITSgeom *geom = GetITSgeom();
1026
1027     const char *all = strstr(opt,"All");
1028     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1029                           strstr(opt,"SSD")};
1030     static Bool_t setDef=kTRUE;
1031     if (setDef) SetDefaultSimulation();
1032     setDef=kFALSE;
1033
1034     Int_t nmodules;
1035     InitModules(size,nmodules);
1036     FillModules(evNumber,bgrev,nmodules,option,filename);
1037
1038     AliITSsimulation *sim      = 0;
1039     AliITSDetType    *iDetType = 0;
1040     AliITSmodule     *mod      = 0;
1041     Int_t id,module;
1042     for(module=0;module<geom->GetIndexMax();module++){
1043         id       = geom->GetModuleType(module);
1044         if (!all && !det[id]) continue;
1045         iDetType = DetType(id);
1046         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1047         if (!sim) {
1048             Error("HitsToSDigits","The simulation class was not "
1049                   "instanciated for module %d type %s!",module,
1050                   geom->GetModuleTypeName(module));
1051             exit(1);
1052         } // end if !sim
1053         mod      = (AliITSmodule *)fITSmodules->At(module);
1054         sim->SDigitiseModule(mod,module,evNumber);
1055         // fills all branches - wasted disk space
1056         fLoader->TreeS()->Fill(); 
1057         ResetSDigits();
1058     } // end for module
1059
1060     ClearModules();
1061
1062     fLoader->TreeS()->GetEntries();
1063     fLoader->TreeS()->AutoSave();
1064     fLoader->WriteSDigits("OVERWRITE");
1065     // reset tree
1066     fLoader->TreeS()->Reset();
1067 }
1068 //______________________________________________________________________
1069 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1070                           Option_t *option,Option_t *opt,
1071                           const char *filename){
1072     //   Keep galice.root for signal and name differently the file for 
1073     // background when add! otherwise the track info for signal will be lost !
1074     // the condition below will disappear when the geom class will be
1075     // initialized for all versions - for the moment it is only for v5 !
1076     // 7 is the SDD beam test version.
1077     // Inputs:
1078     //      Int_t evnt       Event to be processed.
1079     //      Int_t bgrev      Background Hit tree number.
1080     //      Int_t nmodules   Not used.
1081     //      Option_t *option String indicating if merging hits or not. To
1082     //                       merge hits set equal to "Add". Otherwise no
1083     //                       background hits are considered.
1084     //      Test_t *filename File name containing the background hits..
1085     // Outputs:
1086     //      none.
1087     // Return:
1088     //      none.
1089
1090     if(!GetITSgeom()) return; // need transformations to do digitization.
1091     AliITSgeom *geom = GetITSgeom();
1092
1093     const char *all = strstr(opt,"All");
1094     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1095                           strstr(opt,"SSD")};
1096     static Bool_t setDef=kTRUE;
1097     if (setDef) SetDefaultSimulation();
1098     setDef=kFALSE;
1099
1100     Int_t nmodules;
1101     InitModules(size,nmodules);
1102     FillModules(evNumber,bgrev,nmodules,option,filename);
1103
1104     AliITSsimulation *sim      = 0;
1105     AliITSDetType    *iDetType = 0;
1106     AliITSmodule     *mod      = 0;
1107     Int_t id,module;
1108     for(module=0;module<geom->GetIndexMax();module++){
1109         id       = geom->GetModuleType(module);
1110         if (!all && !det[id]) continue;
1111         iDetType = DetType(id);
1112         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1113         if (!sim) {
1114             Error("HitsToDigits","The simulation class was not "
1115                   "instanciated for module %d type %s!",module,
1116                   geom->GetModuleTypeName(module));
1117             exit(1);
1118         } // end if !sim
1119         mod      = (AliITSmodule *)fITSmodules->At(module);
1120         sim->DigitiseModule(mod,module,evNumber);
1121         // fills all branches - wasted disk space
1122         fLoader->TreeD()->Fill(); 
1123         ResetDigits();
1124     } // end for module
1125
1126     ClearModules();
1127
1128     fLoader->TreeD()->GetEntries();
1129     fLoader->TreeD()->AutoSave();
1130     // reset tree
1131     fLoader->TreeD()->Reset();
1132 }
1133 //______________________________________________________________________
1134 void AliITS::ResetDigits(){
1135     // Reset number of digits and the digits array for the ITS detector.
1136     // Inputs:
1137     //      none.
1138     // Outputs:
1139     //      none.
1140
1141     if (!fDtype) return;
1142
1143     Int_t i;
1144     for (i=0;i<kNTYPES;i++ ) {
1145         ResetDigits(i);
1146     } // end for i
1147 }
1148 //______________________________________________________________________
1149 void AliITS::ResetDigits(Int_t i){
1150     // Reset number of digits and the digits array for this branch.
1151     // Inputs:
1152     //      none.
1153     // Outputs:
1154     //      none.
1155
1156     if (fDtype->At(i)) ((TClonesArray*)fDtype->At(i))->Clear();
1157     if (fNdtype)       fNdtype[i]=0;
1158 }
1159 //______________________________________________________________________
1160 void AliITS::AddSumDigit(AliITSpListItem &sdig){
1161     // Adds the a module full of summable digits to the summable digits tree.
1162     // Inputs:
1163     //      AliITSpListItem &sdig   SDigit to be added to SDigits tree.
1164     // Outputs:
1165     //      none.
1166     // Return:
1167     //      none.
1168
1169     TClonesArray &lsdig = *fSDigits;
1170     new(lsdig[fNSDigits++]) AliITSpListItem(sdig);
1171 }
1172 //______________________________________________________________________
1173 void AliITS::AddRealDigit(Int_t id, Int_t *digits){
1174     //   Add a real digit - as coming from data.
1175     // Inputs:
1176     //      Int_t id        Detector type number.
1177     //      Int_t *digits   Integer array containing the digits info. See 
1178     //                      AliITSdigit.h
1179     // Outputs:
1180     //      none.
1181     // Return:
1182     //      none.
1183
1184     TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1185     new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
1186 }
1187 //______________________________________________________________________
1188 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d){
1189     //    Add a simulated digit.
1190     // Inputs:
1191     //      Int_t id        Detector type number.
1192     //      AliITSdigit *d  Digit to be added to the Digits Tree. See 
1193     //                      AliITSdigit.h
1194     // Outputs:
1195     //      none.
1196     // Return:
1197     //      none.
1198
1199     TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1200
1201     switch(id){
1202     case 0:
1203         new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
1204         break;
1205     case 1:
1206         new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
1207         break;
1208     case 2:
1209         new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
1210         break;
1211     } // end switch id
1212 }
1213 //______________________________________________________________________
1214 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,
1215                          Int_t *hits,Float_t *charges){
1216   //   Add a simulated digit to the list.
1217   // Inputs:
1218   //      Int_t id        Detector type number.
1219   //      Float_t phys    Physics indicator. See AliITSdigits.h
1220   //      Int_t *digits   Integer array containing the digits info. See 
1221   //                      AliITSdigit.h
1222   //      Int_t *tracks   Integer array [AliITSdigitS?D::GetNTracks()] 
1223   //                      containing the track numbers that contributed to
1224   //                      this digit.
1225   //      Int_t *hits     Integer array [AliITSdigitS?D::GetNTracks()]
1226   //                      containing the hit numbers, from AliITSmodule, that
1227   //                      contributed to this digit.
1228   //      Float_t *charge Floating point array of the signals contributed
1229   //                      to this digit by each track.
1230   // Outputs:
1231   //      none.
1232   // Return:
1233   //      none.
1234
1235   TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1236   AliITSresponseSDD *resp = 0;
1237   switch(id){
1238   case 0:
1239       new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
1240       break;
1241   case 1:
1242       resp = (AliITSresponseSDD*)DetType(1)->GetResponseModel();
1243       new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,
1244                                                  hits,charges,resp);
1245       break;
1246   case 2:
1247       new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
1248       break;
1249   } // end switch id
1250 }
1251 //______________________________________________________________________
1252 void AliITS::Digits2Raw(){
1253     // convert digits of the current event to raw data
1254
1255   fLoader->LoadDigits();
1256   TTree* digits = fLoader->TreeD();
1257   if (!digits) {
1258       Error("Digits2Raw", "no digits tree");
1259       return;
1260   }
1261   SetTreeAddressD(digits);
1262
1263   AliITSDDLRawData rawWriter;
1264   //Verbose level
1265   // 0: Silent
1266   // 1: cout messages
1267   // 2: txt files with digits 
1268   //BE CAREFUL, verbose level 2 MUST be used only for debugging and
1269   //it is highly suggested to use this mode only for debugging digits files
1270   //reasonably small, because otherwise the size of the txt files can reach
1271   //quickly several MB wasting time and disk space.
1272   rawWriter.SetVerbose(0);
1273     
1274   //SILICON PIXEL DETECTOR
1275   Info("Digits2Raw", "Formatting raw data for SPD");
1276   rawWriter.RawDataSPD(digits->GetBranch("ITSDigitsSPD"));
1277     
1278   //SILICON DRIFT DETECTOR
1279   Info("Digits2Raw", "Formatting raw data for SDD");
1280   rawWriter.RawDataSDD(digits->GetBranch("ITSDigitsSDD"));
1281     
1282   //SILICON STRIP DETECTOR
1283   Info("Digits2Raw", "Formatting raw data for SSD");
1284   rawWriter.RawDataSSD(digits->GetBranch("ITSDigitsSSD"));
1285
1286   fLoader->UnloadDigits();
1287 }
1288
1289 //______________________________________________________________________
1290 void AliITS::MakeTreeC(Option_t *option){
1291   //   Create a separate tree to store the clusters.
1292   // Inputs:
1293   //      Option_t *option  string which must contain "C" otherwise
1294   //                        no Cluster Tree is created.
1295   // Outputs:
1296   //      none.
1297   // Return:
1298   //      none.
1299
1300   AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;    
1301     
1302   if (pITSLoader == 0x0) {
1303       Error("MakeTreeC","fLoader == 0x0 option=%s",option);
1304       return;
1305   }
1306   if (pITSLoader->TreeC() == 0x0) pITSLoader->MakeTree("C");
1307   MakeBranchC();
1308 }
1309 //----------------------------------------------------------------------
1310 void AliITS::MakeBranchC(){
1311     //Makes barnches in treeC
1312     AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;    
1313     if (pITSLoader == 0x0) {
1314         Error("MakeTreeC","fLoader == 0x0");
1315         return;
1316     }
1317     TTree * lTC = pITSLoader->TreeC();
1318     if (lTC == 0x0){
1319         Error("MakeTreeC","Can not get TreeC from Loader");
1320         return;
1321     }
1322
1323     Int_t buffersize = 4000;
1324     char branchname[30];
1325     const char *det[3] = {"SPD","SDD","SSD"};
1326     char digclass[40];
1327     char clclass[40];
1328     
1329     // one branch for Clusters per type of detector
1330     Int_t i;   
1331     for (i=0; i<kNTYPES ;i++) {
1332         AliITSDetType *iDetType=DetType(i); 
1333         iDetType->GetClassNames(digclass,clclass);
1334         // clusters
1335         if (fCtype == 0x0) fCtype  = new TObjArray(fNDetTypes);
1336         if(!ClustersAddress(i)){
1337             fCtype->AddAt(new TClonesArray(clclass,1000),i);
1338         }
1339         if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
1340         else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
1341         if (fCtype  && lTC) {
1342             if (lTC->GetBranch(branchname)){
1343                 Warning("MakeBranchC","Branch %s alread exists in TreeC",
1344                         branchname);
1345             }else{
1346                 Info("MakeBranchC","Creating branch %s in TreeC",branchname);
1347                 lTC->Branch(branchname,&((*fCtype)[i]), buffersize);
1348             }
1349         } // end if fCtype && lTC
1350     } // end for i
1351 }
1352 //______________________________________________________________________
1353 void AliITS::GetTreeC(Int_t event){
1354     //    Get the clusters tree for this event and set the branch address.
1355     // Inputs:
1356     //      Int_t event    Event number for the cluster tree.
1357     // Outputs:
1358     //      none.
1359     // Return:
1360     //      none.
1361     char branchname[30];
1362     const char *det[3] = {"SPD","SDD","SSD"};
1363
1364     AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;
1365     TTree * lTC = pITSLoader->TreeC();
1366
1367     ResetClusters();
1368     if (lTC) {
1369         pITSLoader->CleanRawClusters();
1370     } // end if TreeC()
1371
1372     TBranch *branch;
1373
1374     if (lTC) {
1375         Int_t i;
1376         char digclass[40];
1377         char clclass[40];
1378         for (i=0; i<kNTYPES; i++) {
1379             AliITSDetType *iDetType=DetType(i); 
1380             iDetType->GetClassNames(digclass,clclass);
1381             // clusters
1382             if (fCtype == 0x0) fCtype  = new TObjArray(fNDetTypes);
1383             if(!fCtype->At(i)) fCtype->AddAt(new TClonesArray(clclass,1000),i);
1384             if(kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
1385             else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
1386             if (fCtype) {
1387                 branch = lTC->GetBranch(branchname);
1388                 if (branch) branch->SetAddress(&((*fCtype)[i]));
1389             } // end if fCtype
1390         } // end for i
1391     } else {
1392         Error("GetTreeC","cannot find Clusters Tree for event:%d",event);
1393     } // end if lTC
1394 }
1395 //______________________________________________________________________
1396 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c){
1397     //   Add a cluster to the list.
1398     // Inputs:
1399     //      Int_t id             Detector type number.
1400     //      AliITSRawCluster *c  Cluster class to be added to the tree of
1401     //                           clusters.
1402     // Outputs:
1403     //      none.
1404     // Return:
1405     //      none.
1406
1407     TClonesArray &lc = *((TClonesArray*)fCtype->At(id));
1408
1409     switch(id){
1410     case 0:
1411         new(lc[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
1412         break;
1413     case 1:
1414         new(lc[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
1415         break;
1416     case 2:
1417         new(lc[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
1418         break;
1419     } // end switch id
1420 }
1421 //______________________________________________________________________
1422 void AliITS::ResetClusters(Int_t i){
1423     //    Reset number of clusters and the clusters array for this branch.
1424     // Inputs:
1425     //      Int_t i        Detector type number.
1426     // Outputs:
1427     //      none.
1428     // Return:
1429     //      none.
1430
1431     if (fCtype->At(i))    ((TClonesArray*)fCtype->At(i))->Clear();
1432     if (fNctype)  fNctype[i]=0;
1433 }
1434 //______________________________________________________________________
1435 void AliITS::MakeBranchR(const char *file, Option_t *opt){
1436     // Creates Tree branches for the ITS Reconstructed points.
1437     // Inputs:
1438     //      cont char *file  File name where RecPoints branch is to be written
1439     //                       to. If blank it write the SDigits to the same
1440     //                       file in which the Hits were found.
1441     // Outputs:
1442     //      none.
1443     // Return:
1444     //      none.
1445     Int_t buffsz = 4000;
1446     char branchname[30];
1447
1448     // only one branch for rec points for all detector types
1449     Bool_t oFast= (strstr(opt,"Fast")!=0);
1450     if(oFast){
1451         sprintf(branchname,"%sRecPointsF",GetName());
1452     } else {
1453         sprintf(branchname,"%sRecPoints",GetName());
1454     }
1455
1456     if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1457     if (fLoader->TreeR()) {
1458         if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",
1459                                                           1000);
1460         MakeBranchInTree(fLoader->TreeR(),branchname,&fRecPoints,buffsz,file);
1461     } // end if
1462 }
1463 //______________________________________________________________________
1464 void AliITS::SetTreeAddressR(TTree *treeR){
1465     // Set branch address for the Reconstructed points Trees.
1466     // Inputs:
1467     //      TTree *treeR   Tree containing the RecPoints.
1468     // Outputs:
1469     //      none.
1470     // Return:
1471     //      none.
1472     char branchname[30];
1473
1474     if(!treeR) return;
1475     if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1476     TBranch *branch;
1477     sprintf(branchname,"%sRecPoints",GetName());
1478     branch = treeR->GetBranch(branchname);
1479     if (branch) {
1480         branch->SetAddress(&fRecPoints);
1481     }else {
1482         sprintf(branchname,"%sRecPointsF",GetName());
1483         branch = treeR->GetBranch(branchname);
1484         if (branch) {
1485             branch->SetAddress(&fRecPoints);
1486         }
1487     }
1488 }
1489 //______________________________________________________________________
1490 void AliITS::AddRecPoint(const AliITSRecPoint &r){
1491     // Add a reconstructed space point to the list
1492     // Inputs:
1493     //      const AliITSRecPoint &r RecPoint class to be added to the tree
1494     //                              of reconstructed points TreeR.
1495     // Outputs:
1496     //      none.
1497     // Return:
1498     //      none.
1499
1500     TClonesArray &lrecp = *fRecPoints;
1501     new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
1502 }
1503 //______________________________________________________________________
1504 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1505                                   Option_t *opt0,Option_t *opt1,
1506                                  const char *flnm){
1507     // keep galice.root for signal and name differently the file for 
1508     // background when add! otherwise the track info for signal will be lost !
1509     // the condition below will disappear when the geom class will be
1510     // initialized for all versions - for the moment it is only for v5 !
1511     // Inputs:
1512     //      Int_t evnt       Event to be processed.
1513     //      Int_t bgrev      Background Hit tree number.
1514     //      Int_t size       Size used by InitModules. See InitModules.
1515     //      Option_t *opt0   Option passed to FillModules. See FillModules.
1516     //      Option_t *opt1   String indicating if merging hits or not. To
1517     //                       merge hits set equal to "Add". Otherwise no
1518     //                       background hits are considered.
1519     //      Test_t *flnm     File name containing the background hits..
1520     // Outputs:
1521     //      none.
1522     // Return:
1523     //      none.
1524
1525     if(!GetITSgeom()) return;
1526     AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
1527     AliITSgeom *geom = GetITSgeom();
1528
1529     const char *all = strstr(opt1,"All");
1530     const char *det[3] ={strstr(opt1,"SPD"),strstr(opt1,"SDD"),
1531                          strstr(opt1,"SSD")};
1532     Int_t nmodules;
1533     InitModules(size,nmodules);
1534     FillModules(evNumber,bgrev,nmodules,opt0,flnm);
1535
1536     AliITSsimulation *sim      = 0;
1537     AliITSDetType    *iDetType = 0;
1538     AliITSmodule     *mod      = 0;
1539     Int_t id,module;
1540
1541     //m.b. : this change is nothing but a nice way to make sure
1542     //the CPU goes up !
1543     
1544     if(GetDebug()) cout<<"HitsToFastRecPoints: N mod = "<<
1545                        geom->GetIndexMax()<<endl;
1546     for(module=0;module<geom->GetIndexMax();module++){
1547         id       = geom->GetModuleType(module);
1548         if (!all && !det[id]) continue;
1549         iDetType = DetType(id);
1550         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1551         if (!sim) {
1552             Error("HitsToFastPoints","The simulation class was not "
1553                   "instanciated for module %d type %x!",module,
1554                   geom->GetModuleTypeName(module));
1555             exit(1);
1556         } // end if !sim
1557         mod      = (AliITSmodule *)fITSmodules->At(module);
1558         sim->CreateFastRecPoints(mod,module,gRandom);
1559         cout<<module<<"\r";fflush(0);
1560         //gAlice->TreeR()->Fill();
1561         TTree *lTR = pITSloader->TreeR();
1562         TBranch *br=lTR->GetBranch("ITSRecPointsF");
1563         br->Fill();
1564         ResetRecPoints();
1565     } // end for module
1566
1567     ClearModules();
1568     fLoader->WriteRecPoints("OVERWRITE");
1569 }
1570 //______________________________________________________________________
1571 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt){
1572   // cluster finding and reconstruction of space points
1573   // the condition below will disappear when the geom class will be
1574   // initialized for all versions - for the moment it is only for v5 !
1575   // 7 is the SDD beam test version
1576   // Inputs:
1577   //      Int_t evNumber   Event number to be processed.
1578   //      Int_t lastentry  Offset for module when not all of the modules
1579   //                       are processed.
1580   //      Option_t *opt    String indicating which ITS sub-detectors should
1581   //                       be processed. If ="All" then all of the ITS
1582   //                       sub detectors are processed.
1583   // Outputs:
1584   //      none.
1585   // Return:
1586   //      none.
1587
1588   if(!GetITSgeom()) return;
1589   AliITSgeom *geom = GetITSgeom();
1590     
1591   const char *all = strstr(opt,"All");
1592   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1593                         strstr(opt,"SSD")};
1594   static Bool_t setRec=kTRUE;
1595   if (setRec) SetDefaultClusterFinders();
1596   setRec=kFALSE;
1597
1598   AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
1599   TTree *treeC=pITSloader->TreeC();
1600   AliITSClusterFinder *rec     = 0;
1601   AliITSDetType      *iDetType = 0;
1602   Int_t id,module,first=0;
1603   for(module=0;module<geom->GetIndexMax();module++){
1604       id       = geom->GetModuleType(module);
1605       if (!all && !det[id]) continue;
1606       if(det[id]) first = geom->GetStartDet(id);
1607       iDetType = DetType(id);
1608       rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1609       TClonesArray *itsDigits  = this->DigitsAddress(id);
1610       if (!rec) {
1611           Error("DigitsToRecPoints",
1612                 "The reconstruction class was not instanciated! event=%d",
1613                 evNumber);
1614           exit(1);
1615       } // end if !rec
1616       this->ResetDigits();
1617       TTree *lTD = pITSloader->TreeD();
1618       if (all) {
1619           lTD->GetEvent(lastentry+module);
1620       }else {
1621           lTD->GetEvent(lastentry+(module-first));
1622       }
1623       Int_t ndigits = itsDigits->GetEntriesFast();
1624       if(ndigits>0){
1625           rec->SetDigits(DigitsAddress(id));
1626           rec->SetClusters(ClustersAddress(id));
1627           rec->FindRawClusters(module);
1628       } // end if
1629       pITSloader->TreeR()->Fill(); 
1630       ResetRecPoints();
1631       treeC->Fill();
1632       ResetClusters();
1633   } // end for module
1634
1635   pITSloader->WriteRecPoints("OVERWRITE");
1636   pITSloader->WriteRawClusters("OVERWRITE");
1637 }
1638 //______________________________________________________________________
1639 AliLoader* AliITS::MakeLoader(const char* topfoldername){ 
1640     //builds ITSgetter (AliLoader type)
1641     //if detector wants to use castomized getter, it must overload this method
1642
1643     Info("MakeLoader","Creating AliITSLoader. Top folder is %s.",
1644          topfoldername);
1645     fLoader = new AliITSLoader(GetName(),topfoldername);
1646     return fLoader;
1647 }