]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITS.cxx
Implementation of reconstructors (T.Kuhr)
[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     if(gAlice->GetMCApp()) gAlice->GetMCApp()->AddHitList(fHits);// Not done in AliDetector.
200
201     fEuclidOut  = 0;
202     fITSgeom    = 0;
203     fITSmodules = 0;
204     SetDetectors(); // default to fOpt="All". This variable not written out.
205
206     fIdN        = 0;
207     fIdName     = 0;
208     fIdSens     = 0;
209
210     fNDetTypes  = kNTYPES;
211     fDetTypes   = new TObjArray(fNDetTypes);
212
213     fSDigits    = new TClonesArray("AliITSpListItem",1000);
214     fNSDigits   = 0;
215
216     fNdtype     = new Int_t[fNDetTypes];
217     fDtype      = new TObjArray(fNDetTypes);
218
219     fCtype      = new TObjArray(fNDetTypes);
220     fNctype     = new Int_t[fNDetTypes];
221
222     fRecPoints  = new TClonesArray("AliITSRecPoint",1000);
223     fNRecPoints = 0;
224
225     Int_t i;
226     for(i=0;i<fNDetTypes;i++) {
227         fDetTypes->AddAt(new AliITSDetType(),i); 
228         fNdtype[i] = 0;
229         fNctype[i] = 0;
230     } // end for i
231
232     SetMarkerColor(kRed);
233     SelectVertexer(" ");
234 }
235 //______________________________________________________________________
236 AliITS::~AliITS(){
237     // Default destructor for ITS.
238     //     The default destructor of the AliITS class. In addition to deleting
239     // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
240     // fSDigits, fCtype, fITSmodules, fITSgeom, fRecPoints, fIdSens, fIdName, 
241     // fITSpoints, fDetType and it's contents.
242     // Inputs:
243     //      none.
244     // Outputs:
245     //      none.
246     // Return:
247     //      none.
248
249     if (fHits) {
250       fHits->Delete();
251       delete fHits;
252       fHits=0;
253     }
254     if (fSDigits) {
255       fSDigits->Delete();
256       delete fSDigits;
257       fSDigits=0;
258     }
259     if (fDigits) {
260       fDigits->Delete();
261       delete fDigits;
262       fDigits=0;
263     }
264     if (fRecPoints) {
265       fRecPoints->Delete();
266       delete fRecPoints;
267       fRecPoints=0;
268     }
269     delete[] fIdName;  // Array of TStrings
270     delete[] fIdSens;
271     if(fITSmodules) {
272         this->ClearModules();
273         delete fITSmodules;
274         fITSmodules = 0;
275     }// end if fITSmodules!=0
276
277     if(fDtype) {
278         fDtype->Delete();
279         delete fDtype;
280     } // end if fDtype
281     delete [] fNdtype;
282     if (fCtype) {
283         fCtype->Delete();
284         delete fCtype;
285         fCtype = 0;
286     } // end if fCtype
287     delete [] fNctype;
288
289     if (fDetTypes) {
290         fDetTypes->Delete();
291         delete fDetTypes;
292         fDetTypes = 0;
293     } // end if fDetTypes
294
295
296     if (fITSgeom) delete fITSgeom;
297 }
298 //______________________________________________________________________
299 AliITS::AliITS(const AliITS &source) : AliDetector(source){
300     // Copy constructor. This is a function which is not allowed to be
301     // done to the ITS. It exits with an error.
302     // Inputs:
303     //      AliITS &source  An AliITS class.
304     // Outputs:
305     //      none.
306     // Return:
307     //      none.
308
309     if(this==&source) return;
310     Error("Copy constructor",
311           "You are not allowed to make a copy of the AliITS");
312     exit(1);
313 }
314 //______________________________________________________________________
315 AliITS& AliITS::operator=(AliITS &source){
316     // Assignment operator. This is a function which is not allowed to be
317     // done to the ITS. It exits with an error.
318     // Inputs:
319     //      AliITS &source  An AliITS class.
320     // Outputs:
321     //      none.
322     // Return:
323     //      none.
324
325     if(this==&source) return *this;
326     Error("operator=","You are not allowed to make a copy of the AliITS");
327     exit(1);
328     return *this; //fake return
329 }
330 //______________________________________________________________________
331 Int_t AliITS::DistancetoPrimitive(Int_t,Int_t) const{
332     // Distance from mouse to ITS on the screen. Dummy routine
333     //     A dummy routine used by the ROOT macro display.C to allow for the
334     // use of the mouse (pointing device) in the macro. In general this should
335     // never be called. If it is it returns the number 9999 for any value of
336     // x and y.
337     // Inputs:
338     //      Int_t     Dummy screen coordinate.
339     //      Int_t     Dummy screen coordinate.
340     // Outputs:
341     //      none.
342     // Return:
343     //      Int_t     Dummy = 9999 distance to ITS.
344
345     return 9999;
346 }
347 //______________________________________________________________________
348 void AliITS::Init(){
349     // Initializer ITS after it has been built
350     //     This routine initializes the AliITS class. It is intended to be
351     // called from the Init function in AliITSv?. Besides displaying a banner
352     // indicating that it has been called it initializes the array fIdSens
353     // and sets the default segmentation, response, digit and raw cluster
354     // classes therefore it should be called after a call to CreateGeometry.
355     // Inputs:
356     //      none.
357     // Outputs:
358     //      none.
359     // Return:
360     //      none.
361     Int_t i;
362
363     SetDefaults();
364     // Array of TStrings
365     if(gMC) for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
366 }
367 //______________________________________________________________________
368 void AliITS::SetDefaults(){
369     // sets the default segmentation, response, digit and raw cluster classes.
370     // Inputs:
371     //      none.
372     // Outputs:
373     //      none.
374     // Return:
375     //      none.
376
377     if(fDebug) printf("%s: SetDefaults\n",ClassName());
378
379     AliITSDetType *iDetType;
380
381     //SPD
382     iDetType=DetType(0); 
383     if (!iDetType->GetSegmentationModel()) {
384         AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
385         SetSegmentationModel(0,seg0); 
386     } // end if
387     if (!iDetType->GetResponseModel()) {
388         SetResponseModel(0,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
397     // SDD
398     iDetType=DetType(1); 
399     if (!iDetType->GetResponseModel()) {
400         SetResponseModel(1,new AliITSresponseSDD("simulated")); 
401     } // end if
402     AliITSresponse *resp1=iDetType->GetResponseModel();
403     if (!iDetType->GetSegmentationModel()) {
404         AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
405         SetSegmentationModel(1,seg1); 
406     } // end if
407     const char *kData1=(iDetType->GetResponseModel())->DataType();
408     const char *kopt=iDetType->GetResponseModel()->ZeroSuppOption();
409     if((!strstr(kopt,"2D"))&&(!strstr(kopt,"1D")) || strstr(kData1,"real") ){
410         iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
411     } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
412
413     // SSD
414     iDetType=DetType(2); 
415     if (!iDetType->GetSegmentationModel()) {
416         AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
417         SetSegmentationModel(2,seg2); 
418     } // end if
419     if (!iDetType->GetResponseModel()) {
420         SetResponseModel(2,new AliITSresponseSSD("simulated"));
421     } // end if
422     const char *kData2=(iDetType->GetResponseModel())->DataType();
423     if (strstr(kData2,"real")) {
424         iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
425     } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
426
427     if (kNTYPES>3) {
428         Warning("SetDefaults",
429                 "Only the three basic detector types are initialized!");
430     }  // end if
431 }
432 //______________________________________________________________________
433 void AliITS::SetDefaultSimulation(){
434     // sets the default simulation.
435     // Inputs:
436     //      none.
437     // Outputs:
438     //      none.
439     // Return:
440     //      none.
441
442     AliITSDetType *iDetType;
443     AliITSsimulation *sim;
444     iDetType=DetType(0);
445     sim = iDetType->GetSimulationModel();
446     if (!sim) {
447         AliITSsegmentation *seg0=
448             (AliITSsegmentation*)iDetType->GetSegmentationModel();
449         AliITSresponse *res0 = (AliITSresponse*)iDetType->GetResponseModel();
450         AliITSsimulationSPD *sim0=new AliITSsimulationSPD(seg0,res0);
451         SetSimulationModel(0,sim0);
452     }else{ // simulation exists, make sure it is set up properly.
453         ((AliITSsimulationSPD*)sim)->Init(
454             (AliITSsegmentationSPD*) iDetType->GetSegmentationModel(),
455             (AliITSresponseSPD*) iDetType->GetResponseModel());
456 //        if(sim->GetResponseModel()==0) sim->SetResponseModel(
457 //            (AliITSresponse*)iDetType->GetResponseModel());
458 //        if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
459 //            (AliITSsegmentation*)iDetType->GetSegmentationModel());
460     } // end if
461     iDetType=DetType(1);
462     sim = iDetType->GetSimulationModel();
463     if (!sim) {
464         AliITSsegmentation *seg1=
465             (AliITSsegmentation*)iDetType->GetSegmentationModel();
466         AliITSresponse *res1 = (AliITSresponse*)iDetType->GetResponseModel();
467         AliITSsimulationSDD *sim1=new AliITSsimulationSDD(seg1,res1);
468         SetSimulationModel(1,sim1);
469     }else{ // simulation exists, make sure it is set up properly.
470         ((AliITSsimulationSDD*)sim)->Init(
471             (AliITSsegmentationSDD*) iDetType->GetSegmentationModel(),
472             (AliITSresponseSDD*) iDetType->GetResponseModel());
473 //        if(sim->GetResponseModel()==0) sim->SetResponseModel(
474 //            (AliITSresponse*)iDetType->GetResponseModel());
475 //        if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
476 //            (AliITSsegmentation*)iDetType->GetSegmentationModel());
477     } //end if
478     iDetType=DetType(2);
479     sim = iDetType->GetSimulationModel();
480     if (!sim) {
481         AliITSsegmentation *seg2=
482             (AliITSsegmentation*)iDetType->GetSegmentationModel();
483         AliITSresponse *res2 = (AliITSresponse*)iDetType->GetResponseModel();
484         AliITSsimulationSSD *sim2=new AliITSsimulationSSD(seg2,res2);
485         SetSimulationModel(2,sim2);
486     }else{ // simulation exists, make sure it is set up properly.
487         ((AliITSsimulationSSD*)sim)->Init(
488             (AliITSsegmentationSSD*) iDetType->GetSegmentationModel(),
489             (AliITSresponseSSD*) iDetType->GetResponseModel());
490 //        if(sim->GetResponseModel()==0) sim->SetResponseModel(
491 //            (AliITSresponse*)iDetType->GetResponseModel());
492 //        if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
493 //            (AliITSsegmentation*)iDetType->GetSegmentationModel());
494     } // end if
495 }
496 //______________________________________________________________________
497 void AliITS::SetDefaultClusterFinders(){
498     // Sets the default cluster finders. Used in finding RecPoints.
499     // Inputs:
500     //      none.
501     // Outputs:
502     //      none.
503     // Return:
504     //      none.
505
506     MakeTreeC();
507     AliITSDetType *iDetType;
508
509     // SPD
510     iDetType=DetType(0);
511     if (!iDetType->GetReconstructionModel()) {
512         AliITSsegmentation *seg0 =
513             (AliITSsegmentation*)iDetType->GetSegmentationModel();
514         TClonesArray *dig0=DigitsAddress(0);
515         TClonesArray *recp0=ClustersAddress(0);
516         AliITSClusterFinderSPD *rec0 = new AliITSClusterFinderSPD(seg0,dig0,
517                                                                   recp0);
518         SetReconstructionModel(0,rec0);
519     } // end if
520
521     // SDD
522     iDetType=DetType(1);
523     if (!iDetType->GetReconstructionModel()) {
524         AliITSsegmentation *seg1 =
525             (AliITSsegmentation*)iDetType->GetSegmentationModel();
526         AliITSresponse *res1 = (AliITSresponse*)iDetType->GetResponseModel();
527         TClonesArray *dig1=DigitsAddress(1);
528         TClonesArray *recp1=ClustersAddress(1);
529         AliITSClusterFinderSDD *rec1 =
530             new AliITSClusterFinderSDD(seg1,res1,dig1,recp1);
531       SetReconstructionModel(1,rec1);
532     } // end if
533
534     // SSD
535     iDetType=DetType(2);
536     if (!iDetType->GetReconstructionModel()) {
537         AliITSsegmentation *seg2=
538             (AliITSsegmentation*)iDetType->GetSegmentationModel();
539         TClonesArray *dig2=DigitsAddress(2);
540         AliITSClusterFinderSSD *rec2= new AliITSClusterFinderSSD(seg2,dig2);
541         SetReconstructionModel(2,rec2);
542     } // end if
543 }
544 //______________________________________________________________________
545 void AliITS::MakeBranch(Option_t* option){
546     // Creates Tree branches for the ITS.
547     // Inputs:
548     //      Option_t *option    String of Tree types S,D, and/or R.
549     //      const char *file    String of the file name where these branches
550     //                          are to be stored. If blank then these branches
551     //                          are written to the same tree as the Hits were
552     //                          read from.
553     // Outputs:
554     //      none.
555     // Return:
556     //      none.
557     Bool_t cH = (strstr(option,"H")!=0);
558     Bool_t cS = (strstr(option,"S")!=0);
559     Bool_t cD = (strstr(option,"D")!=0);
560     Bool_t cR = (strstr(option,"R")!=0);
561     Bool_t cRF = (strstr(option,"RF")!=0);
562     
563     if(cRF)cR = kFALSE;
564     if(cH && (fHits == 0x0)) fHits  = new TClonesArray("AliITShit", 1560);
565     
566     AliDetector::MakeBranch(option);
567
568     if(cS) MakeBranchS(0);
569     if(cD) MakeBranchD(0);
570     if(cR) MakeBranchR(0);
571     if(cRF) MakeBranchRF(0);
572 }
573 //______________________________________________________________________
574 void AliITS::SetTreeAddress(){
575     // Set branch address for the Trees.
576     // Inputs:
577     //      none.
578     // Outputs:
579     //      none.
580     // Return:
581     //      none.
582     TTree *treeS = fLoader->TreeS();
583     TTree *treeD = fLoader->TreeD();
584     TTree *treeR = fLoader->TreeR();
585     if (fLoader->TreeH() && (fHits == 0x0)) fHits = new TClonesArray("AliITShit", 1560);
586       
587     AliDetector::SetTreeAddress();
588
589     SetTreeAddressS(treeS);
590     SetTreeAddressD(treeD);
591     SetTreeAddressR(treeR);
592 }
593 //______________________________________________________________________
594 AliITSDetType* AliITS::DetType(Int_t id){
595     // Return pointer to id detector type.
596     // Inputs:
597     //      Int_t id   detector id number.
598     // Outputs:
599     //      none.
600     // Return:
601     //      returned, a pointer to a AliITSDetType.
602
603     return ((AliITSDetType*) fDetTypes->At(id));
604 }
605 //______________________________________________________________________
606 void AliITS::SetResponseModel(Int_t id, AliITSresponse *response){
607     // Set the response model for the id detector type.
608     // Inputs:
609     //      Int_t id        detector id number.
610     //      AliITSresponse* a pointer containing an instance of AliITSresponse
611     //                      to be stored/owned b y AliITSDetType.
612     // Outputs:
613     //      none.
614     // Return:
615     //      none.
616
617     ((AliITSDetType*) fDetTypes->At(id))->ResponseModel(response);
618 }
619 //______________________________________________________________________
620 void AliITS::SetSegmentationModel(Int_t id, AliITSsegmentation *seg){
621     // Set the segmentation model for the id detector type.
622     // Inputs:
623     //      Int_t id            detector id number.
624     //      AliITSsegmentation* a pointer containing an instance of 
625     //                          AliITSsegmentation to be stored/owned b y 
626     //                          AliITSDetType.
627     // Outputs:
628     //      none.
629     // Return:
630     //      none.
631
632     ((AliITSDetType*) fDetTypes->At(id))->SegmentationModel(seg);
633 }
634 //______________________________________________________________________
635 void AliITS::SetSimulationModel(Int_t id, AliITSsimulation *sim){
636     // Set the simulation model for the id detector type.
637     // Inputs:
638     //      Int_t id        detector id number.
639     //      AliITSresponse* a pointer containing an instance of AliITSresponse
640     //                      to be stored/owned b y AliITSDetType.
641     // Outputs:
642     //      none.
643     // Return:
644     //      none.
645
646    ((AliITSDetType*) fDetTypes->At(id))->SimulationModel(sim);
647
648 }
649 //______________________________________________________________________
650 void AliITS::SetReconstructionModel(Int_t id, AliITSClusterFinder *reconst){
651     // Set the cluster finder model for the id detector type.
652     // Inputs:
653     //      Int_t id             detector id number.
654     //      AliITSClusterFinder* a pointer containing an instance of 
655     //                           AliITSClusterFinder to be stored/owned b y 
656     //                           AliITSDetType.
657     // Outputs:
658     //      none.
659     // Return:
660     //      none.
661
662     ((AliITSDetType*) fDetTypes->At(id))->ReconstructionModel(reconst);
663 }
664 //______________________________________________________________________
665 void AliITS::SetClasses(Int_t id, const char *digit, const char *cluster){
666     // Set the digit and cluster classes name to be used for the id detector
667     // type.
668     // Inputs:
669     //      Int_t id            detector id number.
670     //      const char *digit   Digit class name for detector id.
671     //      const char *cluster Cluster class name for detector id.
672     // Outputs:
673     //      none.
674     // Return:
675     //      none.
676
677     ((AliITSDetType*) fDetTypes->At(id))->ClassNames(digit,cluster);
678 }
679 //______________________________________________________________________
680 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
681     // Add an ITS hit
682     //     The function to add information to the AliITShit class. See the
683     // AliITShit class for a full description. This function allocates the
684     // necessary new space for the hit information and passes the variable
685     // track, and the pointers *vol and *hits to the AliITShit constructor
686     // function.
687     // Inputs:
688     //      Int_t   track   Track number which produced this hit.
689     //      Int_t   *vol    Array of Integer Hit information. See AliITShit.h
690     //      Float_t *hits   Array of Floating Hit information.  see AliITShit.h
691     // Outputs:
692     //      none.
693     // Return:
694     //      none.
695
696     TClonesArray &lhits = *fHits;
697     new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
698 }
699 //______________________________________________________________________
700 void AliITS::InitModules(Int_t size,Int_t &nmodules){
701     // Initialize the modules array.
702     // Inputs:
703     //      Int_t size  Size of array of the number of modules to be
704     //                  created. If size <=0 then the number of modules
705     //                  is gotten from AliITSgeom class kept in fITSgeom.
706     // Outputs:
707     //      Int_t &nmodules The number of modules existing.
708     // Return:
709     //      none.
710
711     if(fITSmodules){ 
712         fITSmodules->Delete();
713         delete fITSmodules;
714     } // end fir fITSmoudles
715
716     Int_t nl,indexMAX,index;
717
718     if(size<=0){ // default to using data stored in AliITSgeom
719         if(fITSgeom==0) {
720             Error("InitModules","fITSgeom not defined");
721             return;
722         } // end if fITSgeom==0
723         nl = fITSgeom->GetNlayers();
724         indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
725                                             fITSgeom->GetNdetectors(nl))+1;
726         nmodules = indexMAX;
727         fITSmodules = new TObjArray(indexMAX);
728         for(index=0;index<indexMAX;index++){
729             fITSmodules->AddAt( new AliITSmodule(index),index);
730         } // end for index
731     }else{
732         fITSmodules = new TObjArray(size);
733         for(index=0;index<size;index++) {
734             fITSmodules->AddAt( new AliITSmodule(index),index);
735         } // end for index
736
737         nmodules = size;
738     } // end i size<=0
739 }
740 //______________________________________________________________________
741 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,
742                          Option_t *option, const char *filename){
743     // fill the modules with the sorted by module hits; add hits from
744     // background if option=Add.
745     // Inputs:
746     //      Int_t evnt       Event to be processed.
747     //      Int_t bgrev      Background Hit tree number.
748     //      Int_t nmodules   Not used.
749     //      Option_t *option String indicating if merging hits or not. To
750     //                       merge hits set equal to "Add". Otherwise no
751     //                       background hits are considered.
752     //      Test_t *filename File name containing the background hits..
753     // Outputs:
754     //      none.
755     // Return:
756     //      none.
757     static TTree *trH1;                 //Tree with background hits
758     static Bool_t first=kTRUE;
759     static TFile *file;
760     const char *addBgr = strstr(option,"Add");
761
762     evnt = nmodules; // Dummy use of variables to remove warnings
763     if (addBgr ) {
764         if(first) {
765             file=new TFile(filename);
766         } // end if first
767         first=kFALSE;
768         file->cd();
769         file->ls();
770         // Get Hits Tree header from file
771         if(trH1) delete trH1;
772         trH1=0;
773
774         char treeName[20];
775         sprintf(treeName,"TreeH%d",bgrev);
776         trH1 = (TTree*)gDirectory->Get(treeName);
777         if (!trH1) {
778             Error("FillModules","cannot find Hits Tree for event:%d",bgrev);
779         } // end if !trH1
780         // Set branch addresses
781     } // end if addBgr
782
783     FillModules(fLoader->TreeH(),0); // fill from this file's tree.
784     
785     if (addBgr ) {
786         FillModules(trH1,10000000); // Default mask 10M.
787         TTree *fAli=fLoader->GetRunLoader()->TreeK();
788         TFile *fileAli=0;
789         if (fAli) fileAli =fAli->GetCurrentFile();
790         fileAli->cd();
791     } // end if add
792 }
793 //______________________________________________________________________
794 void AliITS::FillModules(TTree *treeH, Int_t mask) {
795     // fill the modules with the sorted by module hits; 
796     // can be called many times to do a merging
797     // Inputs:
798     //      TTree *treeH  The tree containing the hits to be copied into
799     //                    the modules.
800     //      Int_t mask    The track number mask to indecate which file
801     //                    this hits came from.
802     // Outputs:
803     //      none.
804     // Return:
805     //      none.
806
807     if (treeH == 0x0)
808      {
809        Error("FillModules","Tree is NULL");
810      }
811     Int_t lay,lad,det,index;
812     AliITShit *itsHit=0;
813     AliITSmodule *mod=0;
814     char branchname[20];
815     sprintf(branchname,"%s",GetName());
816     TBranch *branch = treeH->GetBranch(branchname);
817     if (!branch) {
818         Error("FillModules","%s branch in TreeH not found",branchname);
819         return;
820     } // end if !branch
821     branch->SetAddress(&fHits);
822     Int_t nTracks =(Int_t) treeH->GetEntries();
823     Int_t iPrimTrack,h;
824     for(iPrimTrack=0; iPrimTrack<nTracks; iPrimTrack++){
825         ResetHits();
826         Int_t nBytes = treeH->GetEvent(iPrimTrack);
827         if (nBytes <= 0) continue;
828         Int_t nHits = fHits->GetEntriesFast();
829         for(h=0; h<nHits; h++){
830             itsHit = (AliITShit *)fHits->UncheckedAt(h);
831             itsHit->GetDetectorID(lay,lad,det);
832             if (fITSgeom) {
833                 index = fITSgeom->GetModuleIndex(lay,lad,det);
834             } else {
835                 index=det-1; // This should not be used.
836             } // end if [You must have fITSgeom for this to work!]
837             mod = GetModule(index);
838             itsHit->SetTrack(itsHit->GetTrack()+mask); // Set track mask.
839             mod->AddHit(itsHit,iPrimTrack,h);
840         } // end loop over hits 
841     } // end loop over tracks
842 }
843 //______________________________________________________________________
844 void AliITS::ClearModules(){
845     // Clear the modules TObjArray.
846     // Inputs:
847     //      none.
848     // Outputs:
849     //      none.
850
851     if(fITSmodules) fITSmodules->Delete();
852 }
853 //______________________________________________________________________
854 void AliITS::MakeBranchS(const char *fl){
855     // Creates Tree Branch for the ITS summable digits.
856     // Inputs:
857     //      cont char *fl  File name where SDigits branch is to be written
858     //                     to. If blank it write the SDigits to the same
859     //                     file in which the Hits were found.
860     // Outputs:
861     //      none.
862     // Return:
863     //      none.
864     Int_t buffersize = 4000;
865     char branchname[30];
866
867     // only one branch for SDigits.
868     sprintf(branchname,"%s",GetName());
869     
870
871     if(fLoader->TreeS()){
872         if (fSDigits == 0x0)  fSDigits  = new TClonesArray("AliITSpListItem",1000);
873         MakeBranchInTree(fLoader->TreeS(),branchname,&fSDigits,buffersize,fl);
874     } // end if
875 }
876 //______________________________________________________________________
877 void AliITS::SetTreeAddressS(TTree *treeS){
878     // Set branch address for the ITS summable digits Trees.
879     // Inputs:
880     //      TTree *treeS   Tree containing the SDigits.
881     // Outputs:
882     //      none.
883     // Return:
884     //      none.
885     char branchname[30];
886
887     if(!treeS) return;
888     if (fSDigits == 0x0)  fSDigits = new TClonesArray("AliITSpListItem",1000);
889     TBranch *branch;
890     sprintf(branchname,"%s",GetName());
891     branch = treeS->GetBranch(branchname);
892     if (branch) branch->SetAddress(&fSDigits);
893 }
894 //______________________________________________________________________
895 void AliITS::MakeBranchInTreeD(TTree *treeD,const char *file){
896     // Creates Tree branches for the ITS.
897     // Inputs:
898     //      TTree     *treeD Pointer to the Digits Tree.
899     //      cont char *file  File name where Digits branch is to be written
900     //                       to. If blank it write the SDigits to the same
901     //                       file in which the Hits were found.
902     // Outputs:
903     //      none.
904     // Return:
905     //      none.
906     Int_t buffersize = 4000;
907     char branchname[30];
908
909     sprintf(branchname,"%s",GetName());
910     // one branch for digits per type of detector
911     const char *det[3] = {"SPD","SDD","SSD"};
912     char digclass[40];
913     char clclass[40];
914     Int_t i;
915     for (i=0; i<kNTYPES ;i++) {
916         DetType(i)->GetClassNames(digclass,clclass);
917         // digits
918         if (fDtype == 0x0) fDtype = new TObjArray(fNDetTypes);
919         if(!(fDtype->At(i))) fDtype->AddAt(new TClonesArray(digclass,1000),i);
920         else ResetDigits(i);
921     } // end for i
922     for (i=0; i<kNTYPES ;i++) {
923         if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
924         else  sprintf(branchname,"%sDigits%d",GetName(),i+1);      
925         if (fDtype && treeD) {
926             MakeBranchInTree(treeD, branchname, &((*fDtype)[i]),buffersize,file);
927         } // end if
928     } // end for i
929 }
930 //______________________________________________________________________
931 void AliITS::SetTreeAddressD(TTree *treeD){
932     // Set branch address for the Trees.
933     // Inputs:
934     //      TTree *treeD   Tree containing the Digits.
935     // Outputs:
936     //      none.
937     // Return:
938     //      none.
939     char branchname[30];
940     const char *det[3] = {"SPD","SDD","SSD"};
941     TBranch *branch;
942     char digclass[40];
943     char clclass[40];
944     Int_t i;
945
946     if(!treeD) return;
947     for (i=0; i<kNTYPES; i++) {
948         DetType(i)->GetClassNames(digclass,clclass);
949         // digits
950         if (fDtype == 0x0) fDtype = new TObjArray(fNDetTypes);
951         if(!(fDtype->At(i))) fDtype->AddAt(new TClonesArray(digclass,1000),i);
952         else ResetDigits(i);
953         if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
954         else  sprintf(branchname,"%sDigits%d",GetName(),i+1);
955         if (fDtype) {
956             branch = treeD->GetBranch(branchname);
957             if (branch) branch->SetAddress(&((*fDtype)[i]));
958         } // end if fDtype
959     } // end for i
960 }
961 //______________________________________________________________________
962 void AliITS::Hits2SDigits(){
963     // Standard Hits to summable Digits function.
964     // Inputs:
965     //      none.
966     // Outputs:
967     //      none.
968
969 //    return; // Using Hits in place of the larger sDigits.
970     fLoader->LoadHits("read");
971     fLoader->LoadSDigits("recreate");
972     AliRunLoader* rl = fLoader->GetRunLoader(); 
973
974     for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
975     // Do the Hits to Digits operation. Use Standard input values.
976     // Event number from file, no background hit merging , use size from
977     // AliITSgeom class, option="All", input from this file only.
978       rl->GetEvent(iEvent);
979       if (!fLoader->TreeS()) fLoader->MakeTree("S");
980       MakeBranch("S");
981       SetTreeAddress();
982       HitsToSDigits(iEvent,0,-1," ",fOpt," ");
983     }
984     
985     fLoader->UnloadHits();
986     fLoader->UnloadSDigits();
987 }
988 //______________________________________________________________________
989 void AliITS::Hits2PreDigits(){
990     // Standard Hits to summable Digits function.
991     // Inputs:
992     //      none.
993     // Outputs:
994     //      none.
995
996     AliHeader *header=fLoader->GetRunLoader()->GetHeader(); // Get event number from this file.
997     // Do the Hits to Digits operation. Use Standard input values.
998     // Event number from file, no background hit merging , use size from
999     // AliITSgeom class, option="All", input from this file only.
1000     HitsToPreDigits(header->GetEvent(),0,-1," ",fOpt," ");
1001 }
1002 //______________________________________________________________________
1003 AliDigitizer* AliITS::CreateDigitizer(AliRunDigitizer* manager) const
1004 {
1005   return new AliITSDigitizer(manager);
1006 }
1007 //______________________________________________________________________
1008 void AliITS::SDigitsToDigits(Option_t *opt){
1009     // Standard Summable digits to Digits function.
1010     // Inputs:
1011     //      none.
1012     // Outputs:
1013     //      none.
1014     char name[20] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
1015
1016     if(!GetITSgeom()) return; // need transformations to do digitization.
1017     AliITSgeom *geom = GetITSgeom();
1018
1019     const char *all = strstr(opt,"All");
1020     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1021                           strstr(opt,"SSD")};
1022     if( !det[0] && !det[1] && !det[2] ) all = "All";
1023     else all = 0;
1024     static Bool_t setDef=kTRUE;
1025     if (setDef) SetDefaultSimulation();
1026     setDef=kFALSE;
1027
1028     AliITSsimulation *sim      = 0;
1029     AliITSDetType    *iDetType = 0;
1030     TTree            *trees    = fLoader->TreeS();
1031     if( !(trees && this->GetSDigits()) ){
1032         Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
1033         return;
1034     } // end if
1035     sprintf( name, "%s", this->GetName() );
1036     TBranch *brchSDigits = trees->GetBranch( name );
1037     
1038     Int_t id,module;
1039     for(module=0;module<geom->GetIndexMax();module++){
1040         id       = geom->GetModuleType(module);
1041         if (!all && !det[id]) continue;
1042         iDetType = DetType(id);
1043         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1044         if (!sim) {
1045             Error("SDigit2Digits","The simulation class was not "
1046                   "instanciated for module %d type %s!",module,
1047                   geom->GetModuleTypeName(module));
1048             exit(1);
1049         } // end if !sim
1050         sim->InitSimulationModule(module,gAlice->GetEvNumber());
1051 //
1052         // add summable digits to module
1053         this->GetSDigits()->Clear();
1054         brchSDigits->GetEvent(module);
1055         sim->AddSDigitsToModule(GetSDigits(),0);
1056 //
1057         // Digitise current module sum(SDigits)->Digits
1058         sim->FinishSDigitiseModule();
1059
1060         // fills all branches - wasted disk space
1061         fLoader->TreeD()->Fill();
1062         this->ResetDigits();
1063     } // end for module
1064
1065     fLoader->TreeD()->GetEntries();
1066
1067     fLoader->TreeD()->AutoSave();
1068     // reset tree
1069     fLoader->TreeD()->Reset();
1070     
1071 }
1072 //______________________________________________________________________
1073 void AliITS::Hits2Digits(){
1074     // Standard Hits to Digits function.
1075     // Inputs:
1076     //      none.
1077     // Outputs:
1078     //      none.
1079
1080     fLoader->LoadHits("read");
1081     fLoader->LoadDigits("recreate");
1082     AliRunLoader* rl = fLoader->GetRunLoader(); 
1083
1084     for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
1085     // Do the Hits to Digits operation. Use Standard input values.
1086     // Event number from file, no background hit merging , use size from
1087     // AliITSgeom class, option="All", input from this file only.
1088       rl->GetEvent(iEvent);
1089       if (!fLoader->TreeD()) fLoader->MakeTree("D");
1090       MakeBranch("D");
1091       SetTreeAddress();   
1092       HitsToDigits(iEvent,0,-1," ",fOpt," ");
1093     }
1094
1095     fLoader->UnloadHits();
1096     fLoader->UnloadDigits();
1097 }
1098 //______________________________________________________________________
1099 void AliITS::HitsToSDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1100                           Option_t *option, Option_t *opt, const char *filename){
1101     // keep galice.root for signal and name differently the file for 
1102     // background when add! otherwise the track info for signal will be lost !
1103     // the condition below will disappear when the geom class will be
1104     // initialized for all versions - for the moment it is only for v5 !
1105     // 7 is the SDD beam test version. Dummy routine. Hits are ITS's Summable
1106     // Digits.
1107     // Inputs:
1108     //      Int_t evnt       Event to be processed.
1109     //      Int_t bgrev      Background Hit tree number.
1110     //      Int_t nmodules   Not used.
1111     //      Option_t *option String indicating if merging hits or not. To
1112     //                       merge hits set equal to "Add". Otherwise no
1113     //                       background hits are considered.
1114     //      Test_t *filename File name containing the background hits..
1115     // Outputs:
1116     //      none.
1117     // Return:
1118     //      none.
1119 //    return; // using Hits instead of the larger sdigits.
1120
1121     HitsToPreDigits(evNumber,bgrev,size,option,opt,filename);
1122 }
1123 //______________________________________________________________________
1124 void AliITS::HitsToPreDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1125                           Option_t *option, Option_t *opt, const char *filename){
1126     //   Keep galice.root for signal and name differently the file for 
1127     // background when add! otherwise the track info for signal will be lost !
1128     // the condition below will disappear when the geom class will be
1129     // initialized for all versions - for the moment it is only for v5 !
1130     // 7 is the SDD beam test version.
1131     // Inputs:
1132     //      Int_t evnt       Event to be processed.
1133     //      Int_t bgrev      Background Hit tree number.
1134     //      Int_t nmodules   Not used.
1135     //      Option_t *option String indicating if merging hits or not. To
1136     //                       merge hits set equal to "Add". Otherwise no
1137     //                       background hits are considered.
1138     //      Test_t *filename File name containing the background hits..
1139     // Outputs:
1140     //      none.
1141     // Return:
1142     //      none.
1143
1144     if(!GetITSgeom()) return; // need transformations to do digitization.
1145     AliITSgeom *geom = GetITSgeom();
1146
1147     const char *all = strstr(opt,"All");
1148     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1149                           strstr(opt,"SSD")};
1150     static Bool_t setDef=kTRUE;
1151     if (setDef) SetDefaultSimulation();
1152     setDef=kFALSE;
1153
1154     Int_t nmodules;
1155     InitModules(size,nmodules);
1156     FillModules(evNumber,bgrev,nmodules,option,filename);
1157
1158     AliITSsimulation *sim      = 0;
1159     AliITSDetType    *iDetType = 0;
1160     AliITSmodule     *mod      = 0;
1161     Int_t id,module;
1162     for(module=0;module<geom->GetIndexMax();module++){
1163         id       = geom->GetModuleType(module);
1164         if (!all && !det[id]) continue;
1165         iDetType = DetType(id);
1166         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1167         if (!sim) {
1168             Error("HitsToSDigits","The simulation class was not "
1169                   "instanciated for module %d type %s!",module,
1170                   geom->GetModuleTypeName(module));
1171             exit(1);
1172         } // end if !sim
1173         mod      = (AliITSmodule *)fITSmodules->At(module);
1174         sim->SDigitiseModule(mod,module,evNumber);
1175         // fills all branches - wasted disk space
1176         fLoader->TreeS()->Fill(); 
1177         ResetSDigits();
1178     } // end for module
1179
1180     ClearModules();
1181
1182     fLoader->TreeS()->GetEntries();
1183     fLoader->TreeS()->AutoSave();
1184     fLoader->WriteSDigits("OVERWRITE");
1185     // reset tree
1186     fLoader->TreeS()->Reset();
1187 }
1188 //______________________________________________________________________
1189 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1190                           Option_t *option, Option_t *opt, const char *filename){
1191     //   Keep galice.root for signal and name differently the file for 
1192     // background when add! otherwise the track info for signal will be lost !
1193     // the condition below will disappear when the geom class will be
1194     // initialized for all versions - for the moment it is only for v5 !
1195     // 7 is the SDD beam test version.
1196     // Inputs:
1197     //      Int_t evnt       Event to be processed.
1198     //      Int_t bgrev      Background Hit tree number.
1199     //      Int_t nmodules   Not used.
1200     //      Option_t *option String indicating if merging hits or not. To
1201     //                       merge hits set equal to "Add". Otherwise no
1202     //                       background hits are considered.
1203     //      Test_t *filename File name containing the background hits..
1204     // Outputs:
1205     //      none.
1206     // Return:
1207     //      none.
1208
1209     if(!GetITSgeom()) return; // need transformations to do digitization.
1210     AliITSgeom *geom = GetITSgeom();
1211
1212     const char *all = strstr(opt,"All");
1213     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1214                           strstr(opt,"SSD")};
1215     static Bool_t setDef=kTRUE;
1216     if (setDef) SetDefaultSimulation();
1217     setDef=kFALSE;
1218
1219     Int_t nmodules;
1220     InitModules(size,nmodules);
1221     FillModules(evNumber,bgrev,nmodules,option,filename);
1222
1223     AliITSsimulation *sim      = 0;
1224     AliITSDetType    *iDetType = 0;
1225     AliITSmodule     *mod      = 0;
1226     Int_t id,module;
1227     for(module=0;module<geom->GetIndexMax();module++){
1228         id       = geom->GetModuleType(module);
1229         if (!all && !det[id]) continue;
1230         iDetType = DetType(id);
1231         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1232         if (!sim) {
1233             Error("HitsToDigits","The simulation class was not "
1234                   "instanciated for module %d type %s!",module,
1235                   geom->GetModuleTypeName(module));
1236             exit(1);
1237         } // end if !sim
1238         mod      = (AliITSmodule *)fITSmodules->At(module);
1239         sim->DigitiseModule(mod,module,evNumber);
1240         // fills all branches - wasted disk space
1241         fLoader->TreeD()->Fill(); 
1242         ResetDigits();
1243     } // end for module
1244
1245     ClearModules();
1246
1247     fLoader->TreeD()->GetEntries();
1248     fLoader->TreeD()->AutoSave();
1249     // reset tree
1250     fLoader->TreeD()->Reset();
1251 }
1252 //______________________________________________________________________
1253 void AliITS::ResetSDigits(){
1254     // Reset the Summable Digits array.
1255     // Inputs:
1256     //      none.
1257     // Outputs:
1258     //      none.
1259
1260     if (fSDigits) fSDigits->Clear();
1261     fNSDigits = 0;
1262 }
1263 //______________________________________________________________________
1264 void AliITS::ResetDigits(){
1265     // Reset number of digits and the digits array for the ITS detector.
1266     // Inputs:
1267     //      none.
1268     // Outputs:
1269     //      none.
1270
1271     if (!fDtype) return;
1272
1273     Int_t i;
1274     for (i=0;i<kNTYPES;i++ ) {
1275         if (fDtype->At(i))    ((TClonesArray*)fDtype->At(i))->Clear();
1276         if (fNdtype)  fNdtype[i]=0;
1277     } // end for i
1278 }
1279 //______________________________________________________________________
1280 void AliITS::ResetDigits(Int_t i){
1281     // Reset number of digits and the digits array for this branch.
1282     // Inputs:
1283     //      none.
1284     // Outputs:
1285     //      none.
1286
1287     if (fDtype->At(i))    ((TClonesArray*)fDtype->At(i))->Clear();
1288     if (fNdtype)  fNdtype[i]=0;
1289 }
1290 //______________________________________________________________________
1291 void AliITS::AddSumDigit(AliITSpListItem &sdig){
1292     // Adds the a module full of summable digits to the summable digits tree.
1293     // Inputs:
1294     //      AliITSpListItem &sdig   SDigit to be added to SDigits tree.
1295     // Outputs:
1296     //      none.
1297     // Return:
1298     //      none.
1299
1300     TClonesArray &lsdig = *fSDigits;
1301     new(lsdig[fNSDigits++]) AliITSpListItem(sdig);
1302 }
1303 //______________________________________________________________________
1304 void AliITS::AddRealDigit(Int_t id, Int_t *digits){
1305     //   Add a real digit - as coming from data.
1306     // Inputs:
1307     //      Int_t id        Detector type number.
1308     //      Int_t *digits   Integer array containing the digits info. See 
1309     //                      AliITSdigit.h
1310     // Outputs:
1311     //      none.
1312     // Return:
1313     //      none.
1314
1315     TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1316     new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
1317 }
1318 //______________________________________________________________________
1319 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d){
1320     //    Add a simulated digit.
1321     // Inputs:
1322     //      Int_t id        Detector type number.
1323     //      AliITSdigit *d  Digit to be added to the Digits Tree. See 
1324     //                      AliITSdigit.h
1325     // Outputs:
1326     //      none.
1327     // Return:
1328     //      none.
1329
1330     TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1331
1332     switch(id){
1333     case 0:
1334         new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
1335         break;
1336     case 1:
1337         new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
1338         break;
1339     case 2:
1340         new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
1341         break;
1342     } // end switch id
1343 }
1344 //______________________________________________________________________
1345 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,
1346                          Int_t *hits,Float_t *charges){
1347   //   Add a simulated digit to the list.
1348   // Inputs:
1349   //      Int_t id        Detector type number.
1350   //      Float_t phys    Physics indicator. See AliITSdigits.h
1351   //      Int_t *digits   Integer array containing the digits info. See 
1352   //                      AliITSdigit.h
1353   //      Int_t *tracks   Integer array [AliITSdigitS?D::GetNTracks()] 
1354   //                      containing the track numbers that contributed to
1355   //                      this digit.
1356   //      Int_t *hits     Integer array [AliITSdigitS?D::GetNTracks()]
1357   //                      containing the hit numbers, from AliITSmodule, that
1358   //                      contributed to this digit.
1359   //      Float_t *charge Floating point array of the signals contributed
1360   //                      to this digit by each track.
1361   // Outputs:
1362   //      none.
1363   // Return:
1364   //      none.
1365
1366   TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1367   AliITSresponseSDD *resp = 0;
1368   switch(id){
1369   case 0:
1370     new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
1371     break;
1372   case 1:
1373     resp = (AliITSresponseSDD*)DetType(1)->GetResponseModel();
1374     new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,
1375                                                hits,charges,resp);
1376     break;
1377   case 2:
1378     new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
1379     break;
1380   } // end switch id
1381 }
1382
1383 //______________________________________________________________________
1384 void AliITS::Digits2Raw()
1385 {
1386 // convert digits of the current event to raw data
1387
1388   fLoader->LoadDigits();
1389   TTree* digits = fLoader->TreeD();
1390   if (!digits) {
1391     Error("Digits2Raw", "no digits tree");
1392     return;
1393   }
1394   SetTreeAddressD(digits);
1395
1396   AliITSDDLRawData rawWriter;
1397   //Verbose level
1398   // 0: Silent
1399   // 1: cout messages
1400   // 2: txt files with digits 
1401   //BE CAREFUL, verbose level 2 MUST be used only for debugging and
1402   //it is highly suggested to use this mode only for debugging digits files
1403   //reasonably small, because otherwise the size of the txt files can reach
1404   //quickly several MB wasting time and disk space.
1405   rawWriter.SetVerbose(0);
1406     
1407   //SILICON PIXEL DETECTOR
1408   Info("Digits2Raw", "Formatting raw data for SPD");
1409   rawWriter.RawDataSPD(digits->GetBranch("ITSDigitsSPD"));
1410     
1411   //SILICON DRIFT DETECTOR
1412   Info("Digits2Raw", "Formatting raw data for SDD");
1413   rawWriter.RawDataSDD(digits->GetBranch("ITSDigitsSDD"));
1414     
1415   //SILICON STRIP DETECTOR
1416   Info("Digits2Raw", "Formatting raw data for SSD");
1417   rawWriter.RawDataSSD(digits->GetBranch("ITSDigitsSSD"));
1418
1419   fLoader->UnloadDigits();
1420 }
1421
1422 //______________________________________________________________________
1423 void AliITS::MakeTreeC(Option_t *option){
1424   //   Create a separate tree to store the clusters.
1425   // Inputs:
1426   //      Option_t *option  string which must contain "C" otherwise
1427   //                        no Cluster Tree is created.
1428   // Outputs:
1429   //      none.
1430   // Return:
1431   //      none.
1432
1433   AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;    
1434     
1435   if (pITSLoader == 0x0) {
1436     Error("MakeTreeC","fLoader == 0x0 option=%s",option);
1437     return;
1438   }
1439   if (pITSLoader->TreeC() == 0x0) pITSLoader->MakeTree("C");
1440   MakeBranchC();
1441 }
1442
1443 void AliITS::MakeBranchC()
1444 {
1445 //Makes barnches in treeC
1446   AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;    
1447   if (pITSLoader == 0x0) 
1448    {
1449     Error("MakeTreeC","fLoader == 0x0");
1450     return;
1451    }
1452   TTree * lTC = pITSLoader->TreeC();
1453   if (lTC == 0x0)
1454    {
1455      Error("MakeTreeC","Can not get TreeC from Loader");
1456      return;
1457    }
1458
1459   Int_t buffersize = 4000;
1460   char branchname[30];
1461   const char *det[3] = {"SPD","SDD","SSD"};
1462   char digclass[40];
1463   char clclass[40];
1464
1465     // one branch for Clusters per type of detector
1466   Int_t i;   
1467   for (i=0; i<kNTYPES ;i++) 
1468     {
1469         AliITSDetType *iDetType=DetType(i); 
1470         iDetType->GetClassNames(digclass,clclass);
1471         // clusters
1472         if (fCtype == 0x0) fCtype  = new TObjArray(fNDetTypes);
1473         if(!ClustersAddress(i))
1474          {
1475           fCtype->AddAt(new TClonesArray(clclass,1000),i);
1476          }
1477         if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
1478         else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
1479         if (fCtype  && lTC) 
1480          {
1481            if (lTC->GetBranch(branchname))
1482             {
1483               Warning("MakeBranchC","Branch %s alread exists in TreeC",branchname);
1484             }
1485            else
1486             {
1487               Info("MakeBranchC","Creating branch %s in TreeC",branchname);
1488               lTC->Branch(branchname,&((*fCtype)[i]), buffersize);
1489             }
1490          } // end if fCtype && lTC
1491   } // end for i
1492 }
1493
1494 //______________________________________________________________________
1495 void AliITS::GetTreeC(Int_t event){
1496   //    Get the clusters tree for this event and set the branch address.
1497   // Inputs:
1498   //      Int_t event    Event number for the cluster tree.
1499   // Outputs:
1500   //      none.
1501   // Return:
1502   //      none.
1503   char branchname[30];
1504   const char *det[3] = {"SPD","SDD","SSD"};
1505
1506   AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;
1507   TTree * lTC = pITSLoader->TreeC();
1508
1509   ResetClusters();
1510   if (lTC) {
1511     pITSLoader->CleanRawClusters();
1512   } // end if TreeC()
1513
1514
1515   TBranch *branch;
1516
1517   if (lTC) {
1518     Int_t i;
1519         char digclass[40];
1520         char clclass[40];
1521         for (i=0; i<kNTYPES; i++) {
1522       AliITSDetType *iDetType=DetType(i); 
1523       iDetType->GetClassNames(digclass,clclass);
1524       // clusters
1525       if (fCtype == 0x0) fCtype  = new TObjArray(fNDetTypes);
1526       if(!fCtype->At(i)) fCtype->AddAt(new TClonesArray(clclass,1000),i);
1527       if(kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
1528       else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
1529       if (fCtype) {
1530                 branch = lTC->GetBranch(branchname);
1531         if (branch) branch->SetAddress(&((*fCtype)[i]));
1532       } // end if fCtype
1533         } // end for i
1534   } else {
1535         Error("GetTreeC","cannot find Clusters Tree for event:%d",event);
1536   } // end if lTC
1537 }
1538 //______________________________________________________________________
1539 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c){
1540     //   Add a cluster to the list.
1541     // Inputs:
1542     //      Int_t id             Detector type number.
1543     //      AliITSRawCluster *c  Cluster class to be added to the tree of
1544     //                           clusters.
1545     // Outputs:
1546     //      none.
1547     // Return:
1548     //      none.
1549
1550     TClonesArray &lc = *((TClonesArray*)fCtype->At(id));
1551
1552     switch(id){
1553     case 0:
1554         new(lc[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
1555         break;
1556     case 1:
1557         new(lc[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
1558         break;
1559     case 2:
1560         new(lc[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
1561         break;
1562     } // end switch id
1563 }
1564 //______________________________________________________________________
1565 void AliITS::ResetClusters(){
1566     // Reset number of clusters and the clusters array for ITS.
1567     // Inputs:
1568     //      none.
1569     // Outputs:
1570     //      none.
1571
1572     Int_t i;
1573     for (i=0;i<kNTYPES;i++ ) ResetClusters(i);
1574 }
1575 //______________________________________________________________________
1576 void AliITS::ResetClusters(Int_t i){
1577     //    Reset number of clusters and the clusters array for this branch.
1578     // Inputs:
1579     //      Int_t i        Detector type number.
1580     // Outputs:
1581     //      none.
1582     // Return:
1583     //      none.
1584
1585     if (fCtype->At(i))    ((TClonesArray*)fCtype->At(i))->Clear();
1586     if (fNctype)  fNctype[i]=0;
1587 }
1588 //______________________________________________________________________
1589 void AliITS::MakeBranchR(const char *file, Option_t *opt){
1590     // Creates Tree branches for the ITS Reconstructed points.
1591     // Inputs:
1592     //      cont char *file  File name where RecPoints branch is to be written
1593     //                       to. If blank it write the SDigits to the same
1594     //                       file in which the Hits were found.
1595     // Outputs:
1596     //      none.
1597     // Return:
1598     //      none.
1599     Int_t buffsz = 4000;
1600     char branchname[30];
1601
1602     // only one branch for rec points for all detector types
1603     Bool_t oFast= (strstr(opt,"Fast")!=0);
1604     if(oFast){
1605       sprintf(branchname,"%sRecPointsF",GetName());
1606     } else {
1607       sprintf(branchname,"%sRecPoints",GetName());
1608     }
1609     
1610     
1611     if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1612     if (fLoader->TreeR()) {
1613         if (fRecPoints == 0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1614         MakeBranchInTree(fLoader->TreeR(),branchname,&fRecPoints,buffsz,file);
1615     } // end if
1616 }
1617 //______________________________________________________________________
1618 void AliITS::SetTreeAddressR(TTree *treeR){
1619     // Set branch address for the Reconstructed points Trees.
1620     // Inputs:
1621     //      TTree *treeR   Tree containing the RecPoints.
1622     // Outputs:
1623     //      none.
1624     // Return:
1625     //      none.
1626     char branchname[30];
1627
1628     if(!treeR) return;
1629     if (fRecPoints == 0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1630     TBranch *branch;
1631     sprintf(branchname,"%sRecPoints",GetName());
1632     branch = treeR->GetBranch(branchname);
1633     if (branch) {
1634       branch->SetAddress(&fRecPoints);
1635     }
1636     else {
1637       sprintf(branchname,"%sRecPointsF",GetName());
1638       branch = treeR->GetBranch(branchname);
1639       if (branch) {
1640         branch->SetAddress(&fRecPoints);
1641       }
1642     }
1643 }
1644 //______________________________________________________________________
1645 void AliITS::AddRecPoint(const AliITSRecPoint &r){
1646     // Add a reconstructed space point to the list
1647     // Inputs:
1648     //      const AliITSRecPoint &r RecPoint class to be added to the tree
1649     //                              of reconstructed points TreeR.
1650     // Outputs:
1651     //      none.
1652     // Return:
1653     //      none.
1654
1655     TClonesArray &lrecp = *fRecPoints;
1656     new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
1657 }
1658 //______________________________________________________________________
1659 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1660                                   Option_t *opt0,Option_t *opt1, const char *flnm){
1661     // keep galice.root for signal and name differently the file for 
1662     // background when add! otherwise the track info for signal will be lost !
1663     // the condition below will disappear when the geom class will be
1664     // initialized for all versions - for the moment it is only for v5 !
1665     // Inputs:
1666     //      Int_t evnt       Event to be processed.
1667     //      Int_t bgrev      Background Hit tree number.
1668     //      Int_t size       Size used by InitModules. See InitModules.
1669     //      Option_t *opt0   Option passed to FillModules. See FillModules.
1670     //      Option_t *opt1   String indicating if merging hits or not. To
1671     //                       merge hits set equal to "Add". Otherwise no
1672     //                       background hits are considered.
1673     //      Test_t *flnm     File name containing the background hits..
1674     // Outputs:
1675     //      none.
1676     // Return:
1677     //      none.
1678
1679     if(!GetITSgeom()) return;
1680     AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
1681     AliITSgeom *geom = GetITSgeom();
1682
1683     const char *all = strstr(opt1,"All");
1684     const char *det[3] ={strstr(opt1,"SPD"),strstr(opt1,"SDD"),
1685                          strstr(opt1,"SSD")};
1686     Int_t nmodules;
1687     InitModules(size,nmodules);
1688     FillModules(evNumber,bgrev,nmodules,opt0,flnm);
1689
1690     AliITSsimulation *sim      = 0;
1691     AliITSDetType    *iDetType = 0;
1692     AliITSmodule     *mod      = 0;
1693     Int_t id,module;
1694
1695     //m.b. : this change is nothing but a nice way to make sure
1696     //the CPU goes up !
1697     
1698     cout<<"HitsToFastRecPoints: N mod = "<<geom->GetIndexMax()<<endl;
1699     
1700     for(module=0;module<geom->GetIndexMax();module++){
1701         id       = geom->GetModuleType(module);
1702         if (!all && !det[id]) continue;
1703         iDetType = DetType(id);
1704         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1705         if (!sim) 
1706          {
1707            Error("HitsToFastPoints","The simulation class was not "
1708                  "instanciated for module %d type %x!",module,
1709                  geom->GetModuleTypeName(module));
1710            exit(1);
1711          } // end if !sim
1712         mod      = (AliITSmodule *)fITSmodules->At(module);
1713         sim->CreateFastRecPoints(mod,module,gRandom);
1714         cout<<module<<"\r";fflush(0);
1715         //gAlice->TreeR()->Fill();
1716         TTree *lTR = pITSloader->TreeR();
1717         TBranch *br=lTR->GetBranch("ITSRecPointsF");
1718         br->Fill();
1719         ResetRecPoints();
1720     } // end for module
1721
1722     ClearModules();
1723     
1724     fLoader->WriteRecPoints("OVERWRITE");
1725 }
1726 //______________________________________________________________________
1727 void AliITS::Digits2Reco(){
1728     // Find clusters and reconstruct space points.
1729     // Inputs:
1730     //      none.
1731     // Outputs:
1732     //      none.
1733
1734     AliHeader *header=fLoader->GetRunLoader()->GetHeader();
1735     // to Digits to RecPoints for event in file, all digits in file, and
1736     // all ITS detectors.
1737     DigitsToRecPoints(header->GetEvent(),0,fOpt);
1738 }
1739 //______________________________________________________________________
1740 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt){
1741   // cluster finding and reconstruction of space points
1742   // the condition below will disappear when the geom class will be
1743   // initialized for all versions - for the moment it is only for v5 !
1744   // 7 is the SDD beam test version
1745   // Inputs:
1746   //      Int_t evNumber   Event number to be processed.
1747   //      Int_t lastentry  Offset for module when not all of the modules
1748   //                       are processed.
1749   //      Option_t *opt    String indicating which ITS sub-detectors should
1750   //                       be processed. If ="All" then all of the ITS
1751   //                       sub detectors are processed.
1752   // Outputs:
1753   //      none.
1754   // Return:
1755   //      none.
1756
1757   if(!GetITSgeom()) return;
1758   AliITSgeom *geom = GetITSgeom();
1759     
1760   const char *all = strstr(opt,"All");
1761   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1762                         strstr(opt,"SSD")};
1763   static Bool_t setRec=kTRUE;
1764   if (setRec) SetDefaultClusterFinders();
1765   setRec=kFALSE;
1766
1767   AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
1768   TTree *treeC=pITSloader->TreeC();
1769   AliITSClusterFinder *rec     = 0;
1770   AliITSDetType      *iDetType = 0;
1771   Int_t id,module,first=0;
1772   for(module=0;module<geom->GetIndexMax();module++){
1773       id       = geom->GetModuleType(module);
1774       if (!all && !det[id]) continue;
1775       if(det[id]) first = geom->GetStartDet(id);
1776       iDetType = DetType(id);
1777       rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1778       TClonesArray *itsDigits  = this->DigitsAddress(id);
1779       if (!rec) {
1780           Error("DigitsToRecPoints",
1781                 "The reconstruction class was not instanciated! event=%d",
1782                 evNumber);
1783           exit(1);
1784       } // end if !rec
1785       this->ResetDigits();
1786       TTree *lTD = pITSloader->TreeD();
1787       if (all) {
1788           lTD->GetEvent(lastentry+module);
1789       }else {
1790           lTD->GetEvent(lastentry+(module-first));
1791       }
1792       Int_t ndigits = itsDigits->GetEntriesFast();
1793       if (ndigits) rec->FindRawClusters(module);
1794       pITSloader->TreeR()->Fill(); 
1795       ResetRecPoints();
1796       treeC->Fill();
1797       ResetClusters();
1798   } // end for module
1799
1800
1801   pITSloader->WriteRecPoints("OVERWRITE");
1802   pITSloader->WriteRawClusters("OVERWRITE");
1803 }
1804 //______________________________________________________________________
1805 void AliITS::ResetRecPoints(){
1806     // Reset number of rec points and the rec points array.
1807     // Inputs:
1808     //      none.
1809     // Outputs:
1810     //      none.
1811
1812     if (fRecPoints) fRecPoints->Clear();
1813     fNRecPoints = 0;
1814 }
1815 //______________________________________________________________________
1816 AliLoader* AliITS::MakeLoader(const char* topfoldername)
1817
1818   //builds ITSgetter (AliLoader type)
1819   //if detector wants to use castomized getter, it must overload this method
1820
1821   Info("MakeLoader","Creating AliITSLoader. Top folder is %s.",topfoldername);
1822   fLoader = new AliITSLoader(GetName(),topfoldername);
1823   return fLoader;
1824 }
1825