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