]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMFT.cxx
Coordinates of SPD vertex
[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.036),
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.036),
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.036),
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");
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("active"))) return;
289   planeNumber.Remove(0,9);
290   detElemNumber.Remove(0,18);
291   planeNumber.Remove(2);
292   detElemNumber.Remove(3);
293   Int_t detElemID = fSegmentation->GetDetElemGlobalID(planeNumber.Atoi(), detElemNumber.Atoi());
294
295   if (TVirtualMC::GetMC()->IsTrackExiting()) {
296     AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kMFT);
297   }
298
299   static TLorentzVector position, momentum;
300   static AliMFTHit hit;
301   
302   Int_t  status = 0;
303   
304   // Track status
305   if (TVirtualMC::GetMC()->IsTrackInside())      status +=  1;
306   if (TVirtualMC::GetMC()->IsTrackEntering())    status +=  2;
307   if (TVirtualMC::GetMC()->IsTrackExiting())     status +=  4;
308   if (TVirtualMC::GetMC()->IsTrackOut())         status +=  8;
309   if (TVirtualMC::GetMC()->IsTrackDisappeared()) status += 16;
310   if (TVirtualMC::GetMC()->IsTrackStop())        status += 32;
311   if (TVirtualMC::GetMC()->IsTrackAlive())       status += 64;
312
313   // ---------- Fill hit structure
314
315   hit.SetDetElemID(detElemID);
316   hit.SetPlane(planeNumber.Atoi());
317   hit.SetTrack(gAlice->GetMCApp()->GetCurrentTrackNumber());
318     
319   TVirtualMC::GetMC()->TrackPosition(position);
320   TVirtualMC::GetMC()->TrackMomentum(momentum);
321
322   AliDebug(1, Form("AliMFT::StepManager()->%s Hit #%06d (x=%f, y=%f, z=%f) belongs to track %02d\n", 
323                    TVirtualMC::GetMC()->CurrentVolName(), fNhits, position.X(), position.Y(), position.Z(), gAlice->GetMCApp()->GetCurrentTrackNumber())); 
324
325   hit.SetPosition(position);
326   hit.SetTOF(TVirtualMC::GetMC()->TrackTime());
327   hit.SetMomentum(momentum);
328   hit.SetStatus(status);
329   hit.SetEloss(TVirtualMC::GetMC()->Edep());
330   //  hit.SetShunt(GetIshunt());
331 //   if (TVirtualMC::GetMC()->IsTrackEntering()) {
332 //     hit.SetStartPosition(position);
333 //     hit.SetStartTime(TVirtualMC::GetMC()->TrackTime());
334 //     hit.SetStartStatus(status);
335 //     return; // don't save entering hit.
336 //   } 
337
338   // Fill hit structure with this new hit.
339   new ((*fHits)[fNhits++]) AliMFTHit(hit);
340
341   // Save old position... for next hit.
342 //   hit.SetStartPosition(position);
343 //   hit.SetStartTime(TVirtualMC::GetMC()->TrackTime());
344 //   hit.SetStartStatus(status);
345
346   return;
347
348 }
349
350 //====================================================================================================================================================
351
352 TGeoVolumeAssembly* AliMFT::CreateVol() {
353
354   // method to create the MFT Geometry (silicon circular planes)
355
356   if (!fSegmentation) CreateGeometry();
357   
358   TGeoVolumeAssembly *vol = new TGeoVolumeAssembly("MFT");
359   TGeoMedium *silicon = gGeoManager->GetMedium("MFT_Si");
360   TGeoMedium *readout = gGeoManager->GetMedium("MFT_Readout");
361   TGeoMedium *support = gGeoManager->GetMedium("MFT_Support");
362   TGeoMedium *carbon  = gGeoManager->GetMedium("MFT_Carbon");
363   TGeoMedium *alu     = gGeoManager->GetMedium("MFT_Alu");
364   TGeoMedium *water   = gGeoManager->GetMedium("MFT_Water");
365   TGeoMedium *si02    = gGeoManager->GetMedium("MFT_SiO2");
366   TGeoMedium *inox    = gGeoManager->GetMedium("MFT_Inox");
367
368   // ---- Cage & Services Description --------------------------------------------
369   // R. Tieulent - 17/01/2014 - Basic description for ITS/TPC matching studies
370   // -----------------------------------------------------------------------------
371
372   TGeoVolumeAssembly *cageNservices = new TGeoVolumeAssembly("MFT_cageNservices");
373   
374   // cage definition
375   Float_t cageDz   = 150./2.;
376   Float_t cageRMax = 49.5;
377   Float_t cageRMin = cageRMax - 0.120; // 0.64% of X0
378   
379   TGeoVolume *cage = gGeoManager->MakeTube("MFT_cage", carbon, cageRMin, cageRMax, cageDz);
380   cage->SetLineColor(kBlue);
381
382   cageNservices->AddNode(cage,1,new TGeoTranslation(0., 0., 0. ));
383
384   // Services definition
385   TGeoVolumeAssembly *services = new TGeoVolumeAssembly("MFT_services");
386
387   // Aluminum bus-Bar
388   Float_t busBarDz = 150.;
389   Float_t busBarThick = 0.1 ;
390   Float_t busBarWidth = 1.;
391   
392   TGeoVolume *aluBusBar = gGeoManager->MakeBox("MFT_busbar", alu, busBarWidth/2., busBarThick/2., busBarDz/2.);
393   aluBusBar->SetLineColor(kYellow);
394
395   Int_t nBusBar = 30;
396   Float_t dPhiBusBar = 2.*TMath::Pi() / nBusBar;
397   Float_t dShift = cageRMin - busBarThick;
398   
399   TGeoRotation *rot;
400
401   for (Int_t iBusBar=0; iBusBar<nBusBar; iBusBar++) {
402     Float_t phi =  dPhiBusBar*iBusBar;
403     Float_t xp = dShift*TMath::Cos(phi);
404     Float_t yp = dShift*TMath::Sin(phi);
405     rot = new TGeoRotation();
406     rot->RotateZ(phi*TMath::RadToDeg()+90.);
407     services->AddNode(aluBusBar, iBusBar+1, new TGeoCombiTrans(xp,yp,0,rot));
408   }
409   
410   // Cooling Services  definition
411   TGeoVolumeAssembly *cooling = new TGeoVolumeAssembly("MFT_cooling");
412   // Cooling Water
413   Float_t coolingDz = 150.;
414   Float_t coolingR  = 0.3 /2. ; // 3mm in diameter
415   Float_t coolingDR = 0.02 ;    // Thickness of the pipe 0.2mm
416   
417   TGeoVolume *coolingWater = gGeoManager->MakeTube("MFT_coolingWater", water, 0., coolingR, coolingDz/2.);
418   coolingWater->SetLineColor(kCyan);
419   cooling->AddNode(coolingWater, 1, new TGeoTranslation(0,0,0 ));
420
421   // Cooling Pipes
422   TGeoVolume *coolingPipes = gGeoManager->MakeTube("MFT_coolingPipes", inox, coolingR, coolingR+coolingDR, coolingDz/2.);
423   coolingPipes->SetLineColor(kGray);
424   cooling->AddNode(coolingPipes,1,new TGeoTranslation(0,0,0 ));
425   
426   Int_t nCooling = 18;
427   dShift = cageRMin - coolingR ;
428   Float_t phi0 = 0.02;
429   
430   for (Int_t iCooling=0; iCooling<nCooling; iCooling++) {
431     Float_t phi ;
432     if (iCooling<nCooling/2) phi = dPhiBusBar*(iCooling+3) + phi0;
433     else phi = dPhiBusBar*(iCooling+9) + phi0;
434     Float_t xp = dShift*TMath::Cos(phi);
435     Float_t yp = dShift*TMath::Sin(phi);
436     services->AddNode(cooling, iCooling+1, new TGeoTranslation(xp,yp,0 ));
437   }
438   
439   // Optical Fibers
440   Float_t fiberDz = 150.;
441   Float_t fiberRadius = 0.0125 /2. ; // 0.125mm in diameter
442   
443   TGeoVolume *fiber = gGeoManager->MakeTube("MFT_fiber", si02, 0., fiberRadius, fiberDz/2.);
444   fiber->SetLineColor(kCyan);
445
446   Int_t nFiber = 340;
447   dShift = cageRMin - 2*fiberRadius;
448   phi0 = 0.03;
449   
450   for (Int_t iFiber=0; iFiber<nFiber; iFiber++) {
451     Float_t phi = dPhiBusBar*(Int_t)(iFiber/11) - phi0-(iFiber%11)*2.*TMath::ATan(fiberRadius/dShift);
452     Float_t xp  = dShift*TMath::Cos(phi);
453     Float_t yp  = dShift*TMath::Sin(phi);
454     services->AddNode(fiber, iFiber+1, new TGeoTranslation(xp,yp,0 ));
455   }
456
457   cageNservices->AddNode(services, 1, new TGeoTranslation(0., 0., 0. ));
458   
459   vol->AddNode(cageNservices,1,new TGeoTranslation(0., 0., 0. ));
460   
461   // ------------------- Creating volumes for MFT planes --------------------------------
462
463   Double_t origin[3] = {0};
464
465   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
466
467     AliDebug(1, Form("Creating volumes for MFT plane %02d",iPlane));
468
469     AliMFTPlane *plane = fSegmentation->GetPlane(iPlane);
470
471     // --------- support element(s)
472     
473     origin[0] =  0.5*(plane->GetSupportElement(0)->GetAxis(0)->GetXmax() + plane->GetSupportElement(0)->GetAxis(0)->GetXmin());
474     origin[1] =  0.5*(plane->GetSupportElement(0)->GetAxis(1)->GetXmax() + plane->GetSupportElement(0)->GetAxis(1)->GetXmin());
475     origin[2] = -0.5*(plane->GetSupportElement(0)->GetAxis(2)->GetXmax() + plane->GetSupportElement(0)->GetAxis(2)->GetXmin());
476     TGeoVolume *supportElem = gGeoManager->MakeTube(Form("MFT_plane%02d_support", iPlane), support, 
477                                                     plane->GetRMinSupport(), 
478                                                     plane->GetRMaxSupport(),
479                                                     0.5*(plane->GetSupportElement(0)->GetAxis(2)->GetXmax() - 
480                                                          plane->GetSupportElement(0)->GetAxis(2)->GetXmin()) );
481     AliDebug(2, Form("Created vol %s", supportElem->GetName()));
482     supportElem->SetLineColor(kCyan-9);
483     vol -> AddNode(supportElem, 0, new TGeoTranslation(origin[0], origin[1], origin[2]));
484
485     AliDebug(1, "support elements created!");
486     
487     // --------- active elements
488
489     for (Int_t iActive=0; iActive<plane->GetNActiveElements(); iActive++) {
490       
491       Double_t dx = 0.5*TMath::Abs(plane->GetActiveElement(iActive)->GetAxis(0)->GetXmax() - plane->GetActiveElement(iActive)->GetAxis(0)->GetXmin());
492       Double_t dy = 0.5*TMath::Abs(plane->GetActiveElement(iActive)->GetAxis(1)->GetXmax() - plane->GetActiveElement(iActive)->GetAxis(1)->GetXmin());
493       Double_t dz = 0.5*TMath::Abs(plane->GetActiveElement(iActive)->GetAxis(2)->GetXmax() - plane->GetActiveElement(iActive)->GetAxis(2)->GetXmin());
494       dz /= Double_t(fNSlices);
495
496       origin[0] =  0.5*(plane->GetActiveElement(iActive)->GetAxis(0)->GetXmax() + plane->GetActiveElement(iActive)->GetAxis(0)->GetXmin());
497       origin[1] =  0.5*(plane->GetActiveElement(iActive)->GetAxis(1)->GetXmax() + plane->GetActiveElement(iActive)->GetAxis(1)->GetXmin());
498
499       for (Int_t iSlice=0; iSlice<fNSlices; iSlice++) {
500         origin[2] = -0.5*(plane->GetActiveElement(iActive)->GetAxis(2)->GetXmin() + 2*dz*(iSlice+1) + plane->GetActiveElement(iActive)->GetAxis(2)->GetXmin() + 2*dz*(iSlice) );
501         TGeoVolume *activeElem = gGeoManager->MakeBox(Form("MFT_plane%02d_active%03d_slice%02d", iPlane, iActive, iSlice), silicon, dx, dy, dz);
502         AliDebug(2, Form("Created vol %s", activeElem->GetName()));
503         activeElem->SetLineColor(kGreen);
504         vol -> AddNode(activeElem, 0, new TGeoTranslation(origin[0], origin[1], origin[2]));
505       }
506
507     }
508
509     AliDebug(1, "active elements created!");
510
511     // --------- readout elements
512
513     for (Int_t iReadout=0; iReadout<plane->GetNReadoutElements(); iReadout++) {
514       
515       Double_t dx = 0.5*TMath::Abs(plane->GetReadoutElement(iReadout)->GetAxis(0)->GetXmax() - plane->GetReadoutElement(iReadout)->GetAxis(0)->GetXmin());
516       Double_t dy = 0.5*TMath::Abs(plane->GetReadoutElement(iReadout)->GetAxis(1)->GetXmax() - plane->GetReadoutElement(iReadout)->GetAxis(1)->GetXmin());
517       Double_t dz = 0.5*TMath::Abs(plane->GetReadoutElement(iReadout)->GetAxis(2)->GetXmax() - plane->GetReadoutElement(iReadout)->GetAxis(2)->GetXmin());
518
519       origin[0] =  0.5*(plane->GetReadoutElement(iReadout)->GetAxis(0)->GetXmax() + plane->GetReadoutElement(iReadout)->GetAxis(0)->GetXmin());
520       origin[1] =  0.5*(plane->GetReadoutElement(iReadout)->GetAxis(1)->GetXmax() + plane->GetReadoutElement(iReadout)->GetAxis(1)->GetXmin());
521       origin[2] = -0.5*(plane->GetReadoutElement(iReadout)->GetAxis(2)->GetXmax() + plane->GetReadoutElement(iReadout)->GetAxis(2)->GetXmin());
522
523       TGeoVolume *readoutElem = gGeoManager->MakeBox(Form("MFT_plane%02d_readout%03d", iPlane, iReadout), readout, dx, dy, dz);
524       AliDebug(2, Form("Created vol %s", readoutElem->GetName()));
525       readoutElem->SetLineColor(kRed);
526       vol -> AddNode(readoutElem, 0, new TGeoTranslation(origin[0], origin[1], origin[2]));
527       
528     }
529
530     AliDebug(1, "readout elements created!");
531
532   }
533
534   return vol;
535
536 }
537
538 //====================================================================================================================================================
539
540 void AliMFT::Hits2SDigits(){
541   
542   // Interface method invoked from AliSimulation to create a list of sdigits corresponding to list of hits. Every hit generates one sdigit.
543
544   AliDebug(1,"Start Hits2SDigits.");
545   
546   if (!fSegmentation) CreateGeometry();
547  
548   if (!fLoader->TreeH()) fLoader->LoadHits();
549
550   if (!fLoader->TreeS()) {
551
552     for (Int_t iEvt=0;iEvt<fLoader->GetRunLoader()->GetNumberOfEvents(); iEvt++) {
553
554       fLoader->GetRunLoader()->GetEvent(iEvt);
555       fLoader->MakeTree("S");
556       MakeBranch("S");
557       SetTreeAddress();
558
559       AliDebug(1, Form("Event %03d: fLoader->TreeH()->GetEntries() = %2d", iEvt, Int_t(fLoader->TreeH()->GetEntries())));
560
561       for (Int_t iTrack=0; iTrack<fLoader->TreeH()->GetEntries(); iTrack++) {
562         fLoader->TreeH()->GetEntry(iTrack);     
563         Hits2SDigitsLocal(Hits(), GetSDigitsList(), iTrack);    // convert these hits to a list of sdigits  
564       }
565       
566       fLoader->TreeS()->Fill();
567       fLoader->WriteSDigits("OVERWRITE");
568       ResetSDigits();
569
570     }
571   }
572
573   fLoader->UnloadHits();
574   fLoader->UnloadSDigits();
575
576   AliDebug(1,"Stop Hits2SDigits.");
577   
578 }
579
580 //====================================================================================================================================================
581
582 void AliMFT::Hits2SDigitsLocal(TClonesArray *hits, const TObjArray *pSDig, Int_t track) {
583
584   //  Add sdigits of these hits to the list
585   
586   AliDebug(1, "Entering Hits2SDigitsLocal");
587   
588   if (!fSegmentation) CreateGeometry();
589   
590   TClonesArray *pSDigList[fNMaxPlanes];
591   for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) pSDigList[iPlane] = NULL; 
592   for (Int_t iPlane=0; iPlane<fNPlanes;    iPlane++) { 
593     pSDigList[iPlane] = (TClonesArray*) (*pSDig)[iPlane];
594     AliDebug(1,Form("Entries of pSDigList %3d; plane: %02d,",pSDigList[iPlane]->GetEntries(),iPlane));
595     if (!track && pSDigList[iPlane]->GetEntries()!=0) AliErrorClass("Some of sdigits lists is not empty");
596   }
597   
598   for (Int_t iHit=0; iHit<hits->GetEntries(); iHit++) {
599
600     AliMFTHit *hit = (AliMFTHit*) hits->At(iHit);
601
602     // Creating "main digit"
603
604     AliMFTDigit *mainSDigit = new AliMFTDigit();
605     mainSDigit->SetEloss(hit->GetEloss());
606     mainSDigit->SetDetElemID(hit->GetDetElemID());
607     mainSDigit->SetPlane(hit->GetPlane());
608     mainSDigit->AddMCLabel(hit->GetTrack()); 
609
610     Int_t xPixel = -1;
611     Int_t yPixel = -1;
612     if (fSegmentation->Hit2PixelID(hit->X(), hit->Y(), mainSDigit->GetDetElemID(), xPixel, yPixel)) {
613       mainSDigit->SetPixID(xPixel, yPixel, 0);
614       mainSDigit->SetPixWidth(fSegmentation->GetPixelSizeX(mainSDigit->GetDetElemID()), 
615                           fSegmentation->GetPixelSizeY(mainSDigit->GetDetElemID()),
616                           fSegmentation->GetPixelSizeZ(mainSDigit->GetDetElemID()));  
617       mainSDigit->SetPixCenter(fSegmentation->GetPixelCenterX(mainSDigit->GetDetElemID(), xPixel), 
618                            fSegmentation->GetPixelCenterY(mainSDigit->GetDetElemID(), yPixel),
619                            fSegmentation->GetPixelCenterZ(mainSDigit->GetDetElemID(), 0));
620       new ((*fSideDigits)[fSideDigits->GetEntries()]) AliMFTDigit(*mainSDigit);
621       AliDebug(1, Form("Created new sdigit (%f, %f, %f) from hit (%f, %f, %f)",
622                        mainSDigit->GetPixelCenterX(), mainSDigit->GetPixelCenterY(), mainSDigit->GetPixelCenterZ(), hit->X(), hit->Y(), hit->Z()));
623     }
624
625     // creating "side digits" to simulate the effect of charge dispersion
626
627     Double_t pi4 = TMath::Pi()/4.;
628     for (Int_t iStep=0; iStep<fNStepForChargeDispersion; iStep++) {
629       Double_t shift = (iStep+1) * fSingleStepForChargeDispersion;
630       for (Int_t iAngle=0; iAngle<8; iAngle++) {
631         Double_t shiftX = shift*TMath::Cos(iAngle*pi4);
632         Double_t shiftY = shift*TMath::Sin(iAngle*pi4);
633         if (fSegmentation->Hit2PixelID(hit->X()+shiftX, hit->Y()+shiftY, hit->GetDetElemID(), xPixel, yPixel)) {
634           Bool_t digitExists = kFALSE;
635           for (Int_t iSideDigit=0; iSideDigit<fSideDigits->GetEntries(); iSideDigit++) {
636             if (xPixel==((AliMFTDigit*) fSideDigits->At(iSideDigit))->GetPixelX() && 
637                 yPixel==((AliMFTDigit*) fSideDigits->At(iSideDigit))->GetPixelY()) {
638               digitExists = kTRUE;
639               break;
640             }
641           }
642           if (!digitExists) {
643             AliMFTDigit *sideSDigit = new AliMFTDigit();
644             sideSDigit->SetEloss(0.);
645             sideSDigit->SetDetElemID(hit->GetDetElemID());
646             sideSDigit->SetPlane(hit->GetPlane());
647             sideSDigit->AddMCLabel(hit->GetTrack());
648             sideSDigit->SetPixID(xPixel, yPixel, 0);
649             sideSDigit->SetPixWidth(fSegmentation->GetPixelSizeX(sideSDigit->GetDetElemID()), 
650                                     fSegmentation->GetPixelSizeY(sideSDigit->GetDetElemID()),
651                                     fSegmentation->GetPixelSizeZ(sideSDigit->GetDetElemID()));  
652             sideSDigit->SetPixCenter(fSegmentation->GetPixelCenterX(sideSDigit->GetDetElemID(), xPixel), 
653                                      fSegmentation->GetPixelCenterY(sideSDigit->GetDetElemID(), yPixel),
654                                      fSegmentation->GetPixelCenterZ(sideSDigit->GetDetElemID(), 0)); 
655             new ((*fSideDigits)[fSideDigits->GetEntries()]) AliMFTDigit(*sideSDigit);
656           }
657         }
658       }
659     }
660     
661     // ------------ In case we should simulate a rectangular pattern of pixel...
662     
663     if (fSegmentation->GetPlane(mainSDigit->GetPlane())->HasPixelRectangularPatternAlongY()) {
664       for (Int_t iSDigit=0; iSDigit<fSideDigits->GetEntries(); iSDigit++) {
665         AliMFTDigit *mySDig = (AliMFTDigit*) (fSideDigits->At(iSDigit));
666         if (mySDig->GetPixelX()%2 == mySDig->GetPixelY()%2) {   // both pair or both odd
667           xPixel = mySDig->GetPixelX();
668           yPixel = mySDig->GetPixelY()+1;
669           if (fSegmentation->DoesPixelExist(mySDig->GetDetElemID(), xPixel, yPixel)) {
670             AliMFTDigit *newSDigit = new AliMFTDigit();
671             newSDigit->SetEloss(0.);
672             newSDigit->SetDetElemID(mySDig->GetDetElemID());
673             newSDigit->SetPlane(mySDig->GetDetElemID());
674             newSDigit->SetPixID(xPixel, yPixel, 0);
675             newSDigit->SetPixWidth(fSegmentation->GetPixelSizeX(newSDigit->GetDetElemID()), 
676                                    fSegmentation->GetPixelSizeY(newSDigit->GetDetElemID()),
677                                    fSegmentation->GetPixelSizeZ(newSDigit->GetDetElemID()));  
678             newSDigit->SetPixCenter(fSegmentation->GetPixelCenterX(newSDigit->GetDetElemID(), xPixel), 
679                                     fSegmentation->GetPixelCenterY(newSDigit->GetDetElemID(), yPixel),
680                                     fSegmentation->GetPixelCenterZ(newSDigit->GetDetElemID(), 0)); 
681             new ((*fSideDigits)[fSideDigits->GetEntries()]) AliMFTDigit(*newSDigit);
682           }
683         }
684         else {   // pair-odd
685           xPixel = mySDig->GetPixelX();
686           yPixel = mySDig->GetPixelY()-1;
687           if (fSegmentation->DoesPixelExist(mySDig->GetDetElemID(), xPixel, yPixel)) {
688             AliMFTDigit *newSDigit = new AliMFTDigit();
689             newSDigit->SetEloss(0.);
690             newSDigit->SetDetElemID(mySDig->GetDetElemID());
691             newSDigit->SetPlane(mySDig->GetPlane());
692             newSDigit->SetPixID(xPixel, yPixel, 0);
693             newSDigit->SetPixWidth(fSegmentation->GetPixelSizeX(newSDigit->GetDetElemID()), 
694                                    fSegmentation->GetPixelSizeY(newSDigit->GetDetElemID()),
695                                    fSegmentation->GetPixelSizeZ(newSDigit->GetDetElemID()));  
696             newSDigit->SetPixCenter(fSegmentation->GetPixelCenterX(newSDigit->GetDetElemID(), xPixel), 
697                                     fSegmentation->GetPixelCenterY(newSDigit->GetDetElemID(), yPixel),
698                                     fSegmentation->GetPixelCenterZ(newSDigit->GetDetElemID(), 0)); 
699             new ((*fSideDigits)[fSideDigits->GetEntries()]) AliMFTDigit(*newSDigit);
700           }
701         }
702       }
703     }
704
705     // -------- checking which pixels switched on have their diode actually within the charge dispersion radius
706
707     for (Int_t iSDigit=0; iSDigit<fSideDigits->GetEntries(); iSDigit++) {
708       AliMFTDigit *mySDig = (AliMFTDigit*) (fSideDigits->At(iSDigit));
709       Double_t distance = TMath::Sqrt(TMath::Power(mySDig->GetPixelCenterX()-hit->X(),2) + TMath::Power(mySDig->GetPixelCenterY()-hit->Y(),2));
710       if (fSegmentation->GetPlane(mySDig->GetPlane())->HasPixelRectangularPatternAlongY()) {
711         if (mySDig->GetPixelX()%2 == mySDig->GetPixelY()%2) {  // both pair or both odd
712           if (distance<fChargeDispersion) {
713             AliDebug(1, Form("Created new sdigit (%f, %f, %f) from hit (%f, %f, %f)",
714                              mySDig->GetPixelCenterX(), mySDig->GetPixelCenterY(), mySDig->GetPixelCenterZ(), hit->X(), hit->Y(), hit->Z()));
715             new ((*pSDigList[mySDig->GetPlane()])[pSDigList[mySDig->GetPlane()]->GetEntries()]) AliMFTDigit(*mySDig);
716             xPixel = mySDig->GetPixelX();
717             yPixel = mySDig->GetPixelY()+1;
718             if (fSegmentation->DoesPixelExist(mySDig->GetDetElemID(), xPixel, yPixel)) {
719               AliMFTDigit *newSDigit = new AliMFTDigit();
720               newSDigit->SetEloss(0.);
721               newSDigit->SetDetElemID(mySDig->GetDetElemID());
722               newSDigit->SetPlane(mySDig->GetPlane());
723               newSDigit->SetPixID(xPixel, yPixel, 0);
724               newSDigit->SetPixWidth(fSegmentation->GetPixelSizeX(newSDigit->GetDetElemID()), 
725                                      fSegmentation->GetPixelSizeY(newSDigit->GetDetElemID()),
726                                      fSegmentation->GetPixelSizeZ(newSDigit->GetDetElemID()));  
727               newSDigit->SetPixCenter(fSegmentation->GetPixelCenterX(newSDigit->GetDetElemID(), xPixel), 
728                                       fSegmentation->GetPixelCenterY(newSDigit->GetDetElemID(), yPixel),
729                                       fSegmentation->GetPixelCenterZ(newSDigit->GetDetElemID(), 0));
730               AliDebug(1, Form("Created new sdigit (%f, %f, %f) from hit (%f, %f, %f)",
731                                newSDigit->GetPixelCenterX(), newSDigit->GetPixelCenterY(), newSDigit->GetPixelCenterZ(), hit->X(), hit->Y(), hit->Z()));
732               new ((*pSDigList[newSDigit->GetPlane()])[pSDigList[newSDigit->GetPlane()]->GetEntries()]) AliMFTDigit(*newSDigit);
733             }
734           }
735         }
736       }
737       else {
738         if (distance<fChargeDispersion) {
739           AliDebug(1, Form("Created new sdigit (%f, %f, %f) from hit (%f, %f, %f)",
740                            mySDig->GetPixelCenterX(), mySDig->GetPixelCenterY(), mySDig->GetPixelCenterZ(), hit->X(), hit->Y(), hit->Z()));
741           new ((*pSDigList[mySDig->GetPlane()])[pSDigList[mySDig->GetPlane()]->GetEntries()]) AliMFTDigit(*mySDig);
742         }
743       }
744     }
745
746     fSideDigits->Delete(); 
747
748   }
749
750   AliDebug(1,"Exiting Hits2SDigitsLocal");
751
752 }
753
754 //====================================================================================================================================================
755
756 void AliMFT::MakeBranch(Option_t *option) {
757
758   // Create Tree branches 
759   AliDebug(1, Form("Start with option= %s.",option));
760   
761   const Int_t kBufSize = 4000;
762   
763   const Char_t *cH = strstr(option,"H");
764   const Char_t *cD = strstr(option,"D");
765   const Char_t *cS = strstr(option,"S");
766
767   if (cH && fLoader->TreeH()) {
768     CreateHits();
769     MakeBranchInTree(fLoader->TreeH(), "MFT", &fHits, kBufSize, 0);   
770   }
771
772   if (cS && fLoader->TreeS()) {
773     CreateSDigits();
774     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) MakeBranchInTree(fLoader->TreeS(), 
775                                                                     Form("Plane_%02d",iPlane), 
776                                                                     &((*fSDigitsPerPlane)[iPlane]), 
777                                                                     kBufSize, 0);
778   }
779   
780   if (cD && fLoader->TreeD()) {
781     CreateDigits();
782     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) MakeBranchInTree(fLoader->TreeD(), 
783                                                                     Form("Plane_%02d",iPlane), 
784                                                                     &((*fDigitsPerPlane)[iPlane]),
785                                                                     kBufSize, 0);
786   }
787
788   AliDebug(1,"Stop.");
789
790 }
791
792 //====================================================================================================================================================
793
794 void AliMFT::SetTreeAddress() {
795
796   AliDebug(1, "AliMFT::SetTreeAddress()");
797
798   //Set branch address for the Hits and Digits Tree.
799   AliDebug(1, "Start.");
800
801   AliDebug(1, Form("AliMFT::SetTreeAddress Hits  fLoader->TreeH() = %p\n", fLoader->TreeH()));
802   if (fLoader->TreeH() && fLoader->TreeH()->GetBranch("MFT")) {
803     CreateHits();
804     fLoader->TreeH()->SetBranchAddress("MFT", &fHits);
805   }
806     
807   AliDebug(1, Form("AliMFT::SetTreeAddress SDigits  fLoader->TreeS() = %p\n", fLoader->TreeS()));
808   if (fLoader->TreeS() && fLoader->TreeS()->GetBranch("Plane_00")) {
809     CreateSDigits();
810     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
811       fLoader->TreeS()->SetBranchAddress(Form("Plane_%02d",iPlane), &((*fSDigitsPerPlane)[iPlane]));
812     }
813   }
814     
815   AliDebug(1, Form("AliMFT::SetTreeAddress Digits  fLoader->TreeD() = %p\n", fLoader->TreeD()));
816   if (fLoader->TreeD() && fLoader->TreeD()->GetBranch("Plane_00")) {
817     CreateDigits(); 
818     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
819       fLoader->TreeD()->SetBranchAddress(Form("Plane_%02d",iPlane), &((*fDigitsPerPlane)[iPlane]));
820     }
821   }
822
823   AliDebug(1, Form("AliMFT::SetTreeAddress RecPoints  fLoader->TreeR() = %p\n", fLoader->TreeR()));
824   if (fLoader->TreeR() && fLoader->TreeR()->GetBranch("Plane_00")) {
825     CreateRecPoints(); 
826     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
827       fLoader->TreeR()->SetBranchAddress(Form("Plane_%02d",iPlane), &((*fRecPointsPerPlane)[iPlane]));
828     }
829   }
830
831   AliDebug(1,"Stop.");
832
833 }
834
835 //====================================================================================================================================================
836
837 void AliMFT::SetGeometry() {
838
839   AliInfo("AliMFT::SetGeometry\n");
840
841   fSegmentation = new AliMFTSegmentation(fNameGeomFile.Data());
842
843   fNPlanes = fSegmentation->GetNPlanes();
844
845 }
846
847 //====================================================================================================================================================
848
849 void AliMFT::CreateHits() { 
850
851   // create array of hits
852
853   AliDebug(1, "AliMFT::CreateHits()");
854
855   if (fHits) return;    
856   fHits = new TClonesArray("AliMFTHit");  
857
858 }
859
860 //====================================================================================================================================================
861
862 void AliMFT::CreateSDigits() { 
863  
864   // create sdigits list
865
866   AliDebug(1, "AliMFT::CreateSDigits()");
867  
868   if (fSDigitsPerPlane) return; 
869   fSDigitsPerPlane = new TObjArray(fNPlanes); 
870   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) fSDigitsPerPlane->AddAt(new TClonesArray("AliMFTDigit"), iPlane);
871
872   fSideDigits = new TClonesArray("AliMFTDigit");
873
874 }
875
876 //====================================================================================================================================================
877
878 void AliMFT::CreateDigits() {
879
880   // create digits list
881
882   AliDebug(1, "AliMFT::CreateDigits()");
883
884   if (fDigitsPerPlane) return; 
885   fDigitsPerPlane = new TObjArray(fNPlanes);
886   for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) fDigitsPerPlane->AddAt(new TClonesArray("AliMFTDigit"), iPlane);
887
888 }
889
890 //====================================================================================================================================================
891
892 void AliMFT::CreateRecPoints() {
893
894   // create recPoints list
895
896   AliDebug(1, "AliMFT::CreateRecPoints()");
897
898   if (fRecPointsPerPlane) return; 
899   fRecPointsPerPlane = new TObjArray(fNPlanes);
900   for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) fRecPointsPerPlane->AddAt(new TClonesArray("AliMFTCluster"), iPlane);
901
902 }
903
904 //====================================================================================================================================================