]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMFT.cxx
Modifications to MFT classes to manage embedding
[u/mrichter/AliRoot.git] / MFT / AliMFT.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 //====================================================================================================================================================
17 //
18 //      Main class of the ALICE Muon Forward Tracker
19 //
20 //      Contact author: antonio.uras@cern.ch
21 //
22 //====================================================================================================================================================
23
24 #include "AliLog.h"
25 #include "TFile.h"  
26 #include "TGeoManager.h"    
27 #include "TGeoVolume.h"
28 #include "TGeoMatrix.h"
29 #include "TVirtualMC.h"
30 #include "TClonesArray.h"
31 #include "TGeoGlobalMagField.h"
32 #include "AliRun.h"
33 #include "AliLoader.h"
34 #include "AliDetector.h"
35 #include "AliMC.h"
36 #include "AliMagF.h"
37 #include "AliMFT.h"
38 #include "AliMFTHit.h"
39 #include "AliMFTDigit.h"
40 #include "AliMFTCluster.h"
41 #include "AliTrackReference.h"
42 #include "AliMFTSegmentation.h"
43 #include "AliMFTDigitizer.h"
44 #include "AliMFTPlane.h"
45 #include "TString.h"
46 #include "TObjArray.h"
47
48 ClassImp(AliMFT) 
49
50 //====================================================================================================================================================
51
52 AliMFT::AliMFT():
53   AliDetector(),
54   fVersion(1),
55   fNPlanes(0),
56   fNSlices(1),
57   fSDigitsPerPlane(0),
58   fDigitsPerPlane(0),
59   fRecPointsPerPlane(0),
60   fSideDigits(0),
61   fSegmentation(0),
62   fNameGeomFile(0),
63   fChargeDispersion(25.e-4),
64   fSingleStepForChargeDispersion(0),
65   fNStepForChargeDispersion(4),
66   fDensitySupportOverSi(0.15),
67   fFileNameForUnderyingEvent(0), 
68   fFileNameForPileUpEvents(0),
69   fNPileUpEvents(0), 
70   fUnderlyingEventID(-1)
71 {
72
73   // default constructor
74
75   for (Int_t iPileUp=0; iPileUp<AliMFTConstants::fNMaxPileUpEvents; iPileUp++) fPileUpEventsIDs[iPileUp] = -1;
76
77
78 }
79
80 //====================================================================================================================================================
81
82 AliMFT::AliMFT(const Char_t *name, const Char_t *title):
83   AliDetector(name, title),
84   fVersion(1),
85   fNPlanes(0),
86   fNSlices(1),
87   fSDigitsPerPlane(0),
88   fDigitsPerPlane(0),
89   fRecPointsPerPlane(0),
90   fSideDigits(0),
91   fSegmentation(0),
92   fNameGeomFile(0),
93   fChargeDispersion(25.e-4),
94   fSingleStepForChargeDispersion(0),
95   fNStepForChargeDispersion(4),
96   fDensitySupportOverSi(0.15),
97   fFileNameForUnderyingEvent(0), 
98   fFileNameForPileUpEvents(0),
99   fNPileUpEvents(0), 
100   fUnderlyingEventID(-1)
101 {
102
103   for (Int_t iPileUp=0; iPileUp<AliMFTConstants::fNMaxPileUpEvents; iPileUp++) fPileUpEventsIDs[iPileUp] = -1;
104
105   fNameGeomFile = "AliMFTGeometry.root";
106
107   SetGeometry();
108
109   Init();
110
111 }
112
113 //====================================================================================================================================================
114
115 AliMFT::AliMFT(const Char_t *name, const Char_t *title, Char_t *nameGeomFile):
116   AliDetector(name, title),
117   fVersion(1),
118   fNPlanes(0),
119   fNSlices(1),
120   fSDigitsPerPlane(0),
121   fDigitsPerPlane(0),
122   fRecPointsPerPlane(0),
123   fSideDigits(0),
124   fSegmentation(0),
125   fNameGeomFile(0),
126   fChargeDispersion(25.e-4),
127   fSingleStepForChargeDispersion(0),
128   fNStepForChargeDispersion(4),
129   fDensitySupportOverSi(0.15),
130   fFileNameForUnderyingEvent(0), 
131   fFileNameForPileUpEvents(0),
132   fNPileUpEvents(0), 
133   fUnderlyingEventID(-1)
134 {
135
136   for (Int_t iPileUp=0; iPileUp<AliMFTConstants::fNMaxPileUpEvents; iPileUp++) fPileUpEventsIDs[iPileUp] = -1;
137
138   fNameGeomFile = nameGeomFile;
139
140   SetGeometry();
141
142   Init();
143
144 }
145
146 //====================================================================================================================================================
147
148 AliMFT::~AliMFT() {
149
150   if (fSDigitsPerPlane)   { fSDigitsPerPlane->Delete();    delete fSDigitsPerPlane;   }
151   if (fDigitsPerPlane)    { fDigitsPerPlane->Delete();     delete fDigitsPerPlane;    }
152   if (fRecPointsPerPlane) { fRecPointsPerPlane->Delete();  delete fRecPointsPerPlane; }
153
154 }
155
156 //====================================================================================================================================================
157
158 void AliMFT::CreateMaterials() {
159
160   // Definition of MFT materials  - to be updated to the most recent values
161
162   AliInfo("Start MFT materials");
163
164   // data from PDG booklet 2002                 density [gr/cm^3]     rad len [cm]           abs len [cm]    
165   Float_t   aSi = 28.085 ,    zSi   = 14. ,     dSi      =  2.329 ,   radSi   =  21.82/dSi , absSi   = 108.4/dSi  ;    // Silicon
166   Float_t   aCarb = 12.01 ,   zCarb =  6. ,     dCarb    =  2.265 ,   radCarb =  18.8 ,      absCarb = 49.9       ;    // Carbon
167   Float_t   aAlu = 26.98 ,    zAlu  = 13. ,     dAlu     =  2.70  ,   radAlu  =  8.897 ,     absAlu  = 39.70      ;    // Aluminum
168
169   // Air mixture
170   const Int_t nAir = 4;   
171   Float_t   aAir[nAir] = {12, 14, 16, 36} ,  zAir[nAir] = {6, 7, 8, 18} ,   wAir[nAir]={0.000124, 0.755267, 0.231781, 0.012827} , dAir=0.00120479;   
172
173   // Water mixture
174   const Int_t nWater = 2;   
175   Float_t   aWater[nWater] = {1.00794, 15.9994} ,  zWater[nWater] = {1, 8} ,   wWater[nWater] = {0.111894, 0.888106} , dWater=1.;   
176
177   // SiO2 mixture
178   const Int_t nSiO2 = 2; 
179   Float_t   aSiO2[nSiO2] = {15.9994, 28.0855} ,   zSiO2[nSiO2] = {8., 14.} ,   wSiO2[nSiO2] = {0.532565, 0.467435} , dSiO2 = 2.20;  
180
181   // Inox mixture
182   const Int_t nInox = 9;
183   Float_t   aInox[nInox] = {12.0107, 54.9380, 28.0855, 30.9738, 32.0660, 58.6928, 51.9961, 95.9400, 55.8450} ;   
184   Float_t   zInox[nInox] = { 6,      25,      14,      15,      16,      28,      24,      42,      26     } ;   
185   Float_t   wInox[nInox] = {0.0003,  0.02,    0.01,    0.00045, 0.0003,  0.12,    0.17,    0.025,   0.65395} ;
186   Float_t   dInox = 8.03; 
187   
188   Int_t   matId  = 0;                        // tmp material id number
189   Int_t   unsens = 0, sens=1;                // sensitive or unsensitive medium
190   Int_t   itgfld = 3;                        // type of field intergration 0 no field -1 user in guswim 1 Runge Kutta 2 helix 3 const field along z
191   Float_t maxfld = 5.;                       // max field value
192
193   Float_t tmaxfd = -10.0;                    // max deflection angle due to magnetic field in one step
194   Float_t stemax =  0.001;                   // max step allowed [cm]
195   Float_t deemax = -0.2;                     // maximum fractional energy loss in one step 0<deemax<=1  
196   Float_t epsil  =  0.001;                   // tracking precision [cm]   
197   Float_t stmin  = -0.001;                   // minimum step due to continuous processes [cm] (negative value: choose it automatically)
198
199   Float_t tmaxfdSi =  0.1;                   // max deflection angle due to magnetic field in one step
200   Float_t stemaxSi =  5.0e-4;                // maximum step allowed [cm]
201   Float_t deemaxSi =  0.1;                   // maximum fractional energy loss in one step 0<deemax<=1
202   Float_t epsilSi  =  0.5e-4;                // tracking precision [cm]
203   Float_t stminSi  = -0.001;                 // minimum step due to continuous processes [cm] (negative value: choose it automatically)
204
205   Int_t   isxfld = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ();   // from CreateMaterials in STRUCT/AliPIPEv3.cxx
206   Float_t sxmgmx = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();     // from CreateMaterials in STRUCT/AliPIPEv3.cxx
207       
208   AliMixture(++matId,"Air", aAir, zAir, dAir, nAir, wAir); 
209   AliMedium(kAir,    "Air", matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
210   
211   AliMaterial(++matId, "Si", aSi, zSi, dSi, radSi, absSi);  
212   AliMedium(kSi,       "Si", matId, sens, isxfld, sxmgmx, tmaxfdSi, stemaxSi, deemaxSi, epsilSi, stminSi);
213
214   AliMaterial(++matId, "Readout", aSi, zSi, dSi, radSi, absSi);  
215   AliMedium(kReadout,  "Readout", matId, unsens, isxfld, sxmgmx, tmaxfdSi, stemaxSi, deemaxSi, epsilSi, stminSi);
216
217   AliMaterial(++matId, "Support", aSi, zSi, dSi*fDensitySupportOverSi, radSi/fDensitySupportOverSi, absSi/fDensitySupportOverSi);  
218   AliMedium(kSupport,  "Support", matId, unsens, isxfld, sxmgmx, tmaxfdSi, stemaxSi, deemaxSi, epsilSi, stminSi);
219   
220   AliMaterial(++matId, "Carbon", aCarb, zCarb, dCarb, radCarb, absCarb );
221   AliMedium(kCarbon,   "Carbon", matId, unsens, isxfld,  sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
222   
223   AliMaterial(++matId, "Alu", aAlu, zAlu, dAlu, radAlu, absAlu);
224   AliMedium(kAlu,      "Alu", matId, unsens, isxfld,  sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
225   
226   AliMixture(++matId, "Water", aWater, zWater, dWater, nWater, wWater);
227   AliMedium(kWater,   "Water", matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
228   
229   AliMixture(++matId, "SiO2", aSiO2, zSiO2, dSiO2, nSiO2, wSiO2);
230   AliMedium(kSiO2,    "SiO2", matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
231   
232   AliMixture(++matId, "Inox", aInox, zInox, dInox, nInox, wInox);
233   AliMedium(kInox,    "Inox", matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
234
235   AliInfo("End MFT materials");
236           
237 }
238
239 //====================================================================================================================================================
240
241 void AliMFT::CreateGeometry() {
242
243   // Creates detailed geometry simulation (currently GEANT volumes tree)        
244
245   AliInfo("Start MFT preliminary version building");
246   if(!TVirtualMC::GetMC()->IsRootGeometrySupported()) return;                
247   TGeoVolumeAssembly *vol = CreateVol();
248   AliInfo("TGeoVolumeAssembly created!");
249   gGeoManager->GetVolume("ALIC")->AddNode(vol,0);
250   AliInfo("Stop MFT preliminary version building");
251
252   if (fNStepForChargeDispersion) fSingleStepForChargeDispersion = fChargeDispersion/Double_t(fNStepForChargeDispersion);
253
254
255
256 //====================================================================================================================================================
257
258 void AliMFT::AddAlignableVolumes() {
259
260   // Create entries for alignable volumes associating the symbolic volume
261   // name with the corresponding volume path. Needs to be syncronized with
262   // eventual changes in the geometry.
263
264   TString sysName = "MFT";
265   TString volPath = "/ALIC_1/MFT_0";
266   
267   if (!gGeoManager->SetAlignableEntry(sysName.Data(),volPath.Data())) {
268     AliFatal(Form("Alignable entry %s not created. Volume path %s not valid", sysName.Data(), volPath.Data()));
269   }  
270
271 }
272
273 //====================================================================================================================================================
274
275 void AliMFT::StepManager() {
276
277   // Full Step Manager
278
279   AliDebug(2, Form("Entering StepManager: TVirtualMC::GetMC()->CurrentVolName() = %s", TVirtualMC::GetMC()->CurrentVolName()));
280
281   if (!fSegmentation) AliFatal("No segmentation available");    // DO WE HAVE A SEGMENTATION???
282
283   if (!(this->IsActive())) return;
284   if (!(TVirtualMC::GetMC()->TrackCharge())) return;
285
286   TString planeNumber   = TVirtualMC::GetMC()->CurrentVolName();
287   TString detElemNumber = TVirtualMC::GetMC()->CurrentVolName();
288   if (planeNumber.Contains("support")) return;
289   if (planeNumber.Contains("readout")) return;
290   planeNumber.Remove(0,9);
291   detElemNumber.Remove(0,19);
292   planeNumber.Remove(2);
293   detElemNumber.Remove(3);
294   Int_t detElemID = fSegmentation->GetDetElemGlobalID(planeNumber.Atoi(), detElemNumber.Atoi());
295
296   if (TVirtualMC::GetMC()->IsTrackExiting()) {
297     AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kMFT);
298   }
299
300   static TLorentzVector position, momentum;
301   static AliMFTHit hit;
302   
303   Int_t  status = 0;
304   
305   // Track status
306   if (TVirtualMC::GetMC()->IsTrackInside())      status +=  1;
307   if (TVirtualMC::GetMC()->IsTrackEntering())    status +=  2;
308   if (TVirtualMC::GetMC()->IsTrackExiting())     status +=  4;
309   if (TVirtualMC::GetMC()->IsTrackOut())         status +=  8;
310   if (TVirtualMC::GetMC()->IsTrackDisappeared()) status += 16;
311   if (TVirtualMC::GetMC()->IsTrackStop())        status += 32;
312   if (TVirtualMC::GetMC()->IsTrackAlive())       status += 64;
313
314   // ---------- Fill hit structure
315
316   hit.SetDetElemID(detElemID);
317   hit.SetPlane(planeNumber.Atoi());
318   hit.SetTrack(gAlice->GetMCApp()->GetCurrentTrackNumber());
319     
320   TVirtualMC::GetMC()->TrackPosition(position);
321   TVirtualMC::GetMC()->TrackMomentum(momentum);
322
323   AliDebug(1, Form("AliMFT::StepManager()->%s Hit #%06d (x=%f, y=%f, z=%f) belongs to track %02d\n", 
324                    TVirtualMC::GetMC()->CurrentVolName(), fNhits, position.X(), position.Y(), position.Z(), gAlice->GetMCApp()->GetCurrentTrackNumber())); 
325
326   hit.SetPosition(position);
327   hit.SetTOF(TVirtualMC::GetMC()->TrackTime());
328   hit.SetMomentum(momentum);
329   hit.SetStatus(status);
330   hit.SetEloss(TVirtualMC::GetMC()->Edep());
331   //  hit.SetShunt(GetIshunt());
332 //   if (TVirtualMC::GetMC()->IsTrackEntering()) {
333 //     hit.SetStartPosition(position);
334 //     hit.SetStartTime(TVirtualMC::GetMC()->TrackTime());
335 //     hit.SetStartStatus(status);
336 //     return; // don't save entering hit.
337 //   } 
338
339   // Fill hit structure with this new hit.
340   new ((*fHits)[fNhits++]) AliMFTHit(hit);
341
342   // Save old position... for next hit.
343 //   hit.SetStartPosition(position);
344 //   hit.SetStartTime(TVirtualMC::GetMC()->TrackTime());
345 //   hit.SetStartStatus(status);
346
347   return;
348
349 }
350
351 //====================================================================================================================================================
352
353 TGeoVolumeAssembly* AliMFT::CreateVol() {
354
355   // method to create the MFT Geometry (silicon circular planes)
356
357   if (!fSegmentation) CreateGeometry();
358   
359   TGeoVolumeAssembly *vol = new TGeoVolumeAssembly("MFT");
360   TGeoMedium *silicon = gGeoManager->GetMedium("MFT_Si");
361   TGeoMedium *readout = gGeoManager->GetMedium("MFT_Readout");
362   TGeoMedium *support = gGeoManager->GetMedium("MFT_Support");
363   TGeoMedium *carbon = gGeoManager->GetMedium("MFT_Carbon");
364   TGeoMedium *alu = gGeoManager->GetMedium("MFT_Alu");
365   TGeoMedium *water = gGeoManager->GetMedium("MFT_Water");
366   TGeoMedium *si02 = gGeoManager->GetMedium("MFT_SiO2");
367   TGeoMedium *inox = gGeoManager->GetMedium("MFT_Inox");
368
369   // ---- Cage & Services Description --------------------------------------------
370   // R. Tieulent - 17/01/2014 - Basic description for ITS/TPC matching studies
371   // -----------------------------------------------------------------------------
372
373   TGeoVolumeAssembly *cageNservices = new TGeoVolumeAssembly("MFT_cageNservices");
374   
375   // cage definition
376   Float_t cageDz   = 150./2.;
377   Float_t cageRMax = 49.5;
378   Float_t cageRMin = cageRMax - 0.120; // 0.64% of X0
379   
380   TGeoVolume *cage = gGeoManager->MakeTube("MFT_cage", carbon, cageRMin, cageRMax, cageDz);
381   cage->SetLineColor(kBlue);
382
383   cageNservices->AddNode(cage,1,new TGeoTranslation(0., 0., 0. ));
384
385   // Services definition
386   TGeoVolumeAssembly *services = new TGeoVolumeAssembly("MFT_services");
387
388   // Aluminum bus-Bar
389   Float_t busBarDz = 150.;
390   Float_t busBarThick = 0.1 ;
391   Float_t busBarWidth = 1.;
392   
393   TGeoVolume *aluBusBar = gGeoManager->MakeBox("MFT_busbar", alu, busBarWidth/2., busBarThick/2., busBarDz/2.);
394   aluBusBar->SetLineColor(kYellow);
395
396   Int_t nBusBar = 30;
397   Float_t dPhiBusBar = 2.*TMath::Pi() / nBusBar;
398   Float_t dShift = cageRMin - busBarThick;
399   
400   TGeoRotation *rot;
401
402   for (Int_t iBusBar=0; iBusBar<nBusBar; iBusBar++) {
403     Float_t phi =  dPhiBusBar*iBusBar;
404     Float_t xp = dShift*TMath::Cos(phi);
405     Float_t yp = dShift*TMath::Sin(phi);
406     rot = new TGeoRotation();
407     rot->RotateZ(phi*TMath::RadToDeg()+90.);
408     services->AddNode(aluBusBar, iBusBar+1, new TGeoCombiTrans(xp,yp,0,rot));
409   }
410   
411   // Cooling Services  definition
412   TGeoVolumeAssembly *cooling = new TGeoVolumeAssembly("MFT_cooling");
413   // Cooling Water
414   Float_t coolingDz = 150.;
415   Float_t coolingR  = 0.3 /2. ; // 3mm in diameter
416   Float_t coolingDR = 0.02 ;    // Thickness of the pipe 0.2mm
417   
418   TGeoVolume *coolingWater = gGeoManager->MakeTube("MFT_coolingWater", water, 0., coolingR, coolingDz/2.);
419   coolingWater->SetLineColor(kCyan);
420   cooling->AddNode(coolingWater, 1, new TGeoTranslation(0,0,0 ));
421
422   // Cooling Pipes
423   TGeoVolume *coolingPipes = gGeoManager->MakeTube("MFT_coolingPipes", inox, coolingR, coolingR+coolingDR, coolingDz/2.);
424   coolingPipes->SetLineColor(kGray);
425   cooling->AddNode(coolingPipes,1,new TGeoTranslation(0,0,0 ));
426   
427   Int_t nCooling = 18;
428   dShift = cageRMin - coolingR ;
429   Float_t phi0 = 0.02;
430   
431   for (Int_t iCooling=0; iCooling<nCooling; iCooling++) {
432     Float_t phi ;
433     if (iCooling<nCooling/2) phi = dPhiBusBar*(iCooling+3) + phi0;
434     else phi = dPhiBusBar*(iCooling+9) + phi0;
435     Float_t xp = dShift*TMath::Cos(phi);
436     Float_t yp = dShift*TMath::Sin(phi);
437     services->AddNode(cooling, iCooling+1, new TGeoTranslation(xp,yp,0 ));
438   }
439   
440   // Optical Fibers
441   Float_t fiberDz = 150.;
442   Float_t fiberRadius = 0.0125 /2. ; // 0.125mm in diameter
443   
444   TGeoVolume *fiber = gGeoManager->MakeTube("MFT_fiber", si02, 0., fiberRadius, fiberDz/2.);
445   fiber->SetLineColor(kCyan);
446
447   Int_t nFiber = 340;
448   dShift = cageRMin - 2*fiberRadius;
449   phi0 = 0.03;
450   
451   for (Int_t iFiber=0; iFiber<nFiber; iFiber++) {
452     Float_t phi = dPhiBusBar*(Int_t)(iFiber/11) - phi0-(iFiber%11)*2.*TMath::ATan(fiberRadius/dShift);
453     Float_t xp  = dShift*TMath::Cos(phi);
454     Float_t yp  = dShift*TMath::Sin(phi);
455     services->AddNode(fiber, iFiber+1, new TGeoTranslation(xp,yp,0 ));
456   }
457
458   cageNservices->AddNode(services, 1, new TGeoTranslation(0., 0., 0. ));
459   
460   vol->AddNode(cageNservices,1,new TGeoTranslation(0., 0., 0. ));
461   
462   // ------------------- Creating volumes for MFT planes --------------------------------
463
464   Double_t origin[3] = {0};
465
466   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
467
468     AliDebug(1, Form("Creating volumes for MFT plane %02d",iPlane));
469
470     AliMFTPlane *plane = fSegmentation->GetPlane(iPlane);
471
472     // --------- support element(s)
473     
474     origin[0] =  0.5*(plane->GetSupportElement(0)->GetAxis(0)->GetXmax() + plane->GetSupportElement(0)->GetAxis(0)->GetXmin());
475     origin[1] =  0.5*(plane->GetSupportElement(0)->GetAxis(1)->GetXmax() + plane->GetSupportElement(0)->GetAxis(1)->GetXmin());
476     origin[2] = -0.5*(plane->GetSupportElement(0)->GetAxis(2)->GetXmax() + plane->GetSupportElement(0)->GetAxis(2)->GetXmin());
477     TGeoVolume *supportElem = gGeoManager->MakeTube(Form("MFT_plane%02d_support", iPlane), support, 
478                                                     plane->GetRMinSupport(), 
479                                                     plane->GetRMaxSupport(),
480                                                     0.5*(plane->GetSupportElement(0)->GetAxis(2)->GetXmax() - 
481                                                          plane->GetSupportElement(0)->GetAxis(2)->GetXmin()) );
482     AliDebug(2, Form("Created vol %s", supportElem->GetName()));
483     supportElem->SetLineColor(kCyan-9);
484     vol -> AddNode(supportElem, 0, new TGeoTranslation(origin[0], origin[1], origin[2]));
485
486     AliDebug(1, "support elements created!");
487     
488     // --------- active elements
489
490     for (Int_t iActive=0; iActive<plane->GetNActiveElements(); iActive++) {
491       
492       Double_t dx = 0.5*TMath::Abs(plane->GetActiveElement(iActive)->GetAxis(0)->GetXmax() - plane->GetActiveElement(iActive)->GetAxis(0)->GetXmin());
493       Double_t dy = 0.5*TMath::Abs(plane->GetActiveElement(iActive)->GetAxis(1)->GetXmax() - plane->GetActiveElement(iActive)->GetAxis(1)->GetXmin());
494       Double_t dz = 0.5*TMath::Abs(plane->GetActiveElement(iActive)->GetAxis(2)->GetXmax() - plane->GetActiveElement(iActive)->GetAxis(2)->GetXmin());
495       dz /= Double_t(fNSlices);
496
497       origin[0] =  0.5*(plane->GetActiveElement(iActive)->GetAxis(0)->GetXmax() + plane->GetActiveElement(iActive)->GetAxis(0)->GetXmin());
498       origin[1] =  0.5*(plane->GetActiveElement(iActive)->GetAxis(1)->GetXmax() + plane->GetActiveElement(iActive)->GetAxis(1)->GetXmin());
499
500       for (Int_t iSlice=0; iSlice<fNSlices; iSlice++) {
501         origin[2] = -0.5*(plane->GetActiveElement(iActive)->GetAxis(2)->GetXmin() + 2*dz*(iSlice+1) + plane->GetActiveElement(iActive)->GetAxis(2)->GetXmin() + 2*dz*(iSlice) );
502         TGeoVolume *activeElem = gGeoManager->MakeBox(Form("MFT_plane%02d_active%03d_slice%02d", iPlane, iActive, iSlice), silicon, dx, dy, dz);
503         AliDebug(2, Form("Created vol %s", activeElem->GetName()));
504         activeElem->SetLineColor(kGreen);
505         vol -> AddNode(activeElem, 0, new TGeoTranslation(origin[0], origin[1], origin[2]));
506       }
507
508     }
509
510     AliDebug(1, "active elements created!");
511
512     // --------- readout elements
513
514     for (Int_t iReadout=0; iReadout<plane->GetNReadoutElements(); iReadout++) {
515       
516       Double_t dx = 0.5*TMath::Abs(plane->GetReadoutElement(iReadout)->GetAxis(0)->GetXmax() - plane->GetReadoutElement(iReadout)->GetAxis(0)->GetXmin());
517       Double_t dy = 0.5*TMath::Abs(plane->GetReadoutElement(iReadout)->GetAxis(1)->GetXmax() - plane->GetReadoutElement(iReadout)->GetAxis(1)->GetXmin());
518       Double_t dz = 0.5*TMath::Abs(plane->GetReadoutElement(iReadout)->GetAxis(2)->GetXmax() - plane->GetReadoutElement(iReadout)->GetAxis(2)->GetXmin());
519
520       origin[0] =  0.5*(plane->GetReadoutElement(iReadout)->GetAxis(0)->GetXmax() + plane->GetReadoutElement(iReadout)->GetAxis(0)->GetXmin());
521       origin[1] =  0.5*(plane->GetReadoutElement(iReadout)->GetAxis(1)->GetXmax() + plane->GetReadoutElement(iReadout)->GetAxis(1)->GetXmin());
522       origin[2] = -0.5*(plane->GetReadoutElement(iReadout)->GetAxis(2)->GetXmax() + plane->GetReadoutElement(iReadout)->GetAxis(2)->GetXmin());
523
524       TGeoVolume *readoutElem = gGeoManager->MakeBox(Form("MFT_plane%02d_readout%03d", iPlane, iReadout), readout, dx, dy, dz);
525       AliDebug(2, Form("Created vol %s", readoutElem->GetName()));
526       readoutElem->SetLineColor(kRed);
527       vol -> AddNode(readoutElem, 0, new TGeoTranslation(origin[0], origin[1], origin[2]));
528       
529     }
530
531     AliDebug(1, "readout elements created!");
532
533   }
534
535   return vol;
536
537 }
538
539 //====================================================================================================================================================
540
541 void AliMFT::Hits2SDigits(){
542   
543   // Interface method invoked from AliSimulation to create a list of sdigits corresponding to list of hits. Every hit generates one sdigit.
544
545   AliDebug(1,"Start Hits2SDigits.");
546   
547   if (!fSegmentation) CreateGeometry();
548  
549   if (!fLoader->TreeH()) fLoader->LoadHits();
550
551   if (!fLoader->TreeS()) {
552
553     for (Int_t iEvt=0;iEvt<fLoader->GetRunLoader()->GetNumberOfEvents(); iEvt++) {
554
555       fLoader->GetRunLoader()->GetEvent(iEvt);
556       fLoader->MakeTree("S");
557       MakeBranch("S");
558       SetTreeAddress();
559
560       AliDebug(1, Form("Event %03d: fLoader->TreeH()->GetEntries() = %2d", iEvt, Int_t(fLoader->TreeH()->GetEntries())));
561
562       for (Int_t iTrack=0; iTrack<fLoader->TreeH()->GetEntries(); iTrack++) {
563         fLoader->TreeH()->GetEntry(iTrack);     
564         Hits2SDigitsLocal(Hits(), GetSDigitsList(), iTrack);    // convert these hits to a list of sdigits  
565       }
566       
567       fLoader->TreeS()->Fill();
568       fLoader->WriteSDigits("OVERWRITE");
569       ResetSDigits();
570
571     }
572   }
573
574   fLoader->UnloadHits();
575   fLoader->UnloadSDigits();
576
577   AliDebug(1,"Stop Hits2SDigits.");
578   
579 }
580
581 //====================================================================================================================================================
582
583 void AliMFT::Hits2SDigitsLocal(TClonesArray *hits, const TObjArray *pSDig, Int_t track) {
584
585   //  Add sdigits of these hits to the list
586   
587   AliDebug(1, "Entering Hits2SDigitsLocal");
588   
589   if (!fSegmentation) CreateGeometry();
590   
591   TClonesArray *pSDigList[fNMaxPlanes];
592   for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) pSDigList[iPlane] = NULL; 
593   for (Int_t iPlane=0; iPlane<fNPlanes;    iPlane++) { 
594     pSDigList[iPlane] = (TClonesArray*) (*pSDig)[iPlane];
595     AliDebug(1,Form("Entries of pSDigList %3d; plane: %02d,",pSDigList[iPlane]->GetEntries(),iPlane));
596     if (!track && pSDigList[iPlane]->GetEntries()!=0) AliErrorClass("Some of sdigits lists is not empty");
597   }
598   
599   for (Int_t iHit=0; iHit<hits->GetEntries(); iHit++) {
600
601     AliMFTHit *hit = (AliMFTHit*) hits->At(iHit);
602
603     // Creating "main digit"
604
605     AliMFTDigit *mainSDigit = new AliMFTDigit();
606     mainSDigit->SetEloss(hit->GetEloss());
607     mainSDigit->SetDetElemID(hit->GetDetElemID());
608     mainSDigit->SetPlane(hit->GetPlane());
609     mainSDigit->AddMCLabel(hit->GetTrack()); 
610
611     Int_t xPixel = -1;
612     Int_t yPixel = -1;
613     if (fSegmentation->Hit2PixelID(hit->X(), hit->Y(), mainSDigit->GetDetElemID(), xPixel, yPixel)) {
614       mainSDigit->SetPixID(xPixel, yPixel, 0);
615       mainSDigit->SetPixWidth(fSegmentation->GetPixelSizeX(mainSDigit->GetDetElemID()), 
616                           fSegmentation->GetPixelSizeY(mainSDigit->GetDetElemID()),
617                           fSegmentation->GetPixelSizeZ(mainSDigit->GetDetElemID()));  
618       mainSDigit->SetPixCenter(fSegmentation->GetPixelCenterX(mainSDigit->GetDetElemID(), xPixel), 
619                            fSegmentation->GetPixelCenterY(mainSDigit->GetDetElemID(), yPixel),
620                            fSegmentation->GetPixelCenterZ(mainSDigit->GetDetElemID(), 0));
621       new ((*fSideDigits)[fSideDigits->GetEntries()]) AliMFTDigit(*mainSDigit);
622       AliDebug(1, Form("Created new sdigit (%f, %f, %f) from hit (%f, %f, %f)",
623                        mainSDigit->GetPixelCenterX(), mainSDigit->GetPixelCenterY(), mainSDigit->GetPixelCenterZ(), hit->X(), hit->Y(), hit->Z()));
624     }
625
626     // creating "side digits" to simulate the effect of charge dispersion
627
628     Double_t pi4 = TMath::Pi()/4.;
629     for (Int_t iStep=0; iStep<fNStepForChargeDispersion; iStep++) {
630       Double_t shift = (iStep+1) * fSingleStepForChargeDispersion;
631       for (Int_t iAngle=0; iAngle<8; iAngle++) {
632         Double_t shiftX = shift*TMath::Cos(iAngle*pi4);
633         Double_t shiftY = shift*TMath::Sin(iAngle*pi4);
634         if (fSegmentation->Hit2PixelID(hit->X()+shiftX, hit->Y()+shiftY, hit->GetDetElemID(), xPixel, yPixel)) {
635           Bool_t digitExists = kFALSE;
636           for (Int_t iSideDigit=0; iSideDigit<fSideDigits->GetEntries(); iSideDigit++) {
637             if (xPixel==((AliMFTDigit*) fSideDigits->At(iSideDigit))->GetPixelX() && 
638                 yPixel==((AliMFTDigit*) fSideDigits->At(iSideDigit))->GetPixelY()) {
639               digitExists = kTRUE;
640               break;
641             }
642           }
643           if (!digitExists) {
644             AliMFTDigit *sideSDigit = new AliMFTDigit();
645             sideSDigit->SetEloss(0.);
646             sideSDigit->SetDetElemID(hit->GetDetElemID());
647             sideSDigit->SetPlane(hit->GetPlane());
648             sideSDigit->AddMCLabel(hit->GetTrack());
649             sideSDigit->SetPixID(xPixel, yPixel, 0);
650             sideSDigit->SetPixWidth(fSegmentation->GetPixelSizeX(sideSDigit->GetDetElemID()), 
651                                     fSegmentation->GetPixelSizeY(sideSDigit->GetDetElemID()),
652                                     fSegmentation->GetPixelSizeZ(sideSDigit->GetDetElemID()));  
653             sideSDigit->SetPixCenter(fSegmentation->GetPixelCenterX(sideSDigit->GetDetElemID(), xPixel), 
654                                      fSegmentation->GetPixelCenterY(sideSDigit->GetDetElemID(), yPixel),
655                                      fSegmentation->GetPixelCenterZ(sideSDigit->GetDetElemID(), 0)); 
656             new ((*fSideDigits)[fSideDigits->GetEntries()]) AliMFTDigit(*sideSDigit);
657           }
658         }
659       }
660     }
661     
662     // ------------ In case we should simulate a rectangular pattern of pixel...
663     
664     if (fSegmentation->GetPlane(mainSDigit->GetPlane())->HasPixelRectangularPatternAlongY()) {
665       for (Int_t iSDigit=0; iSDigit<fSideDigits->GetEntries(); iSDigit++) {
666         AliMFTDigit *mySDig = (AliMFTDigit*) (fSideDigits->At(iSDigit));
667         if (mySDig->GetPixelX()%2 == mySDig->GetPixelY()%2) {   // both pair or both odd
668           xPixel = mySDig->GetPixelX();
669           yPixel = mySDig->GetPixelY()+1;
670           if (fSegmentation->DoesPixelExist(mySDig->GetDetElemID(), xPixel, yPixel)) {
671             AliMFTDigit *newSDigit = new AliMFTDigit();
672             newSDigit->SetEloss(0.);
673             newSDigit->SetDetElemID(mySDig->GetDetElemID());
674             newSDigit->SetPlane(mySDig->GetDetElemID());
675             newSDigit->SetPixID(xPixel, yPixel, 0);
676             newSDigit->SetPixWidth(fSegmentation->GetPixelSizeX(newSDigit->GetDetElemID()), 
677                                    fSegmentation->GetPixelSizeY(newSDigit->GetDetElemID()),
678                                    fSegmentation->GetPixelSizeZ(newSDigit->GetDetElemID()));  
679             newSDigit->SetPixCenter(fSegmentation->GetPixelCenterX(newSDigit->GetDetElemID(), xPixel), 
680                                     fSegmentation->GetPixelCenterY(newSDigit->GetDetElemID(), yPixel),
681                                     fSegmentation->GetPixelCenterZ(newSDigit->GetDetElemID(), 0)); 
682             new ((*fSideDigits)[fSideDigits->GetEntries()]) AliMFTDigit(*newSDigit);
683           }
684         }
685         else {   // pair-odd
686           xPixel = mySDig->GetPixelX();
687           yPixel = mySDig->GetPixelY()-1;
688           if (fSegmentation->DoesPixelExist(mySDig->GetDetElemID(), xPixel, yPixel)) {
689             AliMFTDigit *newSDigit = new AliMFTDigit();
690             newSDigit->SetEloss(0.);
691             newSDigit->SetDetElemID(mySDig->GetDetElemID());
692             newSDigit->SetPlane(mySDig->GetPlane());
693             newSDigit->SetPixID(xPixel, yPixel, 0);
694             newSDigit->SetPixWidth(fSegmentation->GetPixelSizeX(newSDigit->GetDetElemID()), 
695                                    fSegmentation->GetPixelSizeY(newSDigit->GetDetElemID()),
696                                    fSegmentation->GetPixelSizeZ(newSDigit->GetDetElemID()));  
697             newSDigit->SetPixCenter(fSegmentation->GetPixelCenterX(newSDigit->GetDetElemID(), xPixel), 
698                                     fSegmentation->GetPixelCenterY(newSDigit->GetDetElemID(), yPixel),
699                                     fSegmentation->GetPixelCenterZ(newSDigit->GetDetElemID(), 0)); 
700             new ((*fSideDigits)[fSideDigits->GetEntries()]) AliMFTDigit(*newSDigit);
701           }
702         }
703       }
704     }
705
706     // -------- checking which pixels switched on have their diode actually within the charge dispersion radius
707
708     for (Int_t iSDigit=0; iSDigit<fSideDigits->GetEntries(); iSDigit++) {
709       AliMFTDigit *mySDig = (AliMFTDigit*) (fSideDigits->At(iSDigit));
710       Double_t distance = TMath::Sqrt(TMath::Power(mySDig->GetPixelCenterX()-hit->X(),2) + TMath::Power(mySDig->GetPixelCenterY()-hit->Y(),2));
711       if (fSegmentation->GetPlane(mySDig->GetPlane())->HasPixelRectangularPatternAlongY()) {
712         if (mySDig->GetPixelX()%2 == mySDig->GetPixelY()%2) {  // both pair or both odd
713           if (distance<fChargeDispersion) {
714             AliDebug(1, Form("Created new sdigit (%f, %f, %f) from hit (%f, %f, %f)",
715                              mySDig->GetPixelCenterX(), mySDig->GetPixelCenterY(), mySDig->GetPixelCenterZ(), hit->X(), hit->Y(), hit->Z()));
716             new ((*pSDigList[mySDig->GetPlane()])[pSDigList[mySDig->GetPlane()]->GetEntries()]) AliMFTDigit(*mySDig);
717             xPixel = mySDig->GetPixelX();
718             yPixel = mySDig->GetPixelY()+1;
719             if (fSegmentation->DoesPixelExist(mySDig->GetDetElemID(), xPixel, yPixel)) {
720               AliMFTDigit *newSDigit = new AliMFTDigit();
721               newSDigit->SetEloss(0.);
722               newSDigit->SetDetElemID(mySDig->GetDetElemID());
723               newSDigit->SetPlane(mySDig->GetPlane());
724               newSDigit->SetPixID(xPixel, yPixel, 0);
725               newSDigit->SetPixWidth(fSegmentation->GetPixelSizeX(newSDigit->GetDetElemID()), 
726                                      fSegmentation->GetPixelSizeY(newSDigit->GetDetElemID()),
727                                      fSegmentation->GetPixelSizeZ(newSDigit->GetDetElemID()));  
728               newSDigit->SetPixCenter(fSegmentation->GetPixelCenterX(newSDigit->GetDetElemID(), xPixel), 
729                                       fSegmentation->GetPixelCenterY(newSDigit->GetDetElemID(), yPixel),
730                                       fSegmentation->GetPixelCenterZ(newSDigit->GetDetElemID(), 0));
731               AliDebug(1, Form("Created new sdigit (%f, %f, %f) from hit (%f, %f, %f)",
732                                newSDigit->GetPixelCenterX(), newSDigit->GetPixelCenterY(), newSDigit->GetPixelCenterZ(), hit->X(), hit->Y(), hit->Z()));
733               new ((*pSDigList[newSDigit->GetPlane()])[pSDigList[newSDigit->GetPlane()]->GetEntries()]) AliMFTDigit(*newSDigit);
734             }
735           }
736         }
737       }
738       else {
739         if (distance<fChargeDispersion) {
740           AliDebug(1, Form("Created new sdigit (%f, %f, %f) from hit (%f, %f, %f)",
741                            mySDig->GetPixelCenterX(), mySDig->GetPixelCenterY(), mySDig->GetPixelCenterZ(), hit->X(), hit->Y(), hit->Z()));
742           new ((*pSDigList[mySDig->GetPlane()])[pSDigList[mySDig->GetPlane()]->GetEntries()]) AliMFTDigit(*mySDig);
743         }
744       }
745     }
746
747     fSideDigits->Delete(); 
748
749   }
750
751   AliDebug(1,"Exiting Hits2SDigitsLocal");
752
753 }
754
755 //====================================================================================================================================================
756
757 void AliMFT::MakeBranch(Option_t *option) {
758
759   // Create Tree branches 
760   AliDebug(1, Form("Start with option= %s.",option));
761   
762   const Int_t kBufSize = 4000;
763   
764   const Char_t *cH = strstr(option,"H");
765   const Char_t *cD = strstr(option,"D");
766   const Char_t *cS = strstr(option,"S");
767
768   if (cH && fLoader->TreeH()) {
769     CreateHits();
770     MakeBranchInTree(fLoader->TreeH(), "MFT", &fHits, kBufSize, 0);   
771   }
772
773   if (cS && fLoader->TreeS()) {
774     CreateSDigits();
775     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) MakeBranchInTree(fLoader->TreeS(), 
776                                                                     Form("Plane_%02d",iPlane), 
777                                                                     &((*fSDigitsPerPlane)[iPlane]), 
778                                                                     kBufSize, 0);
779   }
780   
781   if (cD && fLoader->TreeD()) {
782     CreateDigits();
783     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) MakeBranchInTree(fLoader->TreeD(), 
784                                                                     Form("Plane_%02d",iPlane), 
785                                                                     &((*fDigitsPerPlane)[iPlane]),
786                                                                     kBufSize, 0);
787   }
788
789   AliDebug(1,"Stop.");
790
791 }
792
793 //====================================================================================================================================================
794
795 void AliMFT::SetTreeAddress() {
796
797   AliDebug(1, "AliMFT::SetTreeAddress()");
798
799   //Set branch address for the Hits and Digits Tree.
800   AliDebug(1, "Start.");
801
802   AliDebug(1, Form("AliMFT::SetTreeAddress Hits  fLoader->TreeH() = %p\n", fLoader->TreeH()));
803   if (fLoader->TreeH() && fLoader->TreeH()->GetBranch("MFT")) {
804     CreateHits();
805     fLoader->TreeH()->SetBranchAddress("MFT", &fHits);
806   }
807     
808   AliDebug(1, Form("AliMFT::SetTreeAddress SDigits  fLoader->TreeS() = %p\n", fLoader->TreeS()));
809   if (fLoader->TreeS() && fLoader->TreeS()->GetBranch("Plane_00")) {
810     CreateSDigits();
811     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
812       fLoader->TreeS()->SetBranchAddress(Form("Plane_%02d",iPlane), &((*fSDigitsPerPlane)[iPlane]));
813     }
814   }
815     
816   AliDebug(1, Form("AliMFT::SetTreeAddress Digits  fLoader->TreeD() = %p\n", fLoader->TreeD()));
817   if (fLoader->TreeD() && fLoader->TreeD()->GetBranch("Plane_00")) {
818     CreateDigits(); 
819     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
820       fLoader->TreeD()->SetBranchAddress(Form("Plane_%02d",iPlane), &((*fDigitsPerPlane)[iPlane]));
821     }
822   }
823
824   AliDebug(1, Form("AliMFT::SetTreeAddress RecPoints  fLoader->TreeR() = %p\n", fLoader->TreeR()));
825   if (fLoader->TreeR() && fLoader->TreeR()->GetBranch("Plane_00")) {
826     CreateRecPoints(); 
827     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
828       fLoader->TreeR()->SetBranchAddress(Form("Plane_%02d",iPlane), &((*fRecPointsPerPlane)[iPlane]));
829     }
830   }
831
832   AliDebug(1,"Stop.");
833
834 }
835
836 //====================================================================================================================================================
837
838 void AliMFT::SetGeometry() {
839
840   AliInfo("AliMFT::SetGeometry\n");
841
842   fSegmentation = new AliMFTSegmentation(fNameGeomFile.Data());
843
844   fNPlanes = fSegmentation->GetNPlanes();
845
846 }
847
848 //====================================================================================================================================================
849
850 void AliMFT::CreateHits() { 
851
852   // create array of hits
853
854   AliDebug(1, "AliMFT::CreateHits()");
855
856   if (fHits) return;    
857   fHits = new TClonesArray("AliMFTHit");  
858
859 }
860
861 //====================================================================================================================================================
862
863 void AliMFT::CreateSDigits() { 
864  
865   // create sdigits list
866
867   AliDebug(1, "AliMFT::CreateSDigits()");
868  
869   if (fSDigitsPerPlane) return; 
870   fSDigitsPerPlane = new TObjArray(fNPlanes); 
871   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) fSDigitsPerPlane->AddAt(new TClonesArray("AliMFTDigit"), iPlane);
872
873   fSideDigits = new TClonesArray("AliMFTDigit");
874
875 }
876
877 //====================================================================================================================================================
878
879 void AliMFT::CreateDigits() {
880
881   // create digits list
882
883   AliDebug(1, "AliMFT::CreateDigits()");
884
885   if (fDigitsPerPlane) return; 
886   fDigitsPerPlane = new TObjArray(fNPlanes);
887   for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) fDigitsPerPlane->AddAt(new TClonesArray("AliMFTDigit"), iPlane);
888
889 }
890
891 //====================================================================================================================================================
892
893 void AliMFT::CreateRecPoints() {
894
895   // create recPoints list
896
897   AliDebug(1, "AliMFT::CreateRecPoints()");
898
899   if (fRecPointsPerPlane) return; 
900   fRecPointsPerPlane = new TObjArray(fNPlanes);
901   for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) fRecPointsPerPlane->AddAt(new TClonesArray("AliMFTCluster"), iPlane);
902
903 }
904
905 //====================================================================================================================================================