]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITS.cxx
Splitting the variables used to accumulate the sums from the ones used to fill the...
[u/mrichter/AliRoot.git] / ITS / AliITS.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18
19 ///////////////////////////////////////////////////////////////////////////////
20 //                                                                           //
21 //      An overview of the basic philosophy of the ITS code development      //
22 // and analysis is show in the figure below.                                 //
23 //Begin_Html                                                                 //
24 /*                                               
25 <img src="picts/ITS/ITS_Analysis_schema.gif">
26 </pre>
27 <br clear=left>
28 <font size=+2 color=red>
29 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
30 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
31 </font>
32 <pre>
33 */
34 //End_Html
35 //
36 //  AliITS. Inner Traking System base class.
37 //  This class contains the base procedures for the Inner Tracking System
38 //
39 //Begin_Html
40 /*
41 <img src="picts/ITS/AliITS_Class_Diagram.gif">
42 </pre>
43 <br clear=left>
44 <font size=+2 color=red>
45 <p>This show the class diagram of the different elements that are part of
46 the AliITS class.
47 </font>
48 <pre>
49 */
50 //End_Html
51 //
52 // Version: 0
53 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
54 //
55 // Version: 1
56 // Modified and documented by Bjorn S. Nilsen
57 // July 11 1999
58 //
59 // Version: 2
60 // Modified and documented by A. Bologna
61 // October 18 1999
62 //
63 // AliITS is the general base class for the ITS. Also see AliDetector for
64 // futher information.
65 //
66 ///////////////////////////////////////////////////////////////////////////////
67
68 #include <stdlib.h>
69 #include <TClonesArray.h>
70 #include <TFile.h>
71 #include <TParticle.h>
72 #include <TString.h>
73 #include <TTree.h>
74 #include <TVirtualMC.h>
75 #include "AliDetector.h"
76 #include "AliITS.h"
77 #include "AliITSDetTypeSim.h"
78 #include "AliITSDDLRawData.h"
79 #include "AliITSLoader.h"
80 #include "AliITShit.h"
81 #include "AliITSmodule.h"
82 #include "AliITSpListItem.h"
83 #include "AliITSsimulation.h"
84 #include "AliITSsimulationFastPoints.h"
85 #include "AliMC.h"
86 #include "AliITSDigitizer.h"
87 #include "AliITSRecPoint.h"
88 #include "AliITSsegmentationSPD.h"
89 #include "AliITSsegmentationSDD.h"
90 #include "AliITSsegmentationSSD.h"
91 #include "AliITSRawStreamSPD.h"
92 #include "AliITSRawStreamSSD.h"
93 #include "AliITSRawStreamSDD.h"
94 #include "AliRawReader.h"
95 #include "AliRun.h"
96 #include "AliLog.h"
97 #include "AliITSInitGeometry.h"
98
99 ClassImp(AliITS)
100
101 //______________________________________________________________________
102 AliITS::AliITS() : AliDetector(),
103 fDetTypeSim(0),
104 fEuclidOut(0),
105 fOpt("All"),
106 fIdN(0),
107 fIdSens(0),
108 fIdName(0),
109 fITSmodules(0),
110 fTiming(kFALSE)
111 {
112   // Default initializer for ITS
113   //      The default constructor of the AliITS class. In addition to
114   // creating the AliITS class it zeros the variables fIshunt (a member
115   // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
116   // fIdSens, and fIdName. The AliDetector default constructor
117   // is also called.
118   
119 //    SetDetectors(); // default to fOpt="All". This variable not written out.
120 //PH    SetMarkerColor(kRed);
121 }
122 //______________________________________________________________________
123 AliITS::AliITS(const Char_t *title):AliDetector("ITS",title),
124 fDetTypeSim(0),
125 fEuclidOut(0),
126 fOpt("All"),
127 fIdN(0),
128 fIdSens(0),
129 fIdName(0),
130 fITSmodules(0),
131 fTiming(kFALSE)
132 {
133     //     The standard Constructor for the ITS class. 
134     // It also zeros the variables
135     // fIshunt (a member of AliDetector class), fEuclidOut, and zeros
136     // the pointers fIdSens and fIdName. To help in displaying hits via the
137     // ROOT macro display.C AliITS also sets the marker color to red. The
138     // variables passes with this constructor, const char *name and *title,
139     // are used by the constructor of AliDetector class. See AliDetector
140     // class for a description of these parameters and its constructor
141     // functions.
142     // Inputs:
143     //   Char_t *title  Simulation title for the ITS
144     // Outputs:
145     //   none.
146     // Return:
147     //   none.
148   
149     fHits = new TClonesArray("AliITShit",1560); // from AliDetector
150     if(gAlice->GetMCApp()) gAlice->GetMCApp()->AddHitList(fHits);
151     //fNhits=0;  //done in AliDetector(name,title)
152     SetDetectors(); // default to fOpt="All". This variable not written out.
153     fDetTypeSim   = new AliITSDetTypeSim();
154     //PH  SetMarkerColor(kRed);
155     if(!fLoader) MakeLoader(AliConfig::GetDefaultEventFolderName());
156     fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
157 }
158 //______________________________________________________________________
159 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title),
160 fDetTypeSim(0),
161 fEuclidOut(0),
162 fOpt("All"),
163 fIdN(0),
164 fIdSens(0),
165 fIdName(0),
166 fITSmodules(0),
167 fTiming(kFALSE)
168 {
169   //     The standard Constructor for the ITS class. 
170   // It also zeros the variables
171   // fIshunt (a member of AliDetector class), fEuclidOut, and zeros
172   // the pointers fIdSens and fIdName. To help in displaying hits via the
173   // ROOT macro display.C AliITS also sets the marker color to red. The
174   // variables passes with this constructor, const char *name and *title,
175   // are used by the constructor of AliDetector class. See AliDetector
176   // class for a description of these parameters and its constructor
177   // functions.
178   
179   fHits = new TClonesArray("AliITShit",1560);
180   if(gAlice->GetMCApp()) gAlice->GetMCApp()->AddHitList(fHits);
181   //fNhits=0;  //done in AliDetector(name,title)
182
183   SetDetectors(); // default to fOpt="All". This variable not written out.
184     
185   fDetTypeSim   = new AliITSDetTypeSim();
186   
187   //PH  SetMarkerColor(kRed);
188   if(!fLoader) MakeLoader(AliConfig::GetDefaultEventFolderName());
189   fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
190
191 }
192 //______________________________________________________________________
193 AliITS::~AliITS(){
194     // Default destructor for ITS.
195     //     The default destructor of the AliITS class. In addition to deleting
196     // the AliITS class it deletes the memory pointed to by 
197     // fIdSens, fIdName, fDetTypeSim and it's contents.
198     // Inputs:
199     //      none.
200     // Outputs:
201     //      none.
202     // Return:
203     //      none.
204
205     if (fHits) {
206       fHits->Delete();
207       delete fHits;
208       fHits=0;
209     }
210     if(fITSmodules) {
211         this->ClearModules();
212         delete fITSmodules;
213         fITSmodules = 0;
214     }// end if fITSmodules!=0
215
216     delete[] fIdName;  // Array of TStrings
217     delete[] fIdSens;
218
219     if (fDetTypeSim){
220       delete fDetTypeSim;
221       fDetTypeSim = 0;
222     }
223 }
224 //______________________________________________________________________
225 AliDigitizer* AliITS::CreateDigitizer(AliRunDigitizer* manager)const{
226     // Creates the AliITSDigitizer in a standard way for use via AliModule.
227     // This function can not be included in the .h file because of problems
228     // with the order of inclusion (recursive).
229     // Inputs:
230     //    AliRunDigitizer *manager  The Manger class for Digitization
231     // Output:
232     //    none.
233     // Return:
234     //    A new AliITSRunDigitizer (cast as a AliDigitizer).
235
236      return new AliITSDigitizer(manager);
237 }
238 //______________________________________________________________________
239 void AliITS::Init(){
240     // Initializer ITS after it has been built
241     //     This routine initializes the AliITS class. It is intended to be
242     // called from the Init function in AliITSv?. Besides displaying a banner
243     // indicating that it has been called it initializes the array fIdSens
244     // and sets the default segmentation, response, digit and raw cluster
245     // classes therefore it should be called after a call to CreateGeometry.
246     // Inputs:
247     //      none.
248     // Outputs:
249     //      none.
250     // Return:
251     //      none.
252     Int_t i;
253     // Array of TStrings
254     if(gMC) for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
255  
256 }
257
258
259 //______________________________________________________________________
260 void AliITS::SetDefaults(){
261     // sets the default segmentation, response, digit and raw cluster classes.
262     // Inputs:
263     //      none.
264     // Outputs:
265     //      none.
266     // Return:
267     //      none.
268     AliInfoClass("AliITS::Setting Defaults");
269     if(!fDetTypeSim) { 
270      Error("SetDefaults()","fDetTypeSim is 0!"); 
271      return;
272     }
273
274     fDetTypeSim->SetDefaults();
275
276
277 }
278 //______________________________________________________________________
279 void AliITS::SetDefaultSimulation(){
280     // sets the default simulation.
281     // Inputs:
282     //      none.
283     // Outputs:
284     //      none.
285     // Return:
286     //      none.
287     if(!fDetTypeSim) { 
288      Error("SetDefaultSimulation()","fDetTypeSim is 0!"); 
289      return;
290     }
291
292     fDetTypeSim->SetDefaultSimulation();
293
294 }
295
296
297 //______________________________________________________________________
298 void AliITS::MakeBranch(Option_t* option){
299     // Creates Tree branches for the ITS.
300     // Inputs:
301     //      Option_t *option    String of Tree types S,D, and/or R.
302     //      const char *file    String of the file name where these branches
303     //                          are to be stored. If blank then these branches
304     //                          are written to the same tree as the Hits were
305     //                          read from.
306     // Outputs:
307     //      none.
308     // Return:
309     //      none.
310   if(!fDetTypeSim) {
311     Error("MakeBranch","fDetTypeSim is 0!");
312     return;
313   }
314
315   Bool_t cH = (strstr(option,"H")!=0);
316   Bool_t cS = (strstr(option,"S")!=0);
317   Bool_t cD = (strstr(option,"D")!=0);
318   
319   if(cH && (fHits == 0x0)) fHits = new TClonesArray("AliITShit", 1560);
320   AliDetector::MakeBranch(option);
321   
322   if(cS) MakeBranchS(0);
323   if(cD) MakeBranchD(0);
324
325
326 }
327 //___________________________________________________________________
328 void AliITS::MakeBranchS(const char* fl){
329
330   // Creates Tree Branch for the ITS summable digits.
331   // Inputs:
332   //      cont char *fl  File name where SDigits branch is to be written
333   //                     to. If blank it write the SDigits to the same
334   //                     file in which the Hits were found.
335
336   
337   if(!fDetTypeSim){
338     Error("MakeBranchS","fDetTypeSim is 0!");
339   }
340   Int_t buffersize = 4000;
341   char branchname[30];
342
343   // only one branch for SDigits.
344   sprintf(branchname,"%s",GetName());
345
346   if(fLoader->TreeS()){
347     TClonesArray* sdig = (TClonesArray*)fDetTypeSim->GetSDigits();
348     MakeBranchInTree(fLoader->TreeS(),branchname,&sdig,buffersize,fl);
349   } 
350 }
351 //______________________________________________________________________
352 void AliITS::MakeBranchD(const char* file){
353
354   //Make branch for digits
355   if(!fDetTypeSim) {
356     Warning("MakeBranchD","fDetTypeSim is 0!");
357     return;
358   }
359   fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
360   MakeBranchInTreeD(fLoader->TreeD(),file);
361 }
362
363 //___________________________________________________________________
364 void AliITS:: MakeBranchInTreeD(TTree* treeD, const char* file){
365   // Creates Tree branches for the ITS.
366
367   if(!fDetTypeSim){
368     Error("MakeBranchS","fDetTypeSim is 0!");
369   }
370   fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
371
372   const Char_t *det[3] = {"SPD","SDD","SSD"};
373   Char_t* digclass;
374   Int_t buffersize = 4000;
375   Char_t branchname[30];
376   
377   if(!fDetTypeSim->GetDigits()){
378     fDetTypeSim->SetDigits(new TObjArray(fgkNTYPES));
379   }
380   for(Int_t i=0;i<fgkNTYPES;i++){
381     digclass = fDetTypeSim->GetDigitClassName(i);
382     TString classn = digclass;
383     if(!((fDetTypeSim->GetDigits())->At(i))){
384       (fDetTypeSim->GetDigits())->AddAt(new TClonesArray(classn.Data(),1000),i);
385     }
386     else ResetDigits(i);  
387     if(fgkNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
388     else sprintf(branchname,"%sDigits%d",GetName(),i+1);
389     TObjArray* dig = DigitsAddress(i);
390     if(GetDigits() && treeD) AliDetector::MakeBranchInTree(treeD,branchname, &dig,buffersize,file);
391   }
392
393 }
394 //______________________________________________________________________
395 void AliITS::SetTreeAddress(){
396     // Set branch address for the Trees.
397     // Inputs:
398     //      none.
399     // Outputs:
400     //      none.
401     // Return:
402     //      none.
403     
404   if(!fDetTypeSim) {
405     Error("SetTreeAddress","fDetTypeSim is 0!");
406     return;
407   }
408
409   fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
410
411   TTree *treeS = fLoader->TreeS();
412   TTree *treeD = fLoader->TreeD();
413   if (fLoader->TreeH() && (fHits == 0x0)) {
414       fHits = new TClonesArray("AliITShit", 1560);
415   }
416   AliDetector::SetTreeAddress();
417
418   fDetTypeSim->SetTreeAddressS(treeS, (Char_t*)GetName());
419   fDetTypeSim->SetTreeAddressD(treeD, (Char_t*)GetName());
420 }
421 //______________________________________________________________________
422 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
423     // Add an ITS hit
424     //     The function to add information to the AliITShit class. See the
425     // AliITShit class for a full description. This function allocates the
426     // necessary new space for the hit information and passes the variable
427     // track, and the pointers *vol and *hits to the AliITShit constructor
428     // function.
429     // Inputs:
430     //      Int_t   track   Track number which produced this hit.
431     //      Int_t   *vol    Array of Integer Hit information. See AliITShit.h
432     //      Float_t *hits   Array of Floating Hit information.  see AliITShit.h
433     // Outputs:
434     //      none.
435     // Return:
436     //      none.
437   TClonesArray &lhits = *fHits;
438   new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
439 }
440
441 //______________________________________________________________________
442 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,
443                          Option_t *option, const char *filename){
444   // fill the modules with the sorted by module hits; add hits from
445   // background if option=Add.
446
447   static TTree *trH1;                 //Tree with background hits
448   static Bool_t first=kTRUE;
449   static TFile *file;
450   const char *addBgr = strstr(option,"Add");
451   
452   evnt = nmodules; // Dummy use of variables to remove warnings
453   if (addBgr ) {
454     if(first) {
455       file=new TFile(filename);
456     } // end if first
457     first=kFALSE;
458     file->cd();
459     file->ls();
460     // Get Hits Tree header from file
461     if(trH1) delete trH1;
462     trH1=0;
463     
464     char treeName[20];
465     sprintf(treeName,"TreeH%d",bgrev);
466     trH1 = (TTree*)gDirectory->Get(treeName);
467     if (!trH1) {
468       Error("FillModules","cannot find Hits Tree for event:%d",bgrev);
469     } // end if !trH1
470     // Set branch addresses
471   } // end if addBgr
472   
473   FillModules(fLoader->TreeH(),0); // fill from this file's tree.
474     
475   if (addBgr ) {
476     FillModules(trH1,10000000); // Default mask 10M.
477     TTree *fAli=fLoader->GetRunLoader()->TreeK();
478     TFile *fileAli=0;
479     if (fAli) fileAli =fAli->GetCurrentFile();
480     fileAli->cd();
481   } // end if add
482   
483   
484 }
485
486
487 //______________________________________________________________________
488 void AliITS::FillModules(TTree *treeH, Int_t mask) {
489     // fill the modules with the sorted by module hits; 
490     // can be called many times to do a merging
491     // Inputs:
492     //      TTree *treeH  The tree containing the hits to be copied into
493     //                    the modules.
494     //      Int_t mask    The track number mask to indecate which file
495     //                    this hits came from.
496     // Outputs:
497     //      none.
498     // Return:
499     //      none.
500
501     if (treeH == 0x0)
502      {
503        Error("FillModules","Tree is NULL");
504      }
505     Int_t lay,lad,det,index;
506     AliITShit *itsHit=0;
507     AliITSmodule *mod=0;
508     char branchname[20];
509     sprintf(branchname,"%s",GetName());
510     TBranch *branch = treeH->GetBranch(branchname);
511     if (!branch) {
512         Error("FillModules","%s branch in TreeH not found",branchname);
513         return;
514     } // end if !branch
515     branch->SetAddress(&fHits);
516     Int_t nTracks =(Int_t) treeH->GetEntries();
517     Int_t iPrimTrack,h;
518     for(iPrimTrack=0; iPrimTrack<nTracks; iPrimTrack++){
519         ResetHits();
520         Int_t nBytes = treeH->GetEvent(iPrimTrack);
521         if (nBytes <= 0) continue;
522         Int_t nHits = fHits->GetEntriesFast();
523         for(h=0; h<nHits; h++){
524             itsHit = (AliITShit *)fHits->UncheckedAt(h);
525             itsHit->GetDetectorID(lay,lad,det);
526             if (GetITSgeom()) {
527                 index = GetITSgeom()->GetModuleIndex(lay,lad,det);
528             } else {
529                 index=det-1; // This should not be used.
530             } // end if [You must have fITSgeom for this to work!]
531             mod = GetModule(index);
532             itsHit->SetTrack(itsHit->GetTrack()+mask); // Set track mask.
533             mod->AddHit(itsHit,iPrimTrack,h);
534         } // end loop over hits 
535     } // end loop over tracks
536 }
537
538 //______________________________________________________________________
539 Bool_t AliITS::InitModules(Int_t size,Int_t &nmodules){
540     // Initialize the modules array.
541     // Inputs:
542     //      Int_t size  Size of array of the number of modules to be
543     //                  created. If size <=0 then the number of modules
544     //                  is gotten from AliITSgeom class kept in fITSgeom.
545     // Outputs:
546     //      Int_t &nmodules The number of modules existing.
547     // Return:
548     //      none.
549
550     if(fITSmodules){ 
551         fITSmodules->Delete();
552         delete fITSmodules;
553     } // end fir fITSmoudles
554
555     if(!fDetTypeSim) {
556       Error("InitModules","fDetTypeSim is null!");
557       return kFALSE;
558     }
559
560     Int_t nl,indexMAX,index;
561
562     if(size<=0){ // default to using data stored in AliITSgeom
563         if(fDetTypeSim->GetITSgeom()==0) {
564             Error("InitModules","fITSgeom not defined");
565             return kFALSE;
566         } // end if fITSgeom==0
567         nl = fDetTypeSim->GetITSgeom()->GetNlayers();
568         indexMAX = fDetTypeSim->GetITSgeom()->GetIndexMax();
569         nmodules = indexMAX;
570         fITSmodules = new TObjArray(indexMAX);
571         for(index=0;index<indexMAX;index++){
572             fITSmodules->AddAt( new AliITSmodule(index),index);
573         } // end for index
574     }else{
575         fITSmodules = new TObjArray(size);
576         for(index=0;index<size;index++) {
577             fITSmodules->AddAt( new AliITSmodule(index),index);
578         } // end for index
579
580         nmodules = size;
581     } // end i size<=0
582     return kTRUE;
583 }
584 //______________________________________________________________________
585 void AliITS::Hits2SDigits(){
586     // Standard Hits to summable Digits function.
587     // Inputs:
588     //      none.
589     // Outputs:
590     //      none.
591   
592
593    if(!fDetTypeSim) {
594      Error("Hits2SDigits","fDetTypeSim is null!");
595      return;
596     
597   } 
598      
599   SetDefaults();
600   fLoader->LoadHits("read");
601   fLoader->LoadSDigits("recreate");
602   AliRunLoader* rl = fLoader->GetRunLoader(); 
603   fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
604   for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
605         // Do the Hits to Digits operation. Use Standard input values.
606         // Event number from file, no background hit merging , use size from
607         // AliITSgeom class, option="All", input from this file only.
608     rl->GetEvent(iEvent);
609     if (!fLoader->TreeS()) fLoader->MakeTree("S");
610     MakeBranch("S");
611     SetTreeAddress();
612     HitsToPreDigits(iEvent,0,-1," ",fOpt," ");
613   } // end for iEvent
614     
615   fLoader->UnloadHits();
616   fLoader->UnloadSDigits();
617   
618 }
619 //______________________________________________________________________
620 void AliITS::Hits2Digits(){
621
622   //Conversion from hits to digits
623   if(!fDetTypeSim) {
624     Error("Hits2SDigits","fDetTypeSim is 0!");
625     return;
626   }
627    
628   fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
629   SetDefaults();
630
631   fLoader->LoadHits("read");
632   fLoader->LoadDigits("recreate");
633   AliRunLoader* rl = fLoader->GetRunLoader(); 
634   for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
635     rl->GetEvent(iEvent);
636     if (!fLoader->TreeD()) fLoader->MakeTree("D");
637     MakeBranch("D");
638     SetTreeAddress();   
639     HitsToDigits(iEvent,0,-1," ",fOpt," ");
640   } 
641   
642   fLoader->UnloadHits();
643   fLoader->UnloadDigits();
644   
645 }
646
647 //______________________________________________________________________
648 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size,
649                           Option_t *option,Option_t *opt,
650                           const char *filename){
651     //   Keep galice.root for signal and name differently the file for 
652     // background when add! otherwise the track info for signal will be lost !
653     // the condition below will disappear when the geom class will be
654     // initialized for all versions - for the moment it is only for v5 !
655     // 7 is the SDD beam test version.
656     // Inputs:
657     //      Int_t evnt       Event to be processed.
658     //      Int_t bgrev      Background Hit tree number.
659     //      Int_t nmodules   Not used.
660     //      Option_t *option String indicating if merging hits or not. To
661     //                       merge hits set equal to "Add". Otherwise no
662     //                       background hits are considered.
663     //      Test_t *filename File name containing the background hits..
664     // Outputs:
665     //      none.
666     // Return:
667     //      none.
668
669   if(!fDetTypeSim) {
670     Error("HitsToDigits","fDetTypeSim is null!");
671     return;
672   }
673   fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
674   if(!GetITSgeom()) return; // need transformations to do digitization.
675   AliITSgeom *geom = GetITSgeom();
676
677   const char *all = strstr(opt,"All");
678   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
679                         strstr(opt,"SSD")};
680   static Bool_t setDef=kTRUE;
681   if (setDef) SetDefaultSimulation();
682   setDef=kFALSE;
683   
684   Int_t nmodules;
685   InitModules(size,nmodules);
686   FillModules(evNumber,bgrev,nmodules,option,filename);
687
688   AliITSsimulation *sim      = 0;
689   AliITSmodule     *mod      = 0;
690   Int_t id;
691   for(Int_t module=0;module<geom->GetIndexMax();module++){
692     id       = geom->GetModuleType(module);
693     if (!all && !det[id]) continue;
694     sim      = (AliITSsimulation*)fDetTypeSim->GetSimulationModel(id);
695     if (!sim) {
696       Error("HitsToDigits","The simulation class was not "
697             "instanciated for module %d type %s!",module,
698             geom->GetModuleTypeName(module));
699       exit(1);
700     } // end if !sim
701     mod      = (AliITSmodule *)fITSmodules->At(module);
702     sim->DigitiseModule(mod,module,evNumber);
703     // fills all branches - wasted disk space
704     fLoader->TreeD()->Fill(); 
705     ResetDigits();
706   } // end for module
707   
708   ClearModules();
709   
710   fLoader->TreeD()->GetEntries();
711   fLoader->TreeD()->AutoSave();
712   // reset tree
713   fLoader->TreeD()->Reset();
714 }
715 //_____________________________________________________________________
716 void AliITS::Hits2PreDigits(){ 
717   // Turn hits into SDigits
718
719   if(!fDetTypeSim) {
720     Error("Hits2SDigits","fDetTypeSim is 0!");
721     return;
722   }
723    
724   fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
725   SetDefaults();
726   
727   HitsToPreDigits(fLoader->GetRunLoader()->GetEventNumber(),
728                   0,-1," ",fOpt," ");
729 }
730
731 //______________________________________________________________________
732 void AliITS::HitsToPreDigits(Int_t evNumber,Int_t bgrev,Int_t size,
733                              Option_t *option,Option_t *opt,
734                              const char *filename){
735     //   Keep galice.root for signal and name differently the file for 
736     // background when add! otherwise the track info for signal will be lost !
737     // the condition below will disappear when the geom class will be
738     // initialized for all versions - for the moment it is only for v5 !
739     // 7 is the SDD beam test version.
740     // Inputs:
741     //      Int_t evnt       Event to be processed.
742     //      Int_t bgrev      Background Hit tree number.
743     //      Int_t nmodules   Not used.
744     //      Option_t *option String indicating if merging hits or not. To
745     //                       merge hits set equal to "Add". Otherwise no
746     //                       background hits are considered.
747     //      Test_t *filename File name containing the background hits..
748     // Outputs:
749     //      none.
750     // Return:
751     //      none.
752
753  
754   if(!fDetTypeSim) {
755     Error("HitsToPreDigits","fDetTypeSim is null!");
756     return;
757   }
758   fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
759
760   if(!GetITSgeom()){
761     Error("HitsToPreDigits","fGeom is null!");
762     return; // need transformations to do digitization.
763   }
764   AliITSgeom *geom = GetITSgeom();
765
766   const char *all = strstr(opt,"All");
767   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
768                         strstr(opt,"SSD")};
769   static Bool_t setDef=kTRUE;
770   if (setDef) SetDefaultSimulation();
771   setDef=kFALSE;
772   
773   Int_t nmodules;
774   InitModules(size,nmodules);
775   FillModules(evNumber,bgrev,nmodules,option,filename);
776   
777
778   AliITSsimulation *sim      = 0;
779   AliITSmodule     *mod      = 0;
780   Int_t id,module;
781   for(module=0;module<geom->GetIndexMax();module++){
782     id       = geom->GetModuleType(module);
783     if (!all && !det[id]) continue;
784     sim      = (AliITSsimulation*)GetSimulationModel(id);
785     if (!sim) {
786       Error("HitsToPreDigits","The simulation class was not "
787             "instanciated for module %d type %s!",module,
788             geom->GetModuleTypeName(module));
789       exit(1);
790     } // end if !sim
791     mod      = (AliITSmodule *)fITSmodules->At(module);
792     sim->SDigitiseModule(mod,module,evNumber);
793     // fills all branches - wasted disk space
794     fLoader->TreeS()->Fill(); 
795     fDetTypeSim->ResetSDigits();
796   } // end for module
797
798   ClearModules();
799
800   
801   fLoader->TreeS()->GetEntries();
802   fLoader->TreeS()->AutoSave();
803   fLoader->WriteSDigits("OVERWRITE");
804   // reset tree
805   fLoader->TreeS()->Reset();
806 }
807
808 //_____________________________________________________________________
809 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
810                                   Option_t *opt0,Option_t *opt1,
811                                  const char *flnm){
812     // keep galice.root for signal and name differently the file for 
813     // background when add! otherwise the track info for signal will be lost !
814     // the condition below will disappear when the geom class will be
815     // initialized for all versions - for the moment it is only for v5 !
816     // Inputs:
817     //      Int_t evnt       Event to be processed.
818     //      Int_t bgrev      Background Hit tree number.
819     //      Int_t size       Size used by InitModules. See InitModules.
820     //      Option_t *opt0   Option passed to FillModules. See FillModules.
821     //      Option_t *opt1   String indicating if merging hits or not. To
822     //                       merge hits set equal to "Add". Otherwise no
823     //                       background hits are considered.
824     //      Test_t *flnm     File name containing the background hits..
825     // Outputs:
826     //      none.
827     // Return:
828     //      none.
829
830
831
832   if(!GetITSgeom()){
833     Error("HitsToPreDigits","fGeom is null!");
834     return; // need transformations to do digitization.
835   }
836   AliITSgeom *geom = GetITSgeom();
837
838   AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
839
840   const char *all = strstr(opt1,"All");
841   const char *det[3] ={strstr(opt1,"SPD"),strstr(opt1,"SDD"),
842                        strstr(opt1,"SSD")};
843   Int_t nmodules;
844   InitModules(size,nmodules);
845   FillModules(evNumber,bgrev,nmodules,opt0,flnm);
846
847   AliITSsimulation *sim      = 0;
848   AliITSmodule     *mod      = 0;
849   Int_t id,module;
850
851   TTree *lTR = pITSloader->TreeR();
852   if(!lTR) {
853     pITSloader->MakeTree("R");
854     lTR = pITSloader->TreeR();
855   }
856   
857   TClonesArray* ptarray = new TClonesArray("AliITSRecPoint",1000);
858   TBranch* branch = (TBranch*)lTR->Branch("ITSRecPointsF",&ptarray);
859   branch->SetAddress(&ptarray);
860   //m.b. : this change is nothing but a nice way to make sure
861   //the CPU goes up !    
862   for(module=0;module<geom->GetIndexMax();module++){
863     id       = geom->GetModuleType(module);
864     if (!all && !det[id]) continue;
865     sim      = (AliITSsimulation*)GetSimulationModel(id);
866     if (!sim) {
867       Error("HitsToFastPoints","The simulation class was not "
868             "instanciated for module %d type %x!",module,
869             geom->GetModuleTypeName(module));
870       exit(1);
871     } // end if !sim
872     mod      = (AliITSmodule *)fITSmodules->At(module);
873     sim->CreateFastRecPoints(mod,module,gRandom,ptarray);
874     lTR->Fill();
875   } // end for module
876
877   ClearModules();
878   fLoader->WriteRecPoints("OVERWRITE");
879   delete ptarray;
880 }
881 //_____________________________________________________________________
882 Int_t AliITS::Hits2Clusters(TTree *hTree, TTree *cTree) {
883   //------------------------------------------------------------
884   // This function creates ITS clusters
885   //------------------------------------------------------------
886   if(!GetITSgeom()){
887     Error("HitsToPreDigits","fGeom is null!");
888     return 1; // need transformations to do digitization.
889   }
890   AliITSgeom *geom=GetITSgeom();
891   Int_t mmax=geom->GetIndexMax();
892
893   InitModules(-1,mmax);
894   FillModules(hTree,0);
895
896   TClonesArray *points = new TClonesArray("AliITSRecPoint",1000);
897   TBranch *branch=cTree->GetBranch("ITSRecPoints");
898   if (!branch) cTree->Branch("ITSRecPoints",&points);
899   else branch->SetAddress(&points);
900
901   AliITSsimulationFastPoints sim;
902   Int_t ncl=0;
903   for (Int_t m=0; m<mmax; m++) {
904     AliITSmodule *mod=GetModule(m);      
905     sim.CreateFastRecPoints(mod,m,gRandom,points);      
906     ncl+=points->GetEntriesFast();
907     cTree->Fill();
908     points->Clear();
909   }
910
911   Info("Hits2Clusters","Number of found fast clusters : %d",ncl);
912
913   //cTree->Write();
914
915   delete points;
916   return 0;
917 }
918
919 //_____________________________________________________________________
920 void AliITS::CheckLabels(Int_t lab[3]) const {
921   //------------------------------------------------------------
922   // Tries to find mother's labels
923   //------------------------------------------------------------
924
925   if(lab[0]<0 && lab[1]<0 && lab[2]<0) return; // In case of no labels just exit
926
927   Int_t ntracks = gAlice->GetMCApp()->GetNtrack();
928   for (Int_t i=0;i<3;i++){
929     Int_t label = lab[i];
930     if (label>=0 && label<ntracks) {
931       TParticle *part=(TParticle*)gAlice->GetMCApp()->Particle(label);
932       if (part->P() < 0.005) {
933         Int_t m=part->GetFirstMother();
934         if (m<0) {      
935           continue;
936         }
937         if (part->GetStatusCode()>0) {
938           continue;
939         }
940         lab[i]=m;       
941       }
942     }    
943   }
944   
945 }
946
947 //______________________________________________________________________
948 void AliITS::SDigitsToDigits(Option_t *opt){
949     // Standard Summable digits to Digits function.
950     // Inputs:
951     //      none.
952     // Outputs:
953     //      none.
954     if(!fDetTypeSim) {
955       Error("SDigitsToSDigits","fDetTypeSim is 0!");
956       return;
957     }
958    
959     fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
960     SetDefaults();
961     fDetTypeSim->SDigitsToDigits(opt,(Char_t*)GetName());
962
963 }
964
965 //______________________________________________________________________
966 void AliITS::ResetDigits(){
967     // Reset number of digits and the digits array for the ITS detector.
968     // Inputs:
969     //      none.
970     // Outputs:
971     //      none.
972     if(!fDetTypeSim) {
973       Error("ResetDigits","fDetTypeSim is 0!");
974       return;
975     }
976    
977     fDetTypeSim->ResetDigits();
978
979
980 }
981 //______________________________________________________________________
982 void AliITS::ResetDigits(Int_t branch){
983     // Reset number of digits and the digits array for this branch.
984     // Inputs:
985     //      none.
986     // Outputs:
987     //      none.
988
989     if(!fDetTypeSim) {
990       Error("ResetDigits","fDetTypeSim is 0!");
991       return;
992     }
993    
994     fDetTypeSim->ResetDigits(branch);
995
996 }
997 //______________________________________________________________________
998 void AliITS::AddSumDigit(AliITSpListItem &sdig){
999     // Adds the a module full of summable digits to the summable digits tree.
1000     // Inputs:
1001     //      AliITSpListItem &sdig   SDigit to be added to SDigits tree.
1002     // Outputs:
1003     //      none.
1004     // Return:
1005     //      none.
1006
1007     if(!fDetTypeSim) {
1008       Error("AddSumDigit","fDetTypeSim is 0!");
1009       return;
1010     }
1011     fDetTypeSim->AddSumDigit(sdig);
1012     
1013 }
1014 //______________________________________________________________________
1015 void AliITS::AddSimDigit(Int_t branch, AliITSdigit *d){
1016     //    Add a simulated digit.
1017     // Inputs:
1018     //      Int_t id        Detector type number.
1019     //      AliITSdigit *d  Digit to be added to the Digits Tree. See 
1020     //                      AliITSdigit.h
1021     // Outputs:
1022     //      none.
1023     // Return:
1024     //      none.
1025
1026     if(!fDetTypeSim) {
1027       Error("AddSimDigit","fDetTypeSim is 0!");
1028       return;
1029     }
1030     fDetTypeSim->AddSimDigit(branch,d);
1031
1032 }
1033 //______________________________________________________________________
1034 void AliITS::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,Int_t *tracks,
1035                          Int_t *hits,Float_t *charges, Int_t sigexpanded){
1036   //   Add a simulated digit to the list.
1037   // Inputs:
1038   //      Int_t id        Detector type number.
1039   //      Float_t phys    Physics indicator. See AliITSdigits.h
1040   //      Int_t *digits   Integer array containing the digits info. See 
1041   //                      AliITSdigit.h
1042   //      Int_t *tracks   Integer array [AliITSdigitS?D::GetNTracks()] 
1043   //                      containing the track numbers that contributed to
1044   //                      this digit.
1045   //      Int_t *hits     Integer array [AliITSdigitS?D::GetNTracks()]
1046   //                      containing the hit numbers, from AliITSmodule, that
1047   //                      contributed to this digit.
1048   //      Float_t *charge Floating point array of the signals contributed
1049   //                      to this digit by each track.
1050   // Outputs:
1051   //      none.
1052   // Return:
1053   //      none.
1054
1055     if(!fDetTypeSim) {
1056       Error("AddSimDigit","fDetTypeSim is 0!");
1057       return;
1058     }
1059     fDetTypeSim->AddSimDigit(branch,phys,digits,tracks,hits,charges,sigexpanded);
1060
1061 }
1062 //______________________________________________________________________
1063 void AliITS::Digits2Raw(){
1064     // convert digits of the current event to raw data
1065
1066   if(!fDetTypeSim) {
1067     Error("Digits2Raw","fDetTypeSim is 0!");
1068     return;
1069   }
1070   fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
1071   SetDefaults();
1072   fDetTypeSim->GetLoader()->LoadDigits();
1073   TTree* digits = fDetTypeSim->GetLoader()->TreeD();
1074   if (!digits) {
1075       Error("Digits2Raw", "no digits tree");
1076       return;
1077   }
1078   fDetTypeSim->SetTreeAddressD(digits,(Char_t*)GetName());
1079   AliITSDDLModuleMapSDD* ddlsdd=fDetTypeSim->GetDDLModuleMapSDD();
1080   Bool_t isHLTmodeC=fDetTypeSim->IsHLTmodeC();
1081   AliITSDDLRawData rawWriter;
1082   if(isHLTmodeC) rawWriter.SetUseCompressedSDDFormat(kTRUE);
1083   else rawWriter.SetUseCompressedSDDFormat(kFALSE);
1084   //Verbose level
1085   // 0: Silent
1086   // 1: cout messages
1087   // 2: txt files with digits 
1088   //BE CAREFUL, verbose level 2 MUST be used only for debugging and
1089   //it is highly suggested to use this mode only for debugging digits files
1090   //reasonably small, because otherwise the size of the txt files can reach
1091   //quickly several MB wasting time and disk space.
1092   rawWriter.SetVerbose(0);
1093     
1094   //SILICON PIXEL DETECTOR
1095   Info("Digits2Raw", "Formatting raw data for SPD");
1096   rawWriter.RawDataSPD(digits->GetBranch("ITSDigitsSPD"));
1097     
1098   //SILICON DRIFT DETECTOR
1099   Info("Digits2Raw", "Formatting raw data for SDD");
1100   rawWriter.RawDataSDD(digits->GetBranch("ITSDigitsSDD"),ddlsdd);
1101     
1102   //SILICON STRIP DETECTOR
1103   Info("Digits2Raw", "Formatting raw data for SSD");
1104   rawWriter.RawDataSSD(digits->GetBranch("ITSDigitsSSD"));
1105
1106   fLoader->UnloadDigits();
1107 }
1108 //______________________________________________________________________
1109 AliLoader* AliITS::MakeLoader(const char* topfoldername){ 
1110     //builds ITSgetter (AliLoader type)
1111     //if detector wants to use castomized getter, it must overload this method
1112
1113     AliDebug(1,Form("Creating AliITSLoader. Top folder is %s.",
1114          topfoldername));
1115     fLoader = new AliITSLoader(GetName(),topfoldername);
1116     return fLoader;
1117 }
1118
1119 Bool_t AliITS::Raw2SDigits(AliRawReader* rawReader)
1120 {
1121   //
1122   // Converts RAW data to SDigits
1123   //
1124   // Get TreeS
1125   //
1126     Int_t last   = -1;
1127     Int_t size   = GetITSgeom()->GetIndexMax();
1128     TClonesArray** modA = new TClonesArray*[size];
1129     for (Int_t mod = 0; mod < size; mod++) modA[mod] = new TClonesArray("AliITSpListItem", 10000);
1130     
1131     AliLoader* loader =  (gAlice->GetRunLoader())->GetLoader("ITSLoader");
1132     if (!loader)
1133     {
1134         Error("Open","Can not get ITS loader from Run Loader");
1135         return kFALSE;
1136     }
1137
1138     TTree* tree = 0;
1139     tree = loader->TreeS();
1140     if (!tree)
1141     {
1142         loader->MakeTree("S");
1143         tree = loader->TreeS();
1144     }
1145     //
1146     // Array for SDigits
1147     // 
1148     TClonesArray aSDigits("AliITSpListItem",10000), *itsSDigits=&aSDigits;
1149     Int_t bufsize = 32000;
1150     tree->Branch("ITS", &itsSDigits, bufsize);
1151     Int_t npx = 0;
1152     //
1153     // SPD
1154     //
1155     AliITSsegmentationSPD* segSPD = (AliITSsegmentationSPD*) fDetTypeSim->GetSegmentationModel(0);
1156     if(!segSPD){
1157       AliWarning("Set AliITS defaults");
1158       SetDefaults();
1159       segSPD = (AliITSsegmentationSPD*) fDetTypeSim->GetSegmentationModel(0);
1160     }
1161     npx = segSPD->Npx();
1162     Double_t thr, sigma; 
1163     
1164     AliITSRawStreamSPD inputSPD(rawReader);
1165     while(1){
1166         Bool_t next  = inputSPD.Next();
1167         if (!next) break;
1168
1169         Int_t module = inputSPD.GetModuleID();
1170         Int_t column = inputSPD.GetColumn();
1171         Int_t row    = inputSPD.GetRow();
1172         Int_t index  = npx * column + row;
1173
1174         if (module >= size) continue;
1175  
1176         last = (modA[module])->GetEntries();
1177         TClonesArray& dum = *modA[module];
1178         fDetTypeSim->GetSimuParam()->SPDThresholds(module,thr,sigma);
1179         thr += 1.;
1180         new (dum[last]) AliITSpListItem(-1, -1, module, index, thr);
1181     }
1182     rawReader->Reset();
1183
1184     //
1185     // SDD
1186     // 
1187     AliITSsegmentationSDD* segSDD = (AliITSsegmentationSDD*) fDetTypeSim->GetSegmentationModel(1);
1188     npx = segSDD->Npx();
1189     AliITSRawStreamSDD inputSDD(rawReader);
1190     AliITSDDLModuleMapSDD* ddlmap=fDetTypeSim->GetDDLModuleMapSDD();
1191     inputSDD.SetDDLModuleMap(ddlmap);
1192     while(1){
1193         Bool_t next  = inputSDD.Next();
1194         if (!next) break;
1195
1196         Int_t module = inputSDD.GetModuleID();
1197         Int_t anode  = inputSDD.GetAnode();
1198         Int_t time   = inputSDD.GetTime();
1199         Int_t signal10 = inputSDD.GetSignal();
1200         Int_t index  = npx * anode + time;
1201
1202         if (module >= size) continue;
1203         last = modA[module]->GetEntries();
1204         TClonesArray& dum = *modA[module];
1205         new (dum[last]) AliITSpListItem(-1, -1, module, index, Double_t(signal10));
1206         ((AliITSpListItem*) dum.At(last))->AddSignalAfterElect(module, index, Double_t(signal10));
1207         
1208     }
1209     rawReader->Reset();
1210
1211     //
1212     // SSD
1213     // 
1214     AliITSsegmentationSSD* segSSD = (AliITSsegmentationSSD*) fDetTypeSim->GetSegmentationModel(2);
1215     npx = segSSD->Npx();
1216     AliITSRawStreamSSD inputSSD(rawReader);
1217     while(1){
1218         Bool_t next  = inputSSD.Next();
1219         if (!next) break;
1220
1221         Int_t module  = inputSSD.GetModuleID();
1222         Int_t side    = inputSSD.GetSideFlag();
1223         Int_t strip   = inputSSD.GetStrip();
1224         Int_t signal  = inputSSD.GetSignal();
1225         Int_t index  = npx * side + strip;
1226
1227         if (module >= size) continue;
1228         
1229         last = modA[module]->GetEntries();
1230         TClonesArray& dum = *modA[module];
1231         new (dum[last]) AliITSpListItem(-1, -1, module, index, Double_t(signal));
1232     }
1233     rawReader->Reset();
1234      AliITSpListItem* sdig = 0;
1235     
1236     for (Int_t mod = 0; mod < size; mod++)
1237     {
1238         Int_t nsdig =  modA[mod]->GetEntries();
1239         for (Int_t ie = 0; ie < nsdig; ie++) {
1240             sdig = (AliITSpListItem*) (modA[mod]->At(ie));
1241             new (aSDigits[ie]) AliITSpListItem(-1, -1, mod, sdig->GetIndex(), sdig->GetSignal());
1242             Float_t sig = sdig->GetSignalAfterElect();
1243             if (sig > 0.) {
1244                 sdig = (AliITSpListItem*)aSDigits[ie];
1245                 sdig->AddSignalAfterElect(mod, sdig->GetIndex(), Double_t(sig));
1246             }
1247         }
1248         
1249         tree->Fill();
1250         aSDigits.Clear();
1251         modA[mod]->Clear();
1252     }
1253     loader->WriteSDigits("OVERWRITE");    
1254     delete modA;
1255     return kTRUE;
1256 }
1257
1258
1259 //______________________________________________________________________
1260 void AliITS::UpdateInternalGeometry(){
1261
1262   //reads new geometry from TGeo 
1263 //   AliDebug(1,"Delete ITSgeom and create a new one reading TGeo");
1264
1265   AliITSVersion_t version = (AliITSVersion_t)IsVersion();
1266   Int_t minor = 0;
1267   if(version==kvPPRasymmFMD) minor=2;  // default minor version for this geom.
1268   AliITSInitGeometry initgeom;
1269   AliITSgeom* geom = initgeom.CreateAliITSgeom(version,minor);
1270   SetITSgeom(geom);
1271 }
1272