]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/UPGRADE/AliITSupgrade.cxx
bug fix
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSupgrade.cxx
CommitLineData
1d9af2d5 1// **************************************************************************
2// * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3// * *
4// * Author: The ALICE Off-line Project. *
5// * Contributors are mentioned in the code where appropriate. *
6// * *
7// * Permission to use, copy, modify and distribute this software and its *
8// * documentation strictly for non-commercial purposes is hereby granted *
9// * without fee, provided that the above copyright notice appears in all *
10// * copies and that both the copyright notice and this permission notice *
11// * appear in the supporting documentation. The authors make no claims *
12// * about the suitability of this software for any purpose. It is *
13// * provided "as is" without express or implied warranty. *
14// **************************************************************************
15
16/* $Id$ */
17
18#include <TArrayD.h> //new constructor
19#include <TFile.h>
20#include <TGeoManager.h> //CreateGeometry()
21#include <TGeoVolume.h> //CreateGeometry()
22#include <TVirtualMC.h> //->gMC in StepManager
23#include <TPDGCode.h> //StepHistory
24#include <TClonesArray.h>
2accaff0 25#include <TGeoGlobalMagField.h>
1d9af2d5 26#include "AliRun.h" //CreateMaterials()
27#include "AliMC.h" //CreateMaterials()
28#include "AliMagF.h" //CreateMaterials()
29#include "AliCDBManager.h" //CreateMaterials()
30#include "AliCDBEntry.h" //CreateMaterials()
31#include "AliITSupgrade.h" //class header
32#include "AliITShit.h" // StepManager()
33#include "AliITSDigitUpgrade.h"
34#include "AliTrackReference.h" // StepManager()
35#include "AliITSDetTypeSim.h"
36ClassImp(AliITSupgrade)
37//__________________________________________________________________________________________________
38 AliITSupgrade::AliITSupgrade():
39 AliITS(),
40 fWidths(0),
41 fRadii(0),
42 fRadiiCu(0),
43 fWidthsCu(0),
44 fCopper(0),
45 fBeampipe(0),
46 fRadiusBP(0),
47 fWidthBP(0),
48 fHalfLengthBP(0),
49 fNlayers(0),
50 fHalfLength(0),
51 fSdigits(0),
52 fDigits(0),
6ea6235b 53 //fClusters(0),
1d9af2d5 54 fSegmentation(0x0)
55{
56 //
57 // default ctor
387f1a9e 58 //
59}
60
61//__________________________________________________________________________________________________
62AliITSupgrade::AliITSupgrade(const char *name,const char *title, Bool_t isBeamPipe):
63 AliITS(name,title),
64 fWidths(0),
65 fRadii(0),
66 fRadiiCu(0),
67 fWidthsCu(0),
68 fCopper(0),
69 fBeampipe(isBeamPipe),
70 fRadiusBP(0),
71 fWidthBP(0),
72 fHalfLengthBP(0),
73 fNlayers(0),
74 fHalfLength(0),
75 fSdigits(0),
76 fDigits(0),
77 fSegmentation(0x0)
78{
79 fNlayers = 6;
80 fWidths.Set(fNlayers);
81 fRadii.Set(fNlayers);
82 fRadiiCu.Set(fNlayers);
83 fWidthsCu.Set(fNlayers);
84 fCopper.Set(fNlayers);
85 fHalfLength.Set(fNlayers);
86
87 // Default values are used in order to simulate the standard ITS material budget (see The ALICE Collaboration et al 2008 JINST 3 S08002 - Figure 3.2).
88 // The cell segmentation is chosen to achieve the tracking resolution as described in Table 3.2 of The ALICE Collaboration et al 2008 JINST 3 S08002,
89 // apart from SPD0 where the 9 um value is considered : resolution*sqrt(12)
ebaa7604 90 // Cilinder lenghts are chosen according to Table 1 of the Alignment paper : http://iopscience.iop.org/1748-0221/5/03/P03003
387f1a9e 91
92 Double_t xsz[6]={31.18*1e-04,41.6*1e-04,131.6*1e-04,131.6*1e-04,69.3*1e-04,69.3*1e-04};
93 Double_t zsz[6]={416*1e-04,416*1e-04,97*1e-04,97*1e-04,2875*1e-04,2875*1e-04};
ebaa7604 94 Double_t Lhalf[6]={14.1,14.1,22.2,29.7,43.1,48.9};
387f1a9e 95 Double_t r[6]={4.,7.6,14.9,23.8,39.1,43.6};
96 Double_t thick[6]={75.*1e-04,150.*1e-04,150.*1e-04,150.*1e-04,150.*1e-04,150.*1e-04};
97
98 Int_t nlayers =6;
99 TArrayD xsizeSi(nlayers);
100 TArrayD zsizeSi(nlayers);
101
102 // adding beam pipe upgrade
103 if(isBeamPipe){
104 fRadiusBP = 2.9;
105 fWidthBP = 0.08;
106 fHalfLengthBP = Lhalf[0];
107 }
108
109 // setting the geonetry parameters and the segmentation
110 Int_t npixHalf[6], npixR[6];
111 Double_t HalfL[6], xszInt[6];
112 Int_t c[6]={1,1,1,1,1,1};
113
114 for(Int_t i=0; i<nlayers; i++){
115 // recalc length
116 npixHalf[i]=(Int_t)(Lhalf[i]/zsz[i]);
117 HalfL[i]=npixHalf[i]*zsz[i];
118 // recalc segmentation
119 npixR[i] = (Int_t)(2*TMath::Pi()*r[i]/xsz[i]);
120 xszInt[i]= 2*TMath::Pi()*r[i]/npixR[i];
121 xsizeSi.AddAt(xszInt[i],i);
122 zsizeSi.AddAt(zsz[i],i);
123
124 fHalfLength.AddAt(HalfL[i],i);
125
126 fRadii.AddAt(r[i],i);
127 fWidths.AddAt(thick[i],i);
128
129 fRadiiCu.AddAt(r[i]+thick[i],i);
130 fWidthsCu.AddAt(0.015,i);//micron
131 fCopper.AddAt(c[i],i);
132 }
133
134 SetFullSegmentation(xsizeSi,zsizeSi);
135 Init();
1d9af2d5 136}
137//__________________________________________________________________________________________________
138AliITSupgrade::AliITSupgrade(const char *name,const char *title, TArrayD widths, TArrayD radii,TArrayD halfLengths, TArrayD radiiCu, TArrayD widthsCu, TArrayS copper,Bool_t bp,Double_t radiusBP, Double_t widthBP, Double_t halfLengthBP):
139 AliITS(name,title),
140 fWidths(0),
141 fRadii(0),
142 fRadiiCu(radiiCu),
143 fWidthsCu(widthsCu),
144 fCopper(copper),
145 fBeampipe(bp),
387f1a9e 146 fRadiusBP(0),
147 fWidthBP(0),
148 fHalfLengthBP(0),
1d9af2d5 149 fNlayers(widths.GetSize()),
150 fHalfLength(halfLengths),
151 fSdigits(0),
152 fDigits(0),
6ea6235b 153 //fClusters(0),
1d9af2d5 154 fSegmentation(0x0)
155{
156
157 //
158 // constructor
159 for(Int_t i=0;i<fNlayers;i++){
160 fWidths.Set(i+1);
161 fWidths.AddAt(widths.At(i),i);
162 fRadii.Set(i+1);
163 fRadii.AddAt(radii.At(i),i);
164 AliDebug(1,"Creating Volume");
165 }
387f1a9e 166
167 if(bp){
168 fRadiusBP=radiusBP;
169 fWidthBP=widthBP;
170 fHalfLengthBP=halfLengthBP;
171 }
1d9af2d5 172 Init();
173
174}
175
176AliITSupgrade::~AliITSupgrade(){
177
178 if(fSdigits) {fSdigits->Delete(); delete fSdigits;}
179 if(fDigits) {fDigits->Delete(); delete fDigits;}
180}
181
182
183
184//_________________________________________________________________________________________________
185void AliITSupgrade::AddAlignableVolumes()const
186{
187 //not needed
188 return;
189}
190
191//__________________________________________________________________________________________________
192
193void AliITSupgrade::CreateMaterials()
194{
195 //
196 // Definition of ITS materials
197 //
198
199 AliInfo("Start ITS materials");
200 //data from PDG booklet 2002 density [gr/cm^3] rad len [cm] abs len [cm]
201 Float_t aAir[4]={12,14,16,36} , zAir[4]={6,7,8,18} , wAir[4]={0.000124,0.755267,0.231781,0.012827} , dAir=0.00120479; Int_t nAir=4;//mixture
202 Float_t aBe = 9.012 , zBe = 4 , dBe = 1.848 , radBe = 65.19/dBe , absBe = 77.8/dBe ; // UPGRADE -> beryllium beampipe
203 Float_t aSi = 28.085 , zSi = 14 , dSi = 2.329 , radSi = 21.82/dSi , absSi = 108.4/dSi ; // UPGRADE -> Si tube
204 Float_t aCu = 63.54, zCu = 29 , dCu= 8.96 , radCu = 12.86/dCu, absCu= 137.3/dCu;//Upgrade -> Copper Tube
205
206 Int_t matId=0; //tmp material id number
207 Int_t unsens = 0, sens=1; //sensitive or unsensitive medium
208 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
209 Float_t maxfld = 5.; //max field value
210 Float_t tmaxfd = -10.0; //max deflection angle due to magnetic field in one step
211 Float_t deemax = - 0.2; //max fractional energy loss in one step
212 Float_t stemax = - 0.1; //max step allowed [cm]
213 Float_t epsil = 0.001; //abs tracking precision [cm]
214 Float_t stmin = - 0.001; //min step size [cm] in continius process transport, negative value: choose it automatically
215 Float_t tmaxfdSi = 0.1; // .10000E+01; // Degree
216 Float_t stemaxSi = 0.0075; // .10000E+01; // cm
217 Float_t deemaxSi = 0.1; // 0.30000E-02; // Fraction of particle's energy 0<deemax<=1
218 Float_t epsilSi = 1.0E-4;// .10000E+01;
219 Float_t stminSi = 0.0; // cm "Default value used"
220
221 Float_t epsilBe = .001; // Tracking precision,
222 Float_t stemaxBe = -0.01; // Maximum displacement for multiple scat
223 Float_t tmaxfdBe = -20.; // Maximum angle due to field deflection
224 Float_t deemaxBe = -.3; // Maximum fractional energy loss, DLS
225 Float_t stminBe = -.8;
2accaff0 226 Int_t isxfld = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ();//from CreateMaterials in STRUCT/AliPIPEv3.cxx
227 Float_t sxmgmx = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();//from CreateMaterials in STRUCT/AliPIPEv3.cxx
1d9af2d5 228
229
230 AliMixture(++matId,"UpgradeAir" ,aAir ,zAir ,dAir ,nAir ,wAir );
231 AliMedium(kAir ,"UpgradeAir" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
232
233 AliMaterial(++matId,"UpgradeSi" ,aSi ,zSi ,dSi ,radSi ,absSi );
234 AliMedium(kSi ,"UpgradeSi" , matId, sens, itgfld, maxfld, tmaxfdSi, stemaxSi, deemaxSi, epsilSi, stminSi);
235
236 AliMaterial(++matId,"UpgradeBe" ,aBe ,zBe ,dBe ,radBe ,absBe );
2accaff0 237 AliMedium(kBe ,"UpgradeBe" , matId, unsens, isxfld, sxmgmx, tmaxfdBe, stemaxBe, deemaxBe, epsilBe, stminBe);
1d9af2d5 238
239 AliMaterial(++matId, "UpgradeCu", aCu, zCu, dCu, radCu, absCu);
240 AliMedium(kCu, "UpgradeCu", matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
241
242 AliInfo("End ITS materials");
243
244}//void AliITS::CreateMaterials()
245//__________________________________________________________________________________________________
246void AliITSupgrade::CreateGeometry()
247{
248 //
249 //Creates detailed geometry simulation (currently GEANT volumes tree)
250 //
251 AliInfo("Start ITS upgrade preliminary version building");
252 if(!gMC->IsRootGeometrySupported()) return;
253 TGeoVolumeAssembly *vol= CreateVol();
254 gGeoManager->GetVolume("ALIC")->AddNode(vol,0);
255 AliInfo("Stop ITS upgrade preliminary version building");
387f1a9e 256 PrintSummary();
1d9af2d5 257}
258//__________________________________________________________________________________________________
259void AliITSupgrade::Init()
260{
261 // This method defines ID for sensitive volumes, i.e. such geometry volumes for which there are if(gMC->CurrentVolID()==XXX)
262 // statements in StepManager()
263 // Arguments: none
264 // Returns: none
265 AliDebug(1,"Init ITS upgrade preliminary version.");
266}
267//__________________________________________________________________________________________________
268void AliITSupgrade::StepManager()
269{
270 // Full Step Manager.
271 // Arguments: none
272 // Returns: none
273 // StepHistory(); return; //uncomment to print tracks history
274 // StepCount(); return; //uncomment to count photons
275 if(!(this->IsActive())) return;
276 if(!(gMC->TrackCharge())) return;
277 TString volumeName=gMC->CurrentVolName();
278 if(gMC->IsTrackExiting() && !volumeName.Contains("Cu") && !volumeName.Contains("Be")) {
279 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kITS);
280 } // if Outer ITS mother Volume
281
282 static TLorentzVector position, momentum; // Saves on calls to construtors
283 static AliITShit hit;// Saves on calls to constructors
284
285
286 Int_t status = 0;
287
288 //
289 // Track status
290 if(gMC->IsTrackInside()) status += 1;
291 if(gMC->IsTrackEntering()) status += 2;
292 if(gMC->IsTrackExiting()) status += 4;
293 if(gMC->IsTrackOut()) status += 8;
294 if(gMC->IsTrackDisappeared()) status += 16;
295 if(gMC->IsTrackStop()) status += 32;
296 if(gMC->IsTrackAlive()) status += 64;
297
298 //
299 // Fill hit structure.
300 //
301 TString volname = gMC->CurrentVolName();
302 if(volname.Contains("Cu"))return;
303 if(volname.Contains("Be"))return;
b178be37 304 volname.Remove(0,12); // remove letters to get the layer number
1d9af2d5 305 hit.SetModule(volname.Atoi()); // this will be the layer, not the module
306 hit.SetTrack(gAlice->GetMCApp()->GetCurrentTrackNumber());
307
308 gMC->TrackPosition(position);
309 gMC->TrackMomentum(momentum);
310 hit.SetPosition(position);
311
312 hit.SetTime(gMC->TrackTime());
313 hit.SetMomentum(momentum);
314 hit.SetStatus(status);
315 hit.SetEdep(gMC->Edep());
316 hit.SetShunt(GetIshunt());
317 if(gMC->IsTrackEntering()){
318 hit.SetStartPosition(position);
319 hit.SetStartTime(gMC->TrackTime());
320 hit.SetStartStatus(status);
321 return; // don't save entering hit.
322 }
387f1a9e 323 // Fill hit structure with this new hit.
1d9af2d5 324
325 new((*fHits)[fNhits++]) AliITShit(hit); // Use Copy Construtor.
326 // Save old position... for next hit.
327 hit.SetStartPosition(position);
328 hit.SetStartTime(gMC->TrackTime());
329 hit.SetStartStatus(status);
330 return;
331
332}
333//__________________________________________________________________________________________________
334TGeoVolumeAssembly * AliITSupgrade::CreateVol()
335{
336 TGeoVolumeAssembly *vol = new TGeoVolumeAssembly("ITSupgrade");
337 TGeoMedium *si =gGeoManager->GetMedium("ITS_UpgradeSi");
338 TGeoMedium *cu =gGeoManager->GetMedium("ITS_UpgradeCu");
339 TGeoMedium *be =gGeoManager->GetMedium("ITS_UpgradeBe");
340 for(Int_t ivol=0;ivol<fNlayers;ivol++){
341 TGeoVolume *layer=gGeoManager->MakeTube(Form("LayerSilicon%i",ivol),si , fRadii.At(ivol) , fRadii.At(ivol)+fWidths.At(ivol) , fHalfLength.At(ivol)); //upgraded situation
342
343 TGeoVolume *layerCu=gGeoManager->MakeTube(Form("LayerCu%i",ivol),cu , fRadiiCu.At(ivol) , fRadiiCu.At(ivol)+fWidthsCu.At(ivol) , fHalfLength.At(ivol) ); //upgraded situation
344
345
346 vol ->AddNode(layer,ivol);
347 if(fCopper.At(ivol)){
348 vol->AddNode(layerCu,ivol);
349 }
350 }
351 TGeoVolume *beampipe=gGeoManager->MakeTube("BeamPipe", be , fRadiusBP , fRadiusBP+ fWidthBP , fHalfLengthBP ); //upgraded situation
352
353 if(fBeampipe) vol->AddNode(beampipe,0);
354 return vol;
355
356}
357//_________________________________________________________________________________________________
358void AliITSupgrade::SetFullSegmentation(TArrayD xsize,TArrayD zsize){
359 TFile *file= TFile::Open("Segmentation.root","recreate");
360 file->WriteObjectAny(&xsize,"TArrayD","CellSizeX");
361 file->WriteObjectAny(&zsize,"TArrayD","CellSizeZ");
362 file->Close();
363}
364//_________________________________________________________________________________________________
365void AliITSupgrade::StepHistory()
366{
367 // This methode is invoked from StepManager() in order to print out
368 static Int_t iStepN;
369 const char *sParticle;
370 switch(gMC->TrackPid()){
371 case kProton: sParticle="PROTON" ;break;
372 case kNeutron: sParticle="neutron" ;break;
373 case kGamma: sParticle="gamma" ;break;
374 case kPi0: sParticle="Pi0" ;break;
375 case kPiPlus: sParticle="Pi+" ;break;
376 case kPiMinus: sParticle="Pi-" ;break;
377 case kElectron: sParticle="electron" ;break;
378 default: sParticle="not known" ;break;
379 }
380
381 TString flag="fanny combination";
382 if(gMC->IsTrackAlive())
383 if(gMC->IsTrackEntering()) flag="enters to";
384 else if(gMC->IsTrackExiting()) flag="exits from";
385 else if(gMC->IsTrackInside()) flag="inside";
386 else if(gMC->IsTrackStop()) flag="stopped in";
387
388 Int_t vid=0,copy=0;
389 TString path=gMC->CurrentVolName(); path.Prepend("-");path.Prepend(gMC->CurrentVolOffName(1));//current volume and his mother are always there
390 vid=gMC->CurrentVolOffID(2,copy); if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
391 vid=gMC->CurrentVolOffID(3,copy); if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
392
393 AliDebug(10, Form("Step %i: %s (%i) %s %s m=%.6f GeV q=%.1f dEdX=%.4f Etot=%.4f",iStepN,sParticle,gMC->TrackPid(),flag.Data(),path.Data(),gMC->TrackMass(),gMC->TrackCharge(),gMC->Edep()*1e9,gMC->Etot()));
394
395 Double_t gMcTrackPos[3]; gMC->TrackPosition(gMcTrackPos[0],gMcTrackPos[1],gMcTrackPos[2]);
396 Double_t gMcTrackPosLoc[3]; gMC->Gmtod(gMcTrackPos,gMcTrackPosLoc,1);
397 AliDebug(10,Form("gMC Track Position (MARS) x: %5.3lf, y: %5.3lf, z: %5.3lf (r: %5.3lf) ---> (LOC) x: %5.3f, y: %5.3f, z: %5.3f",gMcTrackPos[0],gMcTrackPos[1],gMcTrackPos[2],TMath::Sqrt(gMcTrackPos[0]*gMcTrackPos[0]+gMcTrackPos[1]*gMcTrackPos[1]+gMcTrackPos[2]*gMcTrackPos[2]),gMcTrackPosLoc[0],gMcTrackPosLoc[1],gMcTrackPosLoc[2]));
398
399
400 AliDebug(10,Form("Step %i: tid=%i flags alive=%i disap=%i enter=%i exit=%i inside=%i out=%i stop=%i new=%i",
401 iStepN, gAlice->GetMCApp()->GetCurrentTrackNumber(),
402 gMC->IsTrackAlive(), gMC->IsTrackDisappeared(),gMC->IsTrackEntering(), gMC->IsTrackExiting(),
403 gMC->IsTrackInside(),gMC->IsTrackOut(), gMC->IsTrackStop(), gMC->IsNewTrack()));
404
405 Float_t a,z,den,rad,abs; a=z=den=rad=abs=-1;
406 Int_t mid=gMC->CurrentMaterial(a,z,den,rad,abs);
407 AliDebug(10, Form("Step %i: mid=%i a=%7.2f z=%7.2f den=%9.4f rad=%9.2f abs=%9.2f\n\n",iStepN,mid,a,z,den,rad,abs));
408
409 TArrayI proc; gMC->StepProcesses(proc);
410 AliDebug(1,"Processes in this step:");
411 for ( int i = 0 ; i < proc.GetSize(); i++)
412 {
413 AliDebug(10, Form("%s",TMCProcessName[proc.At(i)]));
414 }
415 AliDebug(1,"End process list");
416 iStepN++;
417}//StepHistory()
418//______________________________________________________
419void AliITSupgrade::Hits2SDigits(){
420
421 // Interface method ivoked from AliSimulation to create a list of sdigits corresponding to list of hits. Every hit generates one or more sdigits.
422 // Arguments: none
423 // Returns: none
424 AliDebug(1,"Start Hits2SDigits.");
425
426 if(!fSegmentation)fSegmentation=new AliITSsegmentationUpgrade();
427
428 if(!GetLoader()->TreeH()) {GetLoader()->LoadHits();}
429 if(!GetLoader()->TreeS())
430
431 for(Int_t iEvt=0;iEvt < GetLoader()->GetRunLoader()->GetNumberOfEvents();iEvt++){//events loop
432 GetLoader()->GetRunLoader()->GetEvent(iEvt); //get next event
433 {
434 GetLoader()->MakeTree("S");
435 MakeBranch("S");
436 }
437 SetTreeAddress();
438 Int_t nSdigit[10]; for(Int_t i=0;i<10;i++) nSdigit[i] =0;
439 for(Int_t iEnt=0;iEnt<GetLoader()->TreeH()->GetEntries();iEnt++){//prims loop
440 GetLoader()->TreeH()->GetEntry(iEnt);
441 Hit2SumDig(Hits(),SDigitsList(),nSdigit);//convert this hit to list of sdigits
442 }//prims loop
443 GetLoader()->TreeS()->Fill();
444 GetLoader()->WriteSDigits("OVERWRITE");
445 SDigitsReset();
446 }//events loop
447 GetLoader()->UnloadHits();
448 GetLoader()->UnloadSDigits();
449 AliDebug(1,"Stop Hits2SDigits.");
450
451}
452//____________________________
453void AliITSupgrade::Hit2SumDig(TClonesArray *hits,TObjArray *pSDig, Int_t *nSdigit)
454{
455 // Adds sdigits of this hit to the list
456 // Returns: none
457
458 AliDebug(1,"Start Hit2SumDig");
459
460
461 if(!fSegmentation){
462 AliDebug(1,"Segmentation Not inizialized!!");
463 return ;
464 }
465
466 TClonesArray *pSdigList[10];
467
468 for(Int_t i=0;i<fNlayers;i++){
469 pSdigList[i]=(TClonesArray*)(*pSDig)[i];
470 if(pSdigList[i]->GetEntries()!=0) AliErrorClass("Some of sdigits lists is not empty"); //in principle those lists should be empty
471 }
472
1d9af2d5 473 for(Int_t iHit=0;iHit<hits->GetEntries();iHit++){ //hits loop
474 AliITShit *hit = (AliITShit*)hits->At(iHit);
475 Double_t xz[2];
476 if(!fSegmentation->GlobalToDet(hit->GetModule(),hit->GetXG(),hit->GetYG(),hit->GetZG(),xz[0],xz[1])) continue;
477 AliITSDigitUpgrade digit;
478 digit.SetSignal(hit->GetIonization());
479 digit.SetNelectrons(hit->GetIonization()/(3.62*1e-09));
480 digit.SetLayer(hit->GetModule());
481 digit.SetTrackID(hit->GetTrack());
482
483 Int_t xpix = 999;
484 if(fSegmentation->GetCellSizeX(hit->GetModule())!=0) xpix =(Int_t) (xz[0]/ fSegmentation->GetCellSizeX(hit->GetModule()));
485
486 Int_t zpix = 999; // shift at the next line to have zpix always positive. Such numbers are used to build the Pixel Id in the layer (> 0!)
487 if(fSegmentation->GetCellSizeZ(hit->GetModule())!=0){
488 zpix =(Int_t)((xz[1]+fHalfLength.At(hit->GetModule())) / fSegmentation->GetCellSizeZ(hit->GetModule()));
489 }
490 digit.SetPixId(xpix,zpix);
491
492
493 new((*pSdigList[hit->GetModule()])[nSdigit[hit->GetModule()]++]) AliITSDigitUpgrade(digit);
494 }
495
496 AliDebug(1,"Stop Hit2SumDig.");
497}
498//_______________________________________________________________________________________________
499void AliITSupgrade::MakeBranch(Option_t *option){
1d9af2d5 500 //Create Tree branches
501 AliDebug(1,Form("Start with option= %s.",option));
502
503 const Int_t kBufSize = 4000;
504
505 const char *cH = strstr(option,"H");
506 const char *cD = strstr(option,"D");
6ea6235b 507 //const char *cR = strstr(option,"R");
1d9af2d5 508 const char *cS = strstr(option,"S");
509
510 if(cH&&fLoader->TreeH()){
511
512 HitCreate();
513 MakeBranchInTree(fLoader->TreeH(),"ITSupgrade",&fHits,kBufSize,0);
514 }
515
516
517 if(cS&&fLoader->TreeS()){
518 SDigitsCreate();
519 for(Int_t i=0;i<fNlayers;i++) MakeBranchInTree(fLoader->TreeS(),Form("Layer%d",i),&((*fSdigits)[i]),kBufSize,0);
520 }
521
522
523 if(cD&&fLoader->TreeD()){
524 DigitsCreate();
525 for(Int_t i=0;i<fNlayers;i++) MakeBranchInTree(fLoader->TreeD(),Form("Layer%d",i),&((*fDigits)[i]),kBufSize,0);
526 }
1d9af2d5 527 AliDebug(1,"Stop.");
528}
529//____________________________________________________________________________________________________
530void AliITSupgrade::SetTreeAddress()
531{
532 //Set branch address for the Hits and Digits Tree.
533 AliDebug(1,"Start.");
534 if(fLoader->TreeH() && fLoader->TreeH()->GetBranch("ITSupgrade" )){
535 HitCreate();
536 fLoader->TreeH()->SetBranchAddress("ITSupgrade",&fHits);
537
538 }
1d9af2d5 539
540 if(fLoader->TreeS() && fLoader->TreeS()->GetBranch("Layer0" )){
541 SDigitsCreate();
542 for(int i=0;i<fNlayers;i++){
543 fLoader->TreeS()->SetBranchAddress(Form("Layer%d",i),&((*fSdigits)[i]));
544 }
545 }
1d9af2d5 546
547 if(fLoader->TreeD() && fLoader->TreeD()->GetBranch("Layer0")){
548 DigitsCreate();
549 for(int i=0;i<fNlayers;i++){
550 fLoader->TreeD()->SetBranchAddress(Form("Layer%d",i),&((*fDigits)[i]));
551 }
552 }
1d9af2d5 553 AliDebug(1,"Stop.");
554}
387f1a9e 555//____________________________________________________________________________________________________
556void AliITSupgrade::PrintSummary()
557{
558
559 // pdg booklet or http://pdg.lbl.gov/2010/reviews/rpp2010-rev-atomic-nuclear-prop.pdf
560
561 // X0 as X0[g/cm^2]/density
562 Double_t beX0= 65.19/1.848;
563 Double_t cuX0= 12.86/8.96;
564 Double_t siX0= 21.82/2.329;
565
566 if(!fSegmentation) fSegmentation= new AliITSsegmentationUpgrade();
567 printf(" %10s %10s %10s %10s %10s %10s %10s %10s\n","Name","R [cm]","x/X0 [%]","Phi res[um]","Z res[um]","X segm[um]","Y segm[um]","widths [um]");
ebaa7604 568
569 if(fBeampipe) printf(" Beampipe %10.3f %10f %10s %10s %10s %10s %10.1f\n", fRadiusBP, 100.*fWidthBP/beX0,"-","-","-","-",fWidthBP*1.e+4);
570 else printf(" * No Beam Pipe Upgrade \n");
387f1a9e 571 for(Int_t i=0; i< fNlayers; i++){
572 printf(" Si Layer %1i %10.3f %10f %10.2f %10.2f %10.2f %10.2f %10.1f\n",i, fRadii.At(i), 100.*fWidths.At(i)/siX0, (fSegmentation->GetCellSizeX(i) *1.e+4)/TMath::Sqrt(12.), (fSegmentation->GetCellSizeZ(i)*1.e+4)/TMath::Sqrt(12.),fSegmentation->GetCellSizeX(i) *1.e+4, fSegmentation->GetCellSizeZ(i)*1.e+4,fWidths.At(i)*1.e+4);
573 printf(" Cu Layer %10.3f %10f %10s %10s %10s %10s %10.1f\n", fRadiiCu.At(i), 100.*fWidthsCu.At(i)/cuX0,"-","-","-","-",fWidthsCu.At(i)*1.e+4);
574 }
575}
576//____________________________________________________________________________________________________
577void AliITSupgrade::SetSegmentationX(Double_t x, Int_t lay){
578TFile *f = TFile::Open("Segmentation.root","UPDATE");
579if(!f) {
580 AliError("\n\n\n Segmentation.root file does not exist. The segmentation is 0! \n\n\n");
581 return;
582 }
583else {
584 TArrayD *a = (TArrayD*)f->Get("CellSizeX");
585 if(!a) {
586 AliError("CessSizeX array does not exist!");
587 return;
588 }
589 a->AddAt(x,lay);
590 f->WriteObjectAny(a,"TArrayD","CellSizeX");
591 f->Close();
592 delete f;
593 }
594}
595//____________________________________________________________________________________________________
596void AliITSupgrade::SetSegmentationZ(Double_t z, Int_t lay){
597
598TFile *f = TFile::Open("Segmentation.root","UPDATE");
599if(!f) {
600 AliError("\n\n\n Segmentation.root file does not exist. The segmentation is 0! \n\n\n");
601 return;
602 }
603else {
604 TArrayD *a = (TArrayD*)f->Get("CellSizeZ");
605 if(!a) {
606 AliError("CellSizeZ array does not exist!");
607 return;
608 }
609 a->AddAt(z,lay);
610 f->WriteObjectAny(a,"TArrayD","CellSizeZ");
611 f->Close();
612 delete f;
613 }
614
615}