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