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