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