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