]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSU.cxx
improved feeddown paramaterisation
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSU.cxx
1 /*
2 */
3
4 /**************************************************************************
5  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
6  *                                                                        *
7  * Author: The ALICE Off-line Project.                                    *
8  * Contributors are mentioned in the code where appropriate.              *
9  *                                                                        *
10  * Permission to use, copy, modify and distribute this software and its   *
11  * documentation strictly for non-commercial purposes is hereby granted   *
12  * without fee, provided that the above copyright notice appears in all   *
13  * copies and that both the copyright notice and this permission notice   *
14  * appear in the supporting documentation. The authors make no claims     *
15  * about the suitability of this software for any purpose. It is          *
16  * provided "as is" without express or implied warranty.                  *
17  **************************************************************************/
18
19 /* $Id: AliITSU.cxx $ */
20
21
22 ///////////////////////////////////////////////////////////////////////////////
23 //                                                                           //
24 //      An overview of the basic philosophy of the ITS code development      //
25 // and analysis is show in the figure below.                                 //
26 //Begin_Html                                                                 //
27 /*                                               
28 <img src="picts/ITS/ITS_Analysis_schema.gif">
29 </pre>
30 <br clear=left>
31 <font size=+2 color=red>
32 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
33 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
34 </font>
35 <pre>
36 */
37 //End_Html
38 //
39 //  AliITSU. Inner Traking System base class.
40 //  This class contains the base procedures for the Inner Tracking System
41 //
42 //Begin_Html
43 /*
44 <img src="picts/ITS/AliITS_Class_Diagram.gif">
45 </pre>
46 <br clear=left>
47 <font size=+2 color=red>
48 <p>This show the class diagram of the different elements that are part of
49 the AliITS class.
50 </font>
51 <pre>
52 */
53 //End_Html
54 //
55 // Version: 0
56 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
57 //
58 // Version: 1
59 // Modified and documented by Bjorn S. Nilsen
60 // July 11 1999
61 //
62 // Version: 2
63 // Modified and documented by A. Bologna
64 // October 18 1999
65 //
66 // AliITSU is the general base class for the ITS. Also see AliDetector for
67 // futher information.
68 //
69 ///////////////////////////////////////////////////////////////////////////////
70
71 #include <stdlib.h>
72 #include <TClonesArray.h>
73 #include <TFile.h>
74 #include <TParticle.h>
75 #include <TString.h>
76 #include <TTree.h>
77 #include <TVirtualMC.h>
78 #include "AliDetector.h"
79 #include "AliITSU.h"
80 #include "AliITSLoader.h"
81 #include "AliITSULoader.h"
82 #include "AliITSUHit.h"
83 #include "AliITSUSDigit.h"
84 #include "AliITSUSimulation.h"
85 #include "AliITSUSimulationPix.h"
86 #include "AliITSsimulationFastPoints.h"
87 #include "AliMC.h"
88 #include "AliITSUDigitizer.h"
89 #include "AliITSRecPoint.h"
90 #include "AliRawReader.h"
91 #include "AliRun.h"
92 #include "AliLog.h"
93 #include "AliITSdigit.h"
94 #include "AliITSUModule.h"
95 #include "AliITSUDigitPix.h"
96 #include "AliITSsegmentation.h"
97 #include "AliITSUSegmentationPix.h"
98 #include "AliITSUSimuParam.h"
99 #include "AliITSFOSignalsSPD.h"
100 #include "AliITSUParamList.h"
101 #include "AliCDBManager.h" // tmp! Later the simuparam should be loaded centrally
102 #include "AliCDBEntry.h"
103
104 ClassImp(AliITSU)
105
106 //______________________________________________________________________
107 AliITSU::AliITSU() : 
108 AliDetector()
109   ,fEuclidOut(0)
110   ,fNLayers(0)
111   ,fIdSens(0)
112   ,fLayerName(0)
113   ,fTiming(kFALSE)
114   ,fGeomTGeo(0)
115   ,fSimuParam(0)
116   ,fModA(0)
117   ,fpSDigits(0)
118   ,fSDigits(0)
119   ,fDetHits(0)
120   ,fModuleHits(0)
121   ,fDetDigits(0)
122   ,fSensMap(0)
123   ,fSimModelLr(0)
124   ,fSegModelLr(0)
125   ,fResponseLr(0)
126   ,fCalibration(0)
127   ,fRunNumber(0)
128   ,fSimInitDone(kFALSE)
129 {
130   // Default initializer for ITS  
131 }
132
133 //______________________________________________________________________
134 AliITSU::AliITSU(const Char_t *title, Int_t nlay) :
135   AliDetector("ITS",title)
136   ,fEuclidOut(0)
137   ,fNLayers(nlay)
138   ,fIdSens(0)
139   ,fLayerName(0)
140   ,fTiming(kFALSE)
141   ,fGeomTGeo(0)
142   ,fSimuParam(0)
143   ,fModA(0)
144   ,fpSDigits(0)
145   ,fSDigits(0)
146   ,fDetHits(0)
147   ,fModuleHits(0)
148   ,fDetDigits(0)
149   ,fSensMap(0)
150   ,fSimModelLr(0)
151   ,fSegModelLr(0)
152   ,fResponseLr(0)
153   ,fCalibration(0)
154   ,fRunNumber(0)
155   ,fSimInitDone(kFALSE)
156 {
157   //     The standard Constructor for the ITS class. 
158   AliMC* mc = gAlice->GetMCApp();
159   if( mc && mc->GetHitLists() ) {
160     fHits = new TClonesArray("AliITSUHit",100); // from AliDetector
161     mc->AddHitList(fHits);  
162   }
163 }
164
165
166 //______________________________________________________________________
167 AliITSU::~AliITSU()
168 {
169   // Default destructor for ITS.
170   //  
171   delete fHits;
172   //  delete fSimuParam; // provided by the CDBManager
173   delete fSensMap;
174   if (fSimModelLr) {
175     for (int i=fNLayers;i--;) { // different layers may use the same simulation model
176       for (int j=i;j--;) if (fSimModelLr[j]==fSimModelLr[i]) fSimModelLr[j] = 0;
177       delete fSimModelLr[i];
178     }
179     delete[] fSimModelLr;
180   }
181   if (fSegModelLr) {
182     for (int i=fNLayers;i--;) { // different layers may use the same simulation model     
183       for (int j=i;j--;) if (fSegModelLr[j]==fSegModelLr[i]) fSegModelLr[j] = 0;
184       delete fSegModelLr[i];
185     }
186     delete[] fSegModelLr;
187   }
188   //
189   delete fResponseLr; // note: the response data is owned by the CDBManager, we don't delete them
190   //
191   delete[] fLayerName;  // Array of TStrings
192   delete[] fIdSens;
193   //
194   int nmod = fGeomTGeo ? fGeomTGeo->GetNModules() : 0;
195   if (fModuleHits) fModuleHits->Delete();
196   delete fModuleHits;
197   if(fModA){
198     for(Int_t j=0; j<nmod; j++){
199       fModA[j]->Delete();
200       delete fModA[j];
201     }
202     delete[] fModA;
203   }
204   if (fpSDigits) { fpSDigits->Delete(); delete fpSDigits; }
205   if (fSDigits)  { fSDigits->Delete(); delete fSDigits; }
206   delete fGeomTGeo;
207   //
208 }
209
210 //______________________________________________________________________
211 AliDigitizer* AliITSU::CreateDigitizer(AliDigitizationInput* manager) const
212 {
213   // Creates the AliITSDigitizer in a standard way for use via AliModule.
214   // This function can not be included in the .h file because of problems
215   // with the order of inclusion (recursive).
216   // Inputs:
217   //    AliDigitizationInput *manager  The Manger class for Digitization
218   // Output:
219   //    none.
220   // Return:
221   //    A new AliITSRunDigitizer (cast as a AliDigitizer).
222   //
223   return new AliITSUDigitizer(manager);
224 }
225
226 //______________________________________________________________________
227 void AliITSU::Init() 
228 {
229   // Initializer ITS after it has been built
230   //     This routine initializes the AliITS class. It is intended to be
231   // called from the Init function in AliITSv?. Besides displaying a banner
232   // indicating that it has been called it initializes the array fIdSens
233   // and sets the default segmentation, response, digit and raw cluster
234   // classes therefore it should be called after a call to CreateGeometry.
235   //
236   if (!fIdSens) fIdSens = new Int_t[fNLayers];
237   for(int i=0;i<fNLayers;i++) fIdSens[i] = gMC ? gMC->VolId(fLayerName[i]) : 0;
238   fGeomTGeo     = new AliITSUGeomTGeo(kTRUE);
239   InitSimulation();
240   //
241 }
242
243 //______________________________________________________________________
244 void AliITSU::MakeBranch(Option_t* option)
245 {
246   // Creates Tree branches for the ITS.
247   // Inputs:
248   //      Option_t *option    String of Tree types S,D, and/or R.
249   //      const char *file    String of the file name where these branches
250   //                          are to be stored. If blank then these branches
251   //                          are written to the same tree as the Hits were
252   //                          read from.
253   // Outputs:
254   //      none.
255   // Return:
256   //      none.
257   
258   Bool_t cH = (strstr(option,"H")!=0);
259   Bool_t cS = (strstr(option,"S")!=0);
260   Bool_t cD = (strstr(option,"D")!=0);
261   
262   if(cH && (fHits == 0x0)) fHits = new TClonesArray("AliITSUHit", 1560);
263   AliDetector::MakeBranch(option);
264   
265   if(cS) MakeBranchS(0);
266   if(cD) MakeBranchD(0);
267   // 
268 }
269
270 //___________________________________________________________________
271 void AliITSU::MakeBranchS(const char* fl) 
272 {
273   // Creates Tree Branch for the ITS summable digits.
274   // Inputs:
275   //      cont char *fl  File name where SDigits branch is to be written
276   //                     to. If blank it write the SDigits to the same
277   //                     file in which the Hits were found.
278   //  
279   Int_t buffersize = 4000;
280   char branchname[31];
281   //
282   // only one branch for SDigits.
283   snprintf(branchname,30,"%s",GetName());
284   if(fLoader->TreeS()) MakeBranchInTree(fLoader->TreeS(),branchname,&fSDigits,buffersize,fl);
285   //
286 }
287
288 //______________________________________________________________________
289 void AliITSU::MakeBranchD(const char* file)
290 {
291   //Make branch for digits
292   MakeBranchInTreeD(fLoader->TreeD(),file);
293 }
294
295 //___________________________________________________________________
296 void AliITSU:: MakeBranchInTreeD(TTree* treeD, const char* file)
297 {
298   // Creates Tree branches for the ITS.
299   //
300   if (!treeD) {AliFatal("No tree provided");}
301   Int_t buffersize = 4000;
302   if (!fDetDigits) InitArrays();
303   //
304   for (Int_t i=0;i<kNDetTypes;i++) {
305     ResetDigits(i);
306     TClonesArray* darr = (TClonesArray*)fDetDigits->At(i);
307     AliDetector::MakeBranchInTree(treeD,Form("%sDigits%s",GetName(),fGeomTGeo->GetDetTypeName(i)),
308                                   &darr,buffersize,file);
309   }
310   //
311 }
312
313 //______________________________________________________________________
314 void AliITSU::InitArrays()
315 {
316   // initialize arrays
317   //
318   if(!fLoader) MakeLoader(AliConfig::GetDefaultEventFolderName());
319   //  
320   fDetDigits = new TObjArray(kNDetTypes);
321   for (Int_t i=0;i<kNDetTypes;i++) fDetDigits->AddAt(new TClonesArray(GetDigitClassName(i),100),i);
322   //
323   fSDigits = new TClonesArray("AliITSUSDigit",100);
324   //
325   fDetHits = new TClonesArray("AliITSUHit",100);
326   //
327   fModuleHits = new TObjArray(fGeomTGeo->GetNModules());
328   for (int i=0;i<fGeomTGeo->GetNModules();i++) fModuleHits->AddLast( new AliITSUModule(i,fGeomTGeo) );
329   //
330 }
331
332 //______________________________________________________________________
333 void AliITSU::SetTreeAddress()
334 {
335   // Set branch address for the Trees.
336   TTree *treeS = fLoader->TreeS();
337   if (treeS) {
338     TBranch* br = treeS->GetBranch(GetName());
339     if (br) br->SetAddress(&fSDigits);
340   }
341   //
342   TTree *treeD = fLoader->TreeD();
343   if (treeD) {
344     if (!fDetDigits) InitArrays();
345     for (int i=0;i<kNDetTypes;i++) {
346       TString brname = Form("%sDigits%s",GetName(),GetDetTypeName(i));
347       TBranch* br = treeD->GetBranch(brname.Data());
348       if (!br) continue;
349       TClonesArray* darr = (TClonesArray*)fDetDigits->At(i);
350       br->SetAddress(&darr);
351     }
352   }
353   if (fLoader->TreeH() && (fHits == 0x0)) fHits = new TClonesArray("AliITSUHit", 1560);
354   AliDetector::SetTreeAddress();
355   //
356 }
357
358 //______________________________________________________________________
359 void AliITSU::AddHit(Int_t track, Int_t *vol, Float_t *hits)
360 {
361   // Add an ITS hit
362   //     The function to add information to the AliITSUHit class. See the
363   // AliITSUHit class for a full description. This function allocates the
364   // necessary new space for the hit information and passes the variable
365   // track, and the pointers *vol and *hits to the AliITSUHit constructor
366   // function.
367   // Inputs:
368   //      Int_t   track   Track number which produced this hit.
369   //      Int_t   *vol    Array of Integer Hit information. See AliITSUHit.h
370   //      Float_t *hits   Array of Floating Hit information.  see AliITSUHit.h
371   TClonesArray &lhits = *fHits;
372   new(lhits[fNhits++]) AliITSUHit(fIshunt,track,vol,hits);
373   //
374 }
375
376 //______________________________________________________________________
377 void AliITSU::FillModules(Int_t bgrev, Option_t *option, const char *filename) 
378 {
379   // fill the modules with the sorted by module hits; add hits from
380   // background if option=Add.
381   //
382   static TTree *trH1=0;                 //Tree with background hits
383   static Bool_t first=kTRUE;
384   static TFile *file = 0;
385   const char *addBgr = strstr(option,"Add");
386   //
387   if (addBgr ) {
388     if(first) {
389       file = new TFile(filename);    
390       first=kFALSE;
391     }
392     file->cd();
393     file->ls();
394     // Get Hits Tree header from file
395     if (trH1) {delete trH1; trH1=0;}
396     //
397     char treeName[21];
398     snprintf(treeName,20,"TreeH%d",bgrev);
399     trH1 = (TTree*)gDirectory->Get(treeName);
400     if (!trH1) Error("FillModules","cannot find Hits Tree for event:%d",bgrev);
401     // Set branch addresses
402   } // end if addBgr
403   
404   FillModules(fLoader->TreeH(),0); // fill from this file's tree.
405   //
406   if (addBgr ) {
407     FillModules(trH1,10000000); // Default mask 10M.
408     TTree *fAli=fLoader->GetRunLoader()->TreeK();
409     TFile *fileAli=0;
410     if (fAli) {
411       fileAli = fAli->GetCurrentFile();
412       fileAli->cd();
413     }
414   } // end if add
415   //  
416 }
417
418 //______________________________________________________________________
419 void AliITSU::FillModules(TTree *treeH, Int_t /*mask*/)
420 {
421   // fill the modules with the sorted by module hits; 
422   // can be called many times to do a merging
423   // Inputs:
424   //      TTree *treeH  The tree containing the hits to be copied into
425   //                    the modules.
426   //      Int_t mask    The track number mask to indecate which file
427   //                    this hits came from.
428   //  
429   if (treeH == 0x0) { AliError("Tree H  is NULL"); return; }
430   //
431   Int_t lay,lad,det,index;
432   AliITSUHit *itsHit=0;
433   char branchname[21];
434   snprintf(branchname,20,"%s",GetName());
435   TBranch *branch = treeH->GetBranch(branchname);
436   if (!branch) {Error("FillModules","%s branch in TreeH not found",branchname); return;} // end if !branch
437   //
438   branch->SetAddress(&fHits);
439   Int_t nTracks =(Int_t) treeH->GetEntries();
440   Int_t iPrimTrack,h;
441   for (iPrimTrack=0; iPrimTrack<nTracks; iPrimTrack++) {
442     ResetHits();
443     Int_t nBytes = treeH->GetEvent(iPrimTrack);
444     if (nBytes <= 0) continue;
445     Int_t nHits = fHits->GetEntriesFast();
446     for (h=0; h<nHits; h++){
447       itsHit = (AliITSUHit *)fHits->UncheckedAt(h);
448       itsHit->GetDetectorID(lay,lad,det);
449       index = fGeomTGeo->GetModuleIndex(lay,lad,det); // !!! AliITSHit counts indices from 1!
450       itsHit = new( (*fDetHits)[fDetHits->GetEntriesFast()] ) AliITSUHit(*itsHit);
451       itsHit->SetUniqueID(h);
452       GetModule(index)->AddHit(itsHit);
453       // do we need to add a mask?
454       // itsHit->SetTrack(itsHit->GetTrack()+mask);
455     } // end loop over hits 
456   } // end loop over tracks
457 }
458
459 //______________________________________________________________________
460 void AliITSU::ClearModules()
461 {
462   // clear accumulated hits
463   if (!fModuleHits || !fDetHits) AliFatal("Hits accumulation arrays are not defined");
464   for (int i=fGeomTGeo->GetNModules();i--;) GetModule(i)->Clear();
465   fDetHits->Clear();
466 }
467
468 //______________________________________________________________________
469 void AliITSU::Hits2SDigits()
470 {
471   // Standard Hits to summable Digits function.
472   if (!IsSimInitDone()) InitSimulation();
473   fLoader->LoadHits("read");
474   fLoader->LoadSDigits("recreate");
475   AliRunLoader* rl = fLoader->GetRunLoader(); 
476   //
477   for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
478     rl->GetEvent(iEvent);
479     if (!fLoader->TreeS()) fLoader->MakeTree("S");
480     MakeBranch("S");
481     SetTreeAddress();
482     Hits2SDigits(iEvent,0," "," ");
483   } // end for iEvent
484     //
485   fLoader->UnloadHits();
486   fLoader->UnloadSDigits();
487   // 
488 }
489
490 //______________________________________________________________________
491 void AliITSU::Hits2SDigits(Int_t evNumber,Int_t bgrev,Option_t *option,const char *filename)
492 {
493   // Keep galice.root for signal and name differently the file for 
494   // background when add! otherwise the track info for signal will be lost !
495   // Inputs:
496   //      Int_t evnt       Event to be processed.
497   //      Int_t bgrev      Background Hit tree number.
498   //      Int_t nmodules   Not used.
499   //      Option_t *option String indicating if merging hits or not. To
500   //                       merge hits set equal to "Add". Otherwise no
501   //                       background hits are considered.
502   //      Test_t *filename File name containing the background hits..
503   //
504   if (!IsSimInitDone()) InitSimulation();
505   FillModules(bgrev,option,filename);
506   //
507   Int_t nmodules = fGeomTGeo->GetNModules();
508   
509   for(int module=0;module<nmodules;module++) {
510     int lr = fGeomTGeo->GetLayer(module);
511     AliITSUSimulation* sim = GetSimulationModel(lr);
512     sim->InitSimulationModule(GetModule(module),evNumber/*,gAlice->GetEvNumber()*/,GetSegmentation(lr),GetResponseParam(lr));
513     sim->SDigitiseModule();
514     fLoader->TreeS()->Fill();      // fills all branches - wasted disk space
515     ResetSDigits();
516   } 
517   //
518   ClearModules();
519   //
520   fLoader->TreeS()->GetEntries();
521   fLoader->TreeS()->AutoSave();
522   fLoader->WriteSDigits("OVERWRITE");
523   fLoader->TreeS()->Reset();
524 }
525
526 //______________________________________________________________________
527 void AliITSU::Hits2Digits()
528 {
529   // Standard Hits to Digits function.
530   if (!IsSimInitDone()) InitSimulation();
531   fLoader->LoadHits("read");
532   fLoader->LoadDigits("recreate");
533   AliRunLoader* rl = fLoader->GetRunLoader(); 
534   //
535   for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
536     rl->GetEvent(iEvent);
537     if (!fLoader->TreeS()) fLoader->MakeTree("S");
538     MakeBranch("D");
539     SetTreeAddress();
540     Hits2Digits(iEvent,0," "," ");
541   } // end for iEvent
542     //
543   fLoader->UnloadHits();
544   fLoader->UnloadSDigits();
545   // 
546 }
547
548 //______________________________________________________________________
549 void AliITSU::Hits2Digits(Int_t evNumber,Int_t bgrev,Option_t *option,const char *filename)
550 {
551   //   Keep galice.root for signal and name differently the file for 
552   // background when add! otherwise the track info for signal will be lost !
553   // Inputs:
554   //      Int_t evnt       Event to be processed.
555   //      Int_t bgrev      Background Hit tree number.
556   //      Option_t *option String indicating if merging hits or not. To
557   //                       merge hits set equal to "Add". Otherwise no
558   //                       background hits are considered.
559   //      Test_t *filename File name containing the background hits..
560   // Outputs:
561   //  
562   if (!IsSimInitDone()) InitSimulation();
563   FillModules(bgrev,option,filename);
564   // 
565   Int_t nmodules = fGeomTGeo->GetNModules();
566   for (Int_t module=0;module<nmodules;module++) {
567     int lr = fGeomTGeo->GetLayer(module);
568     AliITSUSimulation* sim = GetSimulationModel(lr);
569     sim->InitSimulationModule(GetModule(module),evNumber/*gAlice->GetEvNumber()*/,GetSegmentation(lr),GetResponseParam(lr));
570     sim->DigitiseModule();
571     // fills all branches - wasted disk space
572     fLoader->TreeD()->Fill(); 
573     ResetDigits();
574   } // end for module
575   //
576   ClearModules();
577   //
578   //    WriteFOSignals(); // Add Fast-OR signals to event (only one object per event)
579   fLoader->TreeD()->GetEntries();
580   fLoader->TreeD()->AutoSave();
581   fLoader->TreeD()->Reset(); 
582   //
583 }
584
585 //_____________________________________________________________________
586 void AliITSU::Hits2FastRecPoints(Int_t bgrev,Option_t *opt,const char *flnm)
587 {
588   // keep galice.root for signal and name differently the file for 
589   // background when add! otherwise the track info for signal will be lost !
590   // Inputs:
591   //      Int_t evnt       Event to be processed.
592   //      Int_t bgrev      Background Hit tree number.
593   //      Option_t *opt    Option passed to FillModules. See FillModules.
594   //      Test_t *flnm     File name containing the background hits..
595   // Outputs:
596   //      none.
597   // Return:
598   //      none.
599   if (!IsSimInitDone()) InitSimulation();
600   AliITSULoader *pITSloader = (AliITSULoader*)fLoader;
601   Int_t nmodules = fGeomTGeo->GetNModules();
602   FillModules(bgrev,opt,flnm);
603   //
604   TTree *lTR = pITSloader->TreeR();
605   if(!lTR) {
606     pITSloader->MakeTree("R");
607     lTR = pITSloader->TreeR();
608   }
609   //
610   TClonesArray* ptarray = new TClonesArray("AliITSRecPoint",1000);
611   TBranch* branch = (TBranch*)lTR->Branch("ITSRecPointsF",&ptarray);
612   branch->SetAddress(&ptarray);
613   for (int module=0;module<nmodules;module++){
614     int id = fGeomTGeo->GetModuleDetTypeID(module);
615     AliITSUSimulation* sim = GetSimulationModel(id);
616     if (!sim) AliFatal(Form("The sim.class for module %d of DetTypeID %d is missing",module,id));
617     sim->CreateFastRecPoints( GetModule(module) ,module,gRandom,ptarray);
618     lTR->Fill();
619     ptarray->Clear();
620   } // end for module
621   //
622   ClearModules();
623   fLoader->WriteRecPoints("OVERWRITE");
624   delete ptarray;
625 }
626
627 //_____________________________________________________________________
628 Int_t AliITSU::Hits2Clusters(TTree */*hTree*/, TTree */*cTree*/)
629 {
630   /* RS: TODO
631   // This function creates ITS clusters
632   if (!IsSimInitDone()) InitSimulation();
633   Int_t mmax = 0;
634   FillModules(hTree,0);
635   //
636   TClonesArray *points = new TClonesArray("AliITSRecPoint",1000);
637   TBranch *branch=cTree->GetBranch("ITSRecPoints");
638   if (!branch) cTree->Branch("ITSRecPoints",&points);
639   else branch->SetAddress(&points);
640   //
641   AliITSsimulationFastPoints sim;
642   Int_t ncl=0;
643   for (Int_t m=0; m<mmax; m++) {
644     sim.CreateFastRecPoints(GetModule(m),m,gRandom,points);      
645     ncl+=points->GetEntriesFast();
646     cTree->Fill();
647     points->Clear();
648   }
649   //
650   ClearModules();
651   //
652   AliDebug(1,Form("Number of found fast clusters : %d",ncl));
653   //cTree->Write();
654   delete points;
655   */
656   return 0;
657 }
658
659 //_____________________________________________________________________
660 void AliITSU::CheckLabels(Int_t lab[3]) const //RSDONE
661 {
662   // Tries to find mother's labels
663   //
664   if(lab[0]<0 && lab[1]<0 && lab[2]<0) return; // In case of no labels just exit
665   //
666   Int_t ntracks = gAlice->GetMCApp()->GetNtrack();
667   for (Int_t i=0;i<3;i++){
668     Int_t label = lab[i];
669     if (label>=0 && label<ntracks) {
670       TParticle *part=(TParticle*)gAlice->GetMCApp()->Particle(label);
671       if (part->P() < 0.005) {
672         Int_t m=part->GetFirstMother();
673         if (m<0) continue;
674         if (part->GetStatusCode()>0) continue;
675         lab[i]=m;       
676       }
677     }    
678   }
679   //
680 }
681
682 //______________________________________________________________________
683 void AliITSU::ResetDigits() //RSDONE?
684 {
685   // Reset number of digits and the digits array for the ITS detector.
686   if (fDetDigits) for (int i=kNDetTypes;i--;) ResetDigits(i);
687   //
688 }
689
690 //______________________________________________________________________
691 void AliITSU::ResetDigits(Int_t branch)
692 {
693   // Reset number of digits and the digits array for this branch.
694   if (fDetDigits) ((TClonesArray*)fDetDigits->At(branch))->Clear();
695   //
696 }
697
698 //______________________________________________________________________
699 void AliITSU::AddSumDigit(AliITSUSDigit &sdig)
700 {
701   // Adds the module summable digits to the summable digits tree.
702   new( (*fSDigits)[fSDigits->GetEntriesFast()]) AliITSUSDigit(sdig);
703   //  
704 }
705
706 //______________________________________________________________________
707 void AliITSU::AddSimDigit(Int_t branch, AliITSdigit *d)
708 {
709   //    Add a simulated digit.
710   // Inputs:
711   //      Int_t id        Detector type number.
712   //      AliITSdigit *d  Digit to be added to the Digits Tree. See 
713   //                      AliITSdigit.h
714   TClonesArray &ldigits = *((TClonesArray*)fDetDigits->At(branch));
715   int nd = ldigits.GetEntriesFast();
716   switch(branch){
717   case AliITSUGeomTGeo::kDetTypePix:
718     new(ldigits[nd]) AliITSUDigitPix(*((AliITSUDigitPix*)d));
719     break;
720   default:
721     AliFatal(Form("Unknown digits branch %d",branch));
722   }
723 }
724
725 //______________________________________________________________________
726 void AliITSU::AddSimDigit(Int_t branch,Float_t /*phys*/,Int_t *digits,Int_t *tracks,
727                           Int_t *hits,Float_t */*charges*/, Int_t /*sigexpanded*/)
728 {
729   // Add a simulated digit to the list.
730   // Inputs:
731   //      Int_t id        Detector type number.
732   //      Float_t phys    Physics indicator. See AliITSdigits.h
733   //      Int_t *digits   Integer array containing the digits info. See 
734   //                      AliITSdigit.h
735   //      Int_t *tracks   Integer array [AliITSdigitS?D::GetNTracks()] 
736   //                      containing the track numbers that contributed to
737   //                      this digit.
738   //      Int_t *hits     Integer array [AliITSdigitS?D::GetNTracks()]
739   //                      containing the hit numbers, from AliITSmodule, that
740   //                      contributed to this digit.
741   //      Float_t *charge Floating point array of the signals contributed
742   //                      to this digit by each track.
743   TClonesArray &ldigits = *((TClonesArray*)fDetDigits->At(branch));
744   int nd = ldigits.GetEntriesFast();
745   switch(branch){
746   case AliITSUGeomTGeo::kDetTypePix:
747     new(ldigits[nd]) AliITSUDigitPix(digits,tracks,hits);
748     break;
749   default:
750     AliFatal(Form("Unknown digits branch %d",branch));
751   }  
752   //
753 }
754
755 //______________________________________________________________________
756 void AliITSU::Digits2Raw()
757 {
758   AliError("Not ready");
759 }
760
761 //______________________________________________________________________
762 AliLoader* AliITSU::MakeLoader(const char* topfoldername)
763
764   //builds ITSgetter (AliLoader type)
765   //if detector wants to use castomized getter, it must overload this method
766   
767   AliDebug(1,Form("Creating AliITSULoader. Top folder is %s.",topfoldername));
768   fLoader = new AliITSULoader(GetName(),topfoldername);
769   return fLoader;
770 }
771
772 //______________________________________________________________________
773 Bool_t AliITSU::Raw2SDigits(AliRawReader* /*rawReader*/)
774 {
775   AliError("Not ready");
776   return kFALSE;
777 }
778
779 //______________________________________________________________________
780 /*
781 AliTriggerDetector* AliITSU::CreateTriggerDetector() const 
782 {
783   // create an AliITSTrigger object (and set trigger conditions as input)
784   return new AliITSTrigger(fDetTypeSim->GetTriggerConditions());
785 }
786 */
787
788 //______________________________________________________________________
789 void AliITSU::WriteFOSignals()
790 {
791   // This method write FO signals in Digits tree both in Hits2Digits
792   // or SDigits2Digits
793   AliError("Not ready");
794   //  fDetTypeSim->ProcessNoiseForFastOr();
795 }
796
797 //_______________________________________________________________________
798 void AliITSU::SDigits2Digits()
799 {
800   // Standard Summable digits to Digits function.
801   //
802   if (!IsSimInitDone()) InitSimulation();
803   TTree* trees = fLoader->TreeS();
804   if( !(trees && fSDigits) ) AliFatal("Error: No trees or SDigits.");
805   TBranch* brchSDigits = trees->GetBranch(GetName());
806   //
807   int nmodules = fGeomTGeo->GetNModules();
808   for (int module=0;module<nmodules;module++) {
809     int lr = fGeomTGeo->GetLayer(module);
810     AliITSUSimulation* sim = GetSimulationModel(lr);
811     sim->InitSimulationModule(GetModule(module),gAlice->GetEvNumber(),GetSegmentation(lr),GetResponseParam(lr));
812     fSDigits->Clear();
813     brchSDigits->GetEvent(module);
814     sim->AddSDigitsToModule(fSDigits,0);
815     sim->FinishSDigitiseModule();
816     fLoader->TreeD()->Fill();
817     ResetDigits();
818   }
819   //  WriteFOSignals(); 
820   fLoader->TreeD()->GetEntries();
821   fLoader->TreeD()->AutoSave();
822   fLoader->TreeD()->Reset();
823 }
824
825 //_______________________________________________________________________
826 void AliITSU::InitSimulation()
827 {
828   // Initialize arrays, segmentations ets, needed for simulation
829   // Equivalent of old AliITSDetTypeSim construction
830   //
831   if (fSimInitDone) {AliInfo("Already done"); return;}
832   //
833   AliCDBEntry* cdbEnt = AliCDBManager::Instance()->Get("ITS/Calib/SimuParam"); // tmp: load it centrally
834   if (!cdbEnt) {AliFatal("Failed to find ITS/Calib/SimuParam on CDB"); exit(1);}
835   fSimuParam    = (AliITSUSimuParam*)cdbEnt->GetObject();
836   //
837   fSensMap      = new AliITSUSensMap("AliITSUSDigit",0,0);
838   fSimModelLr   = new AliITSUSimulation*[fNLayers];
839   fSegModelLr   = new AliITSsegmentation*[fNLayers];
840   fResponseLr   = new AliITSUParamList*[fNLayers];
841   //
842   TObjArray arrSeg;
843   AliITSUSegmentationPix::LoadSegmentations(&arrSeg, AliITSUGeomTGeo::GetITSsegmentationFileName());
844   //
845   // add known simulation types used in the setup
846   for (int i=fNLayers;i--;) {
847     fSimModelLr[i] = 0;
848     fSegModelLr[i] = 0;
849     fResponseLr[i] = 0;
850     int dType = fGeomTGeo->GetLayerDetTypeID(i);           // fine detector type: class + segmentation
851     int sType = dType/AliITSUGeomTGeo::kMaxSegmPerDetType; // detector simulation class
852     //
853     // check if the simulation of this sType was already created for preceeding layers
854     AliITSUSimulation* simUpg = 0;
855     for (int j=fNLayers-1;j>i;j--) {
856       simUpg = GetSimulationModel(j);
857       if (simUpg && int(simUpg->GetUniqueID())==sType) break;
858       else simUpg = 0;
859     }
860     //
861     if (!simUpg) { // need to create simulation for detector class sType
862       switch (sType) 
863         {
864         case AliITSUGeomTGeo::kDetTypePix : 
865           simUpg = new AliITSUSimulationPix(fSimuParam,fSensMap); 
866           break;
867         default: AliFatal(Form("No %d detector type is defined",sType));
868         }
869     }
870     fSimModelLr[i] = simUpg;
871     //
872     // add segmentations used in the setup
873     if (!(fSegModelLr[i]=(AliITSsegmentation*)arrSeg[dType])) {AliFatal(Form("Segmentation for DetType#%d is not found",dType)); exit(1);}
874     //
875     // add response function for the detectors of this layer
876     if ( !(fResponseLr[i]=(AliITSUParamList*)fSimuParam->FindRespFunParams(dType)) ) {AliFatal(Form("Response for DetType#%d is not found in SimuParams",dType)); exit(1);}
877   }
878   // delete non needed segmentations
879   for (int i=fNLayers;i--;) arrSeg.Remove(fSegModelLr[i]);
880   arrSeg.Delete();
881   //
882   InitArrays();
883   //
884   fSimInitDone = kTRUE;
885   //
886 }