]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITS.cxx
Patches and updates for fixes to this and other routines.
[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 /*
17 $Log$
18 Revision 1.18  2000/07/12 05:32:20  fca
19 Correcting several syntax problem with static members
20
21 Revision 1.17  2000/07/10 16:07:18  fca
22 Release version of ITS code
23
24 Revision 1.9.2.3  2000/02/02 13:42:09  barbera
25 fixed AliITS.cxx for new AliRun structure. Added ITS hits list to list of hits which will have their track numbers updated
26
27 Revision 1.9.2.2  2000/01/23 03:03:13  nilsen
28 //fixed FillModule. Removed fi(fabs(xl)<dx....
29
30 Revision 1.9.2.1  2000/01/12 19:03:32  nilsen
31 This is the version of the files after the merging done in December 1999.
32 See the ReadMe110100.txt file for details
33
34 Revision 1.9  1999/11/14 14:33:25  fca
35 Correct problems with distructors and pointers, thanks to I.Hrivnacova
36
37 Revision 1.8  1999/09/29 09:24:19  fca
38 Introduction of the Copyright and cvs Log
39
40 */
41
42 ///////////////////////////////////////////////////////////////////////////////
43 //
44 //      An overview of the basic philosophy of the ITS code development
45 // and analysis is show in the figure below.
46 //Begin_Html
47 /*
48 <img src="picts/ITS/ITS_Analysis_schema.gif">
49 </pre>
50 <br clear=left>
51 <font size=+2 color=red>
52 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
53 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
54 </font>
55 <pre>
56 */
57 //End_Html
58 //
59 //  AliITS. Inner Traking System base class.
60 //  This class contains the base procedures for the Inner Tracking System
61 //
62 //Begin_Html
63 /*
64 <img src="picts/ITS/AliITS_Class_Diagram.gif">
65 </pre>
66 <br clear=left>
67 <font size=+2 color=red>
68 <p>This show the class diagram of the different elements that are part of
69 the AliITS class.
70 </font>
71 <pre>
72 */
73 //End_Html
74 //
75 // Version: 0
76 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
77 //
78 // Version: 1
79 // Modified and documented by Bjorn S. Nilsen
80 // July 11 1999
81 //
82 // Version: 2
83 // Modified and documented by A. Bologna
84 // October 18 1999
85 //
86 // AliITS is the general base class for the ITS. Also see AliDetector for
87 // futher information.
88 //
89 ///////////////////////////////////////////////////////////////////////////////
90  
91 #include <TMath.h>
92 #include <TRandom.h>
93 #include <TBranch.h>
94 #include <TVector.h>
95 #include <TObjArray.h>
96 #include <TROOT.h>
97 #include <TObjectTable.h>
98 #include <TFile.h>
99 #include <TTree.h>
100
101
102
103 #include "AliRun.h"
104 #include "AliITS.h"
105 #include "AliITSMap.h"
106 #include "AliITSDetType.h"
107 #include "AliITSClusterFinder.h"
108 #include "AliITSsimulation.h"
109 #include "AliITSsegmentationSPD.h"
110 #include "AliITSresponseSPD.h"
111 #include "AliITSsegmentationSDD.h"
112 #include "AliITSresponseSDD.h"
113 #include "AliITSsegmentationSSD.h"
114 #include "AliITSresponseSSD.h"
115 #include "AliITShit.h"
116 #include "AliITSgeom.h"
117 #include "AliITSdigit.h"
118 #include "AliITSmodule.h"
119 #include "AliITSRecPoint.h"
120 #include "AliITSRawCluster.h"
121
122 const Int_t AliITS::fgkNTYPES=3;
123
124 ClassImp(AliITS)
125  
126 //_____________________________________________________________________________
127 AliITS::AliITS() : AliDetector() {
128   //
129   // Default initialiser for ITS
130   //     The default constructor of the AliITS class. In addition to
131   // creating the AliITS class it zeros the variables fIshunt (a member
132   // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
133   // fITSpoints, fIdSens, and fIdName. The AliDetector default constructor
134   // is also called.
135   //
136
137
138   fIshunt     = 0;
139   fEuclidOut  = 0;
140
141   //fNDetTypes = fgkNTYPES;
142   fIdN        = 0;
143   fIdName     = 0;
144   fIdSens     = 0;
145   fITSmodules = 0;
146   //
147   fDetTypes   = 0;
148   //
149   fDtype  = 0;
150   fNdtype = 0;
151   fCtype  = 0;
152   fNctype = 0;
153   fRecPoints = 0;
154   fNRecPoints = 0;
155   fTreeC = 0;
156   //
157   fITSgeom=0;
158 }
159
160 //_____________________________________________________________________________
161 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title){
162   //
163   // Default initialiser for ITS
164   //     The constructor of the AliITS class. In addition to creating the
165   // AliITS class, it allocates memory for the TClonesArrays fHits and
166   // fDigits, and for the TObjArray fITSpoints. It also zeros the variables
167   // fIshunt (a member of AliDetector class), fEuclidOut, and fIdN, and zeros
168   // the pointers fIdSens and fIdName. To help in displaying hits via the ROOT
169   // macro display.C AliITS also sets the marker color to red. The variables
170   // passes with this constructor, const char *name and *title, are used by
171   // the constructor of AliDetector class. See AliDetector class for a
172   // description of these parameters and its constructor functions.
173   //
174
175
176   fHits       = new TClonesArray("AliITShit", 1560);
177   gAlice->AddHitList(fHits);
178
179   //fNDetTypes = fgkNTYPES;
180
181   fNdtype = new Int_t[fgkNTYPES];
182   fDtype = new TObjArray(fgkNTYPES);
183
184   fNctype = new Int_t[fgkNTYPES];
185   fCtype = new TObjArray(fgkNTYPES);
186
187
188   fRecPoints = 0;
189   fNRecPoints = 0;
190
191   fTreeC = 0;
192
193   fITSmodules = 0; 
194
195   fIshunt     = 0;
196   fEuclidOut  = 0;
197   fIdN        = 0;
198   fIdName     = 0;
199   fIdSens     = 0;
200  
201   fDetTypes = new TObjArray(fgkNTYPES);  
202
203   Int_t i;
204   for(i=0;i<fgkNTYPES;i++) {
205     (*fDetTypes)[i]=new AliITSDetType(); 
206     fNdtype[i]=0;
207     fNctype[i]=0;
208    }
209   //
210
211   SetMarkerColor(kRed);
212
213   fITSgeom=0;
214 }
215 //___________________________________________________________________________
216 AliITS::AliITS(AliITS &source){
217   // copy constructor
218   if(this==&source) return;
219   printf("Error: You are not allowed to make a copy of the AliITS\n");
220   exit(1);
221 }
222 //____________________________________________________________________________
223 AliITS& AliITS::operator=(AliITS &source){
224   // assignment operator
225   if(this==&source) return *this;
226   printf("Error: You are not allowed to make a copy of the AliITS\n");
227   exit(1);
228   return *this; //fake return
229 }
230 //____________________________________________________________________________
231 void AliITS::ClearModules(){
232   //clear the modules TObjArray
233
234   if(fITSmodules) fITSmodules->Delete();
235
236 }
237 //_____________________________________________________________________________
238 AliITS::~AliITS(){
239   //
240   // Default distructor for ITS
241   //     The default destructor of the AliITS class. In addition to deleting
242   // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
243   // fIdSens, fIdName, and fITSpoints.
244   //
245
246
247   delete fHits;
248   delete fDigits;
249   delete fRecPoints;
250   if(fIdName!=0) delete[] fIdName;
251   if(fIdSens!=0) delete[] fIdSens;
252   if(fITSmodules!=0) {
253       this->ClearModules();
254       delete fITSmodules;
255   }// end if fITSmodules!=0
256
257   //
258   if(fDtype) {
259     fDtype->Delete();
260     delete fDtype;
261   }
262   delete [] fNdtype;
263   if (fCtype) {
264     fCtype->Delete();
265     delete fCtype;
266   }
267   delete [] fNctype;
268   //
269
270   if (fDetTypes) {
271     fDetTypes->Delete();
272     delete fDetTypes;
273   }
274
275   if (fTreeC) delete fTreeC;
276
277   if (fITSgeom) delete fITSgeom;
278
279 }
280
281 //___________________________________________
282 AliITSDetType* AliITS::DetType(Int_t id)
283 {
284   //return pointer to id detector type
285     return ((AliITSDetType*) (*fDetTypes)[id]);
286
287 }
288 //___________________________________________
289 void AliITS::SetClasses(Int_t id, const char *digit, const char *cluster)
290 {
291   //set the digit and cluster classes to be used for the id detector type
292     ((AliITSDetType*) (*fDetTypes)[id])->ClassNames(digit,cluster);
293
294 }
295 //___________________________________________
296 void AliITS::SetResponseModel(Int_t id, AliITSresponse *response)
297 {
298   //set the response model for the id detector type
299
300     ((AliITSDetType*) (*fDetTypes)[id])->ResponseModel(response);
301
302 }
303
304 //___________________________________________
305 void AliITS::SetSegmentationModel(Int_t id, AliITSsegmentation *seg)
306 {
307   //set the segmentation model for the id detector type
308
309     ((AliITSDetType*) (*fDetTypes)[id])->SegmentationModel(seg);
310
311 }
312
313 //___________________________________________
314 void AliITS::SetSimulationModel(Int_t id, AliITSsimulation *sim)
315 {
316   //set the simulation model for the id detector type
317
318    ((AliITSDetType*) (*fDetTypes)[id])->SimulationModel(sim);
319
320 }
321 //___________________________________________
322 void AliITS::SetReconstructionModel(Int_t id, AliITSClusterFinder *reconst)
323 {
324   //set the cluster finder model for the id detector type
325
326    ((AliITSDetType*) (*fDetTypes)[id])->ReconstructionModel(reconst);
327
328 }
329
330 //_____________________________________________________________________________
331 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
332   //
333   // Add an ITS hit
334   //     The function to add information to the AliITShit class. See the
335   // AliITShit class for a full description. This function allocates the
336   // necessary new space for the hit information and passes the variable
337   // track, and the pointers *vol and *hits to the AliITShit constructor
338   // function.
339   //
340   TClonesArray &lhits = *fHits;
341   new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
342 }
343 //_____________________________________________________________________________
344 void AliITS::AddRealDigit(Int_t id, Int_t *digits) 
345 {
346   // add a real digit - as coming from data
347
348   TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
349   new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
350
351 }
352 //_____________________________________________________________________________
353 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d) 
354 {
355
356   // add a simulated digit
357
358   TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
359
360   switch(id)
361   {
362   case 0:
363      new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
364      break;
365   case 1:
366      new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
367      break;
368   case 2:
369      new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
370      break;
371   }
372
373 }
374
375 //_____________________________________________________________________________
376 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,Int_t *hits,Float_t *charges){
377
378   // add a simulated digit to the list
379
380   TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
381   switch(id)
382   {
383   case 0:
384      new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
385      break;
386   case 1:
387      new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,hits,charges);
388      break;
389   case 2:
390      new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
391      break;
392   }
393  
394 }
395
396 //_____________________________________________________________________________
397 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c) 
398 {
399
400   // add a cluster to the list
401
402   TClonesArray &lcl = *((TClonesArray*)(*fCtype)[id]);
403
404   switch(id)
405   {
406   case 0:
407      new(lcl[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
408      break;
409   case 1:
410      new(lcl[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
411      break;
412   case 2:
413      new(lcl[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
414      break;
415   }
416
417 }
418
419
420 //_____________________________________________________________________________
421 void AliITS::AddRecPoint(const AliITSRecPoint &r)
422 {
423   //
424   // Add a reconstructed space point to the list
425   //
426   TClonesArray &lrecp = *fRecPoints;
427   new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
428 }
429  
430
431 //____________________________________________
432 void AliITS::ResetDigits()
433 {
434     //
435     // Reset number of digits and the digits array for the ITS detector
436     //
437
438     if (!fDtype) return;
439
440     Int_t i;
441     for (i=0;i<fgkNTYPES;i++ ) {
442         if ((*fDtype)[i])    ((TClonesArray*)(*fDtype)[i])->Clear();
443         if (fNdtype)  fNdtype[i]=0;
444     }
445 }
446
447 //____________________________________________
448 void AliITS::ResetDigits(Int_t i)
449 {
450     //
451     // Reset number of digits and the digits array for this branch
452     //
453   if ((*fDtype)[i])    ((TClonesArray*)(*fDtype)[i])->Clear();
454   if (fNdtype)  fNdtype[i]=0;
455 }
456
457
458 //____________________________________________
459 void AliITS::ResetClusters()
460 {
461     //
462     // Reset number of clusters and the clusters array for ITS
463     //
464
465     Int_t i;
466     for (i=0;i<fgkNTYPES;i++ ) {
467         if ((*fCtype)[i])    ((TClonesArray*)(*fCtype)[i])->Clear();
468         if (fNctype)  fNctype[i]=0;
469     }
470
471 }
472
473 //____________________________________________
474 void AliITS::ResetClusters(Int_t i)
475 {
476     //
477     // Reset number of clusters and the clusters array for this branch
478     //
479         if ((*fCtype)[i])    ((TClonesArray*)(*fCtype)[i])->Clear();
480         if (fNctype)  fNctype[i]=0;
481
482 }
483
484
485 //____________________________________________
486 void AliITS::ResetRecPoints()
487 {
488     //
489     // Reset number of rec points and the rec points array 
490     //
491     if (fRecPoints) fRecPoints->Clear();
492     fNRecPoints = 0;
493
494 }
495
496 //_____________________________________________________________________________
497 Int_t AliITS::DistancetoPrimitive(Int_t , Int_t ){
498   //
499   // Distance from mouse to ITS on the screen. Dummy routine
500   //     A dummy routine used by the ROOT macro display.C to allow for the
501   // use of the mouse (pointing device) in the macro. In general this should
502   // never be called. If it is it returns the number 9999 for any value of
503   // x and y.
504   //
505   return 9999;
506 }
507
508 //_____________________________________________________________________________
509 void AliITS::Init(){
510   //
511   // Initialise ITS after it has been built
512   //     This routine initializes the AliITS class. It is intended to be called
513   // from the Init function in AliITSv?. Besides displaying a banner
514   // indicating that it has been called it initializes the array fIdSens
515   // and sets the default segmentation, response, digit and raw cluster classes
516   // Therefore it should be called after a call to CreateGeometry.
517   //
518
519
520   SetDefaults();
521
522   Int_t i;
523   printf("\n");
524   for(i=0;i<35;i++) printf("*");
525   printf(" ITS_INIT ");
526   for(i=0;i<35;i++) printf("*");
527   printf("\n");
528   //
529   //
530   for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
531   //
532   for(i=0;i<80;i++) printf("*");
533   printf("\n");
534 }
535
536 //_____________________________________________________________________________
537 void AliITS::SetDefaults()
538 {
539   // sets the default segmentation, response, digit and raw cluster classes
540
541   printf("SetDefaults\n");
542
543   AliITSDetType *iDetType;
544
545   //SPD 
546
547   AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
548   AliITSresponseSPD *resp0=new AliITSresponseSPD();
549   iDetType=DetType(0); 
550   if (!iDetType->GetSegmentationModel()) SetSegmentationModel(0,seg0); 
551   if (!iDetType->GetResponseModel()) SetResponseModel(0,resp0); 
552   // set digit and raw cluster classes to be used
553   const char *kData0=resp0->DataType();
554   if (strstr(kData0,"real")) {
555       iDetType->ClassNames("AliITSdigit","AliITSRawClusterSPD");
556   } else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
557
558   // SDD                                          //
559   AliITSresponseSDD *resp1=new AliITSresponseSDD();
560   AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
561   iDetType=DetType(1); 
562   printf("SetDefaults: iDetType %p\n",iDetType);
563   if (!iDetType->GetSegmentationModel()) SetSegmentationModel(1,seg1); 
564   printf("SetDefaults: segm %p\n",iDetType->GetSegmentationModel());
565   if (!iDetType->GetResponseModel()) SetResponseModel(1,resp1); 
566   printf("SetDefaults: resp %p\n",iDetType->GetResponseModel());
567   const char *kData1=resp1->DataType();
568   const char *kopt=resp1->ZeroSuppOption();
569   if ((!strstr(kopt,"2D")) && (!strstr(kopt,"1D")) || strstr(kData1,"real") ) {
570       iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
571   } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
572
573   // SSD
574   AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
575   AliITSresponseSSD *resp2=new AliITSresponseSSD();
576   iDetType=DetType(2); 
577   if (!iDetType->GetSegmentationModel()) SetSegmentationModel(2,seg2); 
578   if (!iDetType->GetResponseModel()) SetResponseModel(2,resp2); 
579   const char *kData2=resp2->DataType();
580   if (strstr(kData2,"real")) {
581       iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
582   } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
583
584   if (fgkNTYPES>3) {
585     Warning("SetDefaults","Only the three basic detector types are initialised!");
586   } 
587
588 }
589
590
591 //_____________________________________________________________________________
592 void AliITS::SetDefaultSimulation()
593 {
594   // to be written
595
596 }
597 //_____________________________________________________________________________
598 void AliITS::SetDefaultClusterFinders()
599 {
600   // to be written
601
602 }
603 //_____________________________________________________________________________
604
605 void AliITS::MakeTreeC(Option_t *option)
606 {
607   // create a separate tree to store the clusters
608
609   printf("MakeTreeC \n");
610
611      char *optC = strstr(option,"C");
612      if (optC && !fTreeC) fTreeC = new TTree("TC","Clusters in ITS");
613      else return;
614
615      Int_t buffersize = 4000;
616      char branchname[30];
617
618      char *det[3] = {"SPD","SDD","SSD"};
619
620      // one branch for Clusters per type of detector
621      Int_t i;
622      for (i=0; i<fgkNTYPES ;i++) {
623         if (fgkNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
624         else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
625         if (fCtype   && fTreeC) {
626            TreeC()->Branch(branchname,&((*fCtype)[i]), buffersize);
627            printf("Making Branch %s for Clusters of detector type %d\n",branchname,i+1);
628         }       
629      }
630
631 }
632
633 //_____________________________________________________________________________
634 void AliITS::GetTreeC(Int_t event)
635 {
636
637   printf("GetTreeC \n");
638
639   // get the clusters tree for this event and set the branch address
640     char treeName[20];
641     char branchname[30];
642
643     char *det[3] = {"SPD","SDD","SSD"};
644
645     ResetClusters();
646     if (fTreeC) {
647           delete fTreeC;
648     }
649
650     sprintf(treeName,"TreeC%d",event);
651     fTreeC = (TTree*)gDirectory->Get(treeName);
652
653
654     TBranch *branch;
655     if (fTreeC) {
656         Int_t i;
657         for (i=0; i<fgkNTYPES; i++) {
658            if (fgkNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
659            else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
660            if (fCtype) {
661                 branch = fTreeC->GetBranch(branchname);
662                 if (branch) branch->SetAddress(&((*fCtype)[i]));
663             }
664         }
665     } else {
666         printf("ERROR: cannot find Clusters Tree for event:%d\n",event);
667     }
668
669 }
670 //_____________________________________________________________________________
671 void AliITS::MakeBranch(Option_t* option){
672   //
673   // Creates Tree branches for the ITS.
674   //
675
676
677   Int_t buffersize = 4000;
678   char branchname[30];
679   sprintf(branchname,"%s",GetName());
680
681   AliDetector::MakeBranch(option);
682
683
684 // one branch for digits per type of detector
685   
686    char *det[3] = {"SPD","SDD","SSD"};
687
688    char digclass[40];
689    char clclass[40];
690
691    Int_t i;
692    for (i=0; i<fgkNTYPES ;i++) {
693        AliITSDetType *iDetType=DetType(i); 
694        iDetType->GetClassNames(digclass,clclass);
695        //printf("i, digclass, recclass %d %s %s\n",i,digclass,clclass); 
696        // digits
697        (*fDtype)[i] = new TClonesArray(digclass,10000); 
698        // clusters
699        (*fCtype)[i] = new TClonesArray(clclass,10000); 
700    }
701
702
703   for (i=0; i<fgkNTYPES ;i++) {
704       if (fgkNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
705       else  sprintf(branchname,"%sDigits%d",GetName(),i+1);
706       
707       if (fDtype   && gAlice->TreeD()) {
708           gAlice->TreeD()->Branch(branchname,&((*fDtype)[i]), buffersize);
709           printf("Making Branch %s for digits of type %d\n",branchname,i+1);
710       } 
711   }
712
713   // only one branch for rec points for all detector types
714   sprintf(branchname,"%sRecPoints",GetName());
715
716   fRecPoints=new TClonesArray("AliITSRecPoint",10000);
717
718   if (fRecPoints && gAlice->TreeR()) {
719     gAlice->TreeR()->Branch(branchname,&fRecPoints, buffersize);
720     printf("Making Branch %s for reconstructed space points\n",branchname);
721   }     
722
723
724 }
725
726 //___________________________________________
727 void AliITS::SetTreeAddress()
728 {
729
730   // Set branch address for the Trees.
731
732   char branchname[30];
733   AliDetector::SetTreeAddress();
734
735   char *det[3] = {"SPD","SDD","SSD"};
736
737   TBranch *branch;
738   TTree *treeD = gAlice->TreeD();
739   TTree *treeR = gAlice->TreeR();
740
741   Int_t i;
742   if (treeD) {
743       for (i=0; i<fgkNTYPES; i++) {
744           if (fgkNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
745           else  sprintf(branchname,"%sDigits%d",GetName(),i+1);
746           if (fDtype) {
747               branch = treeD->GetBranch(branchname);
748               if (branch) branch->SetAddress(&((*fDtype)[i]));
749           }
750       }
751   }
752
753  
754   if (treeR) {
755     sprintf(branchname,"%sRecPoints",GetName());
756     if (fRecPoints) {
757       branch = treeR->GetBranch(branchname);
758       if (branch) branch->SetAddress(&fRecPoints);
759     }
760   }
761   
762
763 }
764
765 //____________________________________________________________________________
766 void AliITS::InitModules(Int_t size,Int_t &nmodules){
767
768   //initialize the modules array
769
770   if(fITSmodules){ 
771       fITSmodules->Delete();
772       delete fITSmodules;
773   }
774
775     Int_t nl,indexMAX,index;
776
777     if(size<=0){ // default to using data stored in AliITSgeom
778         if(fITSgeom==0) {
779             printf("Error in AliITS::InitModule fITSgeom not defined\n");
780             return;
781         } // end if fITSgeom==0
782         nl = fITSgeom->GetNlayers();
783         indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
784                                             fITSgeom->GetNdetectors(nl))+1;
785         nmodules = indexMAX;
786         fITSmodules = new TObjArray(indexMAX);
787         for(index=0;index<indexMAX;index++){
788                 fITSmodules->AddAt( new AliITSmodule(index),index);
789         } // end for index
790     }else{
791         fITSmodules = new TObjArray(size);
792         for(index=0;index<size;index++) {
793             fITSmodules->AddAt( new AliITSmodule(index),index);
794         }
795
796         nmodules = size;
797     } // end i size<=0
798 }
799
800 //____________________________________________________________________________
801 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,Option_t *option,Text_t *filename){
802
803   // fill the modules with the sorted by module hits; add hits from background
804   // if option=Add
805
806
807     static TTree *trH1;                 //Tree with background hits
808     static TClonesArray *fHits2;        //List of hits for one track only
809
810     static Bool_t first=kTRUE;
811     static TFile *file;
812     char *addBgr = strstr(option,"Add");
813
814
815     if (addBgr ) {
816         if(first) {
817             cout<<"filename "<<filename<<endl;
818             file=new TFile(filename);
819             cout<<"I have opened "<<filename<<" file "<<endl;
820             fHits2     = new TClonesArray("AliITShit",1000  );
821         }           
822         first=kFALSE;
823         file->cd();
824         file->ls();
825         // Get Hits Tree header from file
826         if(fHits2) fHits2->Clear();
827         if(trH1) delete trH1;
828         trH1=0;
829         
830         char treeName[20];
831         sprintf(treeName,"TreeH%d",bgrev);
832         trH1 = (TTree*)gDirectory->Get(treeName);
833         //printf("TrH1 %p of treename %s for event %d \n",trH1,treeName,bgrev);
834         
835         if (!trH1) {
836             printf("ERROR: cannot find Hits Tree for event:%d\n",bgrev);
837         }
838         // Set branch addresses
839         TBranch *branch;
840         char branchname[20];
841         sprintf(branchname,"%s",GetName());
842         if (trH1 && fHits2) {
843             branch = trH1->GetBranch(branchname);
844             if (branch) branch->SetAddress(&fHits2);
845         }
846
847         // test
848         //Int_t ntracks1 =(Int_t)TrH1->GetEntries();
849         //printf("background - ntracks1 - %d\n",ntracks1);
850    }
851
852     Int_t npart = gAlice->GetEvent(evnt);
853     if(npart<=0) return;
854     TClonesArray *itsHits = this->Hits();
855     Int_t lay,lad,det,index;
856     AliITShit *itsHit=0;
857     AliITSmodule *mod=0;
858
859     TTree *iTH = gAlice->TreeH();
860     Int_t ntracks =(Int_t) iTH->GetEntries();
861
862     Int_t t,h;
863     for(t=0; t<ntracks; t++){
864         gAlice->ResetHits();
865         iTH->GetEvent(t);
866         Int_t nhits = itsHits->GetEntriesFast();
867         //printf("nhits %d\n",nhits);
868         if (!nhits) continue;
869         for(h=0; h<nhits; h++){
870             itsHit = (AliITShit *)itsHits->UncheckedAt(h);
871             itsHit->GetDetectorID(lay,lad,det);
872             // temporarily index=det-1 !!!
873             if(fITSgeom) index = fITSgeom->GetModuleIndex(lay,lad,det);
874             else index=det-1;
875             //
876             mod = this->GetModule(index);
877             mod->AddHit(itsHit,t,h);
878         } // end loop over hits 
879     } // end loop over tracks
880
881     // open the file with background
882     
883     if (addBgr ) {
884           Int_t track,i;
885           ntracks =(Int_t)trH1->GetEntries();
886             //printf("background - ntracks1 %d\n",ntracks);
887             //printf("background - Start loop over tracks \n");     
888             //   Loop over tracks
889
890             for (track=0; track<ntracks; track++) {
891
892                 if (fHits2)       fHits2->Clear();
893                 trH1->GetEvent(track);
894                 //   Loop over hits
895                 for(i=0;i<fHits2->GetEntriesFast();++i) {
896
897                     itsHit=(AliITShit*) (*fHits2)[i];
898                     itsHit->GetDetectorID(lay,lad,det);
899                     // temporarily index=det-1 !!!
900                     if(fITSgeom) index = fITSgeom->GetModuleIndex(lay,lad,det);
901                     else index=det-1;
902                     //
903                     mod = this->GetModule(index);
904                     mod->AddHit(itsHit,track,i);
905                }  // end loop over hits
906             } // end loop over tracks
907
908             TTree *fAli=gAlice->TreeK();
909             TFile *fileAli=0;
910             
911             if (fAli) fileAli =fAli->GetCurrentFile();
912             fileAli->cd();
913
914     } // end if add
915
916     //gObjectTable->Print();
917
918 }
919
920
921 //____________________________________________________________________________
922 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size, Option_t *option, Option_t *opt,Text_t *filename)
923 {
924     // keep galice.root for signal and name differently the file for 
925     // background when add! otherwise the track info for signal will be lost !
926   
927    // the condition below will disappear when the geom class will be
928    // initialised for all versions - for the moment it is only for v5 !
929    // 7 is the SDD beam test version  
930    Int_t ver = this->IsVersion(); 
931    if(ver!=5 && ver!=7) return; 
932
933    char *all = strstr(opt,"All");
934    char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
935
936    Int_t nmodules;
937    InitModules(size,nmodules); 
938    FillModules(evNumber,bgrev,nmodules,option,filename);
939
940    //TBranch *branch;
941    AliITSsimulation* sim;
942    //TObjArray *branches=gAlice->TreeD()->GetListOfBranches();
943    AliITSgeom *geom = GetITSgeom();
944
945    Int_t id,module;
946    Int_t first,last;
947    for (id=0;id<fgkNTYPES;id++) {
948         if (!all && !det[id]) continue;
949         //branch = (TBranch*)branches->UncheckedAt(id);
950         AliITSDetType *iDetType=DetType(id); 
951         sim = (AliITSsimulation*)iDetType->GetSimulationModel();
952         if (!sim) {
953            Error("HitsToDigits","The simulation class was not instantiated!");
954            exit(1);
955            // or SetDefaultSimulation();
956         }
957         if(geom) {
958           first = geom->GetStartDet(id);
959           last = geom->GetLastDet(id);
960         } else first=last=0;
961         printf("det type %d first, last %d %d \n",id,first,last);
962         for(module=first;module<=last;module++) {
963             AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
964             sim->DigitiseModule(mod,module,evNumber);
965             // fills all branches - wasted disk space
966             gAlice->TreeD()->Fill(); 
967             ResetDigits();
968             // try and fill only the branch 
969             //branch->Fill();
970             //ResetDigits(id);
971         } // loop over modules
972    } // loop over detector types
973
974    ClearModules();
975
976    Int_t nentries=(Int_t)gAlice->TreeD()->GetEntries();
977    printf("nentries in TreeD %d\n",nentries);
978
979    char hname[30];
980    sprintf(hname,"TreeD%d",evNumber);
981    gAlice->TreeD()->Write(hname);
982    // reset tree
983    gAlice->TreeD()->Reset();
984
985 }
986
987
988 //____________________________________________________________________________
989 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt)
990 {
991   // cluster finding and reconstruction of space points
992   
993    // the condition below will disappear when the geom class will be
994    // initialised for all versions - for the moment it is only for v5 !
995    // 7 is the SDD beam test version  
996    Int_t ver = this->IsVersion(); 
997    if(ver!=5 && ver!=7) return; 
998
999    char *all = strstr(opt,"All");
1000    char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1001
1002    static Bool_t first=kTRUE;
1003    if (first) {
1004        MakeTreeC("C");
1005        first=kFALSE;
1006    }
1007  
1008    TTree *iTC=TreeC();
1009
1010    //TBranch *branch;
1011    AliITSClusterFinder* rec;
1012
1013    //TObjArray *branches=gAlice->TreeR()->GetListOfBranches();
1014    AliITSgeom *geom = GetITSgeom();
1015
1016    Int_t id,module;
1017    for (id=0;id<fgkNTYPES;id++) {
1018         if (!all && !det[id]) continue;
1019         //branch = (TBranch*)branches->UncheckedAt(id);
1020         AliITSDetType *iDetType=DetType(id); 
1021         rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1022         if (!rec) {
1023            Error("DigitsToRecPoints","The cluster finder class was not instantiated!");
1024            exit(1);
1025            // or SetDefaultClusterFinders();
1026         }
1027         TClonesArray *itsDigits  = this->DigitsAddress(id);
1028
1029         Int_t first,last;
1030         if(geom) {
1031           first = geom->GetStartDet(id);
1032           last = geom->GetLastDet(id);
1033         } else first=last=0;
1034         //printf("first last %d %d\n",first,last);
1035         for(module=first;module<=last;module++) {
1036               this->ResetDigits();
1037               if (all) gAlice->TreeD()->GetEvent(lastentry+module);
1038               else gAlice->TreeD()->GetEvent(lastentry+(module-first));
1039               Int_t ndigits = itsDigits->GetEntriesFast();
1040               if (ndigits) rec->FindRawClusters();
1041               gAlice->TreeR()->Fill(); 
1042               ResetRecPoints();
1043               iTC->Fill();
1044               ResetClusters();
1045               // try and fill only the branch 
1046               //branch->Fill();
1047               //ResetRecPoints(id);
1048         } // loop over modules
1049    } // loop over detector types
1050
1051
1052    Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1053    Int_t ncentries=(Int_t)iTC->GetEntries();
1054    printf(" nentries ncentries %d %d\n", nentries, ncentries);
1055
1056    char hname[30];
1057    sprintf(hname,"TreeR%d",evNumber);
1058    gAlice->TreeR()->Write(hname);
1059    // reset tree
1060    gAlice->TreeR()->Reset();
1061
1062    sprintf(hname,"TreeC%d",evNumber);
1063    iTC->Write(hname);
1064    iTC->Reset();
1065 }
1066
1067
1068 //____________________________________________________________________________
1069 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1070 Option_t *option,Option_t *opt,Text_t *filename)
1071 {
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   
1075
1076    // the condition below will disappear when the geom class will be
1077    // initialised for all versions - for the moment it is only for v5 !  
1078    Int_t ver = this->IsVersion(); 
1079    if(ver!=5) return; 
1080
1081    char *all = strstr(opt,"All");
1082    char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1083
1084    Int_t nmodules;
1085    InitModules(size,nmodules);
1086    FillModules(evNumber,bgrev,nmodules,option,filename);
1087
1088
1089    AliITSsimulation* sim;
1090    AliITSgeom *geom = GetITSgeom();
1091
1092    TRandom *random=new TRandom[9];
1093    random[0].SetSeed(111);
1094    random[1].SetSeed(222);
1095    random[2].SetSeed(333);              
1096    random[3].SetSeed(444);
1097    random[4].SetSeed(555);
1098    random[5].SetSeed(666);              
1099    random[6].SetSeed(777);
1100    random[7].SetSeed(888);
1101    random[8].SetSeed(999);              
1102
1103
1104    Int_t id,module;
1105    for (id=0;id<fgkNTYPES;id++) {
1106         if (!all && !det[id]) continue;
1107         AliITSDetType *iDetType=DetType(id); 
1108         sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1109         if (!sim) {
1110            Error("HitsToFastPoints","The simulation class was not instantiated!");
1111            exit(1);
1112            // or SetDefaultSimulation();
1113         }
1114         Int_t first = geom->GetStartDet(id);
1115         Int_t last = geom->GetLastDet(id);
1116         for(module=first;module<=last;module++) {
1117             AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1118             sim->CreateFastRecPoints(mod,module,random);
1119             gAlice->TreeR()->Fill(); 
1120             ResetRecPoints();
1121         } // loop over modules
1122    } // loop over detector types
1123
1124
1125    ClearModules();
1126
1127    //Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1128
1129    char hname[30];
1130    sprintf(hname,"TreeR%d",evNumber);
1131    gAlice->TreeR()->Write(hname);
1132    // reset tree
1133    gAlice->TreeR()->Reset();
1134
1135    delete [] random;
1136
1137 }
1138
1139 //____________________________________________________________________________
1140 void AliITS::Streamer(TBuffer &R__b){
1141    // Stream an object of class AliITS.
1142
1143    Int_t i;
1144
1145    if (R__b.IsReading()) {
1146       Version_t R__v = R__b.ReadVersion();
1147       if (R__v) {
1148           AliDetector::Streamer(R__b);
1149           R__b >> fIdN;
1150           R__b.ReadArray(fIdSens); 
1151           for(i=0;i<fIdN;i++) fIdName[i].Streamer(R__b);
1152           R__b >> fITSgeom;
1153           R__b >> fITSmodules;
1154           R__b >> fEuclidOut;
1155           R__b >> fMajorVersion;
1156           R__b >> fMinorVersion;
1157           R__b >> fDetTypes;
1158           R__b >> fDtype;
1159           delete []fNdtype; 
1160           fNdtype = new Int_t[fgkNTYPES];
1161           R__b.ReadFastArray(fNdtype,fgkNTYPES);
1162           R__b >> fCtype;
1163           delete []fNctype; 
1164           fNctype = new Int_t[fgkNTYPES];
1165           R__b.ReadFastArray(fNctype,fgkNTYPES);
1166           R__b >> fRecPoints;
1167           R__b >> fNRecPoints;
1168           R__b >> fTreeC;
1169       } // end if R__v
1170    } else { // writing
1171       R__b.WriteVersion(AliITS::IsA());
1172       AliDetector::Streamer(R__b);
1173       R__b << fIdN;
1174       R__b.WriteArray(fIdSens,fIdN); 
1175       for(i=0;i<fIdN;i++) fIdName[i].Streamer(R__b);
1176       R__b << fITSgeom;
1177       R__b << fITSmodules;
1178       R__b << fEuclidOut;
1179       R__b << fMajorVersion;
1180       R__b << fMinorVersion;
1181       R__b << fDetTypes;
1182       R__b << fDtype;
1183       R__b.WriteFastArray(fNdtype,fgkNTYPES);
1184       R__b << fCtype;
1185       R__b.WriteFastArray(fNctype,fgkNTYPES);
1186       R__b << fRecPoints;
1187       R__b << fNRecPoints;
1188       R__b << fTreeC;
1189    } // end if
1190
1191 }
1192
1193