1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.8 1999/09/29 09:24:23 fca
19 Introduction of the Copyright and cvs Log
23 ////////////////////////////////////////////////
24 // Manager and hits classes for set:PHOS //
25 ////////////////////////////////////////////////
27 // --- ROOT system ---
36 // --- Standard library ---
42 // --- galice header files ---
46 //______________________________________________________________________________
51 //______________________________________________________________________________
53 AliPHOS::~AliPHOS(void)
55 delete fHits; // 28.12.1998
56 delete fTreePHOS; // 28.12.1998
61 //______________________________________________________________________________
66 fBranchNameOfCradles ("AliPHOSCradles"),
71 if( NULL==(fCradles=new TObjArray) )
73 Error("AliPHOS","Can not create fCradles");
79 //______________________________________________________________________________
81 AliPHOS::AliPHOS(const char *name, const char *title)
82 : AliDetector (name,title),
85 fBranchNameOfCradles ("AliPHOSCradles"),
90 <img src="picts/aliphos.gif">
94 fHits = new TClonesArray("AliPHOShit", 405);
98 SetMarkerColor(kGreen);
102 if( NULL==(fCradles=new TObjArray) ) {
103 Error("AliPHOS","Can not create fCradles");
109 //______________________________________________________________________________
111 void AliPHOS::DefPars()
156 //______________________________________________________________________________
158 void AliPHOS::AddHit(Int_t track, Int_t *vol, Float_t *hits)
160 TClonesArray &lhits = *fHits;
161 new(lhits[fNhits++]) AliPHOShit(fIshunt,track,vol,hits);
164 //___________________________________________
165 void AliPHOS::BuildGeometry()
170 const int kColorPHOS = kRed;
172 Top=gAlice->GetGeometry()->GetNode("alice");
176 Float_t pphi=12.9399462;
177 new TRotMatrix("rot988","rot988",90,-3*pphi,90,90-3*pphi,0,0);
178 new TRotMatrix("rot989","rot989",90,- pphi,90,90- pphi,0,0);
179 new TRotMatrix("rot990","rot990",90, pphi,90,90+ pphi,0,0);
180 new TRotMatrix("rot991","rot991",90, 3*pphi,90,90+3*pphi,0,0);
181 new TBRIK("S_PHOS","PHOS box","void",107.3,40,130);
183 Node = new TNode("PHOS1","PHOS1","S_PHOS",-317.824921,-395.014343,0,"rot988");
184 Node->SetLineColor(kColorPHOS);
187 Node = new TNode("PHOS2","PHOS2","S_PHOS",-113.532333,-494.124908,0,"rot989");
189 Node->SetLineColor(kColorPHOS);
191 Node = new TNode("PHOS3","PHOS3","S_PHOS", 113.532333,-494.124908,0,"rot990");
192 Node->SetLineColor(kColorPHOS);
195 Node = new TNode("PHOS4","PHOS4","S_PHOS", 317.824921,-395.014343,0,"rot991");
196 Node->SetLineColor(kColorPHOS);
200 //___________________________________________
201 void AliPHOS::CreateMaterials()
203 // *** DEFINITION OF AVAILABLE PHOS MATERIALS ***
205 // CALLED BY : PHOS_MEDIA
206 // ORIGIN : NICK VAN EIJNDHOVEN
210 Int_t ISXFLD = gAlice->Field()->Integ();
211 Float_t SXMGMX = gAlice->Field()->Max();
213 // --- The PbWO4 crystals ---
214 Float_t ax[3] = { 207.19,183.85,16. };
215 Float_t zx[3] = { 82.,74.,8. };
216 Float_t wx[3] = { 1.,1.,4. };
218 // --- Stainless Steel ---
219 Float_t as[5] = { 55.847,12.011,51.9961,58.69,28.0855 };
220 Float_t zs[5] = { 26.,6.,24.,28.,14. };
221 Float_t ws[5] = { .6392,8e-4,.2,.14,.02 };
223 // --- The polysterene scintillator (CH) ---
224 Float_t ap[2] = { 12.011,1.00794 };
225 Float_t zp[2] = { 6.,1. };
226 Float_t wp[2] = { 1.,1. };
229 Float_t at[2] = { 12.011,1.00794 };
230 Float_t zt[2] = { 6.,1. };
231 Float_t wt[2] = { 1.,2. };
233 // --- Polystyrene foam ---
234 Float_t af[2] = { 12.011,1.00794 };
235 Float_t zf[2] = { 6.,1. };
236 Float_t wf[2] = { 1.,1. };
238 //--- Foam thermo insulation (actual chemical composition unknown yet!) ---
239 Float_t ati[2] = { 12.011,1.00794 };
240 Float_t zti[2] = { 6.,1. };
241 Float_t wti[2] = { 1.,1. };
243 // --- Textolit (actual chemical composition unknown yet!) ---
244 Float_t atx[2] = { 12.011,1.00794 };
245 Float_t ztx[2] = { 6.,1. };
246 Float_t wtx[2] = { 1.,1. };
249 Int_t *idtmed = fIdtmed->GetArray()-699;
251 AliMixture( 0, "PbWO4$", ax, zx, dx, -3, wx);
252 AliMixture( 1, "Polystyrene$", ap, zp, dp, -2, wp);
253 AliMaterial( 2, "Al$", 26.98, 13., 2.7, 8.9, 999);
254 // --- Absorption length^ is ignored ---
255 AliMixture( 3, "Tyvek$", at, zt, dt, -2, wt);
256 AliMixture( 4, "Foam$", af, zf, df, -2, wf);
257 AliMixture( 5, "Stainless Steel$", as, zs, ds, 5, ws);
258 AliMaterial( 6, "Si$", 28.09, 14., 2.33, 9.36, 42.3);
259 AliMixture( 7, "Thermo Insul.$", ati, zti, dti, -2, wti);
260 AliMixture( 8, "Textolit$", atx, ztx, dtx, -2, wtx);
261 AliMaterial(99, "Air$", 14.61, 7.3, .001205, 30420., 67500);
263 AliMedium(0, "PHOS Xtal $", 0, 1, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
264 AliMedium(2, "Al parts $", 2, 0, ISXFLD, SXMGMX, 10., .1, .1, .001, .001);
265 AliMedium(3, "Tyvek wrapper$", 3, 0, ISXFLD, SXMGMX, 10., .1, .1, .001, .001);
266 AliMedium(4, "Polyst. foam $", 4, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
267 AliMedium(5, "Steel cover $", 5, 0, ISXFLD, SXMGMX, 10., .1, .1, 1e-4, 1e-4);
268 AliMedium(6, "Si PIN $", 6, 0, ISXFLD, SXMGMX, 10., .1, .1, .01, .01);
269 AliMedium(7, "Thermo Insul.$", 7, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
270 AliMedium(8, "Textolit $", 8, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
271 AliMedium(99, "Air $",99, 0, ISXFLD, SXMGMX, 10., 1., .1, .1, 10);
273 // --- Generate explicitly delta rays in the steel cover ---
274 gMC->Gstpar(idtmed[704], "LOSS", 3.);
275 gMC->Gstpar(idtmed[704], "DRAY", 1.);
276 // --- and in aluminium parts ---
277 gMC->Gstpar(idtmed[701], "LOSS", 3.);
278 gMC->Gstpar(idtmed[701], "DRAY", 1.);
281 //______________________________________________________________________________
283 void AliPHOS::AddPHOSCradles()
286 for(i=0;i<GetCradlesAmount();i++) {
288 int n = fCradles->GetEntries();
289 fCradles->Add(new AliPHOSCradle( IsVersion(), // geometry.
290 GetCrystalSideSize (),
299 GetCradleAngle (i)));
301 if( n+1 != fCradles->GetEntries() || NULL == fCradles->At(n) )
303 cout << " Can not create or add AliPHOSCradle.\n";
309 //______________________________________________________________________________
311 Int_t AliPHOS::DistancetoPrimitive(Int_t , Int_t )
316 //___________________________________________
322 for(i=0;i<35;i++) printf("*");
323 printf(" PHOS_INIT ");
324 for(i=0;i<35;i++) printf("*");
327 // Here the ABSO initialisation code (if any!)
328 for(i=0;i<80;i++) printf("*");
332 //______________________________________________________________________________
334 void AliPHOS::MakeBranch(Option_t *)
336 // ROOT output initialization to ROOT file.
338 // AliDetector::MakeBranch() is always called.
340 // There will be also special tree "PHOS" with one branch "AliPHOSCradles"
341 // if it was set next flag in the galice card file:
342 // * PHOSflags: YES: X<>0 NO: X=0
343 // * PHOSflags(1) : -----X. Create branch for TObjArray of AliPHOSCradle
347 // In that case special bit CradlesBranch_Bit will be set for AliPHOS
349 AliDetector::MakeBranch();
352 float t = GetPHOS_flag(0)/10;
354 i = (int) ((t-i)*10);
358 SetBit(CradlesBranch_Bit);
360 if( NULL==(fTreePHOS=new TTree(fTreeName.Data(),"PHOS events tree")) )
362 Error("MakeBranch","Can not create TTree");
366 if( NULL==fTreePHOS->GetCurrentFile() )
368 Error("MakeBranch","There is no opened ROOT file");
372 // Create a new branch in the current Root Tree.
374 if( NULL==fTreePHOS->Branch(fBranchNameOfCradles.Data(),"TObjArray",&fCradles,4000,0) )
376 Error("MakeBranch","Can not create branch");
380 printf("The branch %s has been created\n",fBranchNameOfCradles.Data());
383 //______________________________________________________________________________
385 void AliPHOS::SetTreeAddress(void)
387 // ROOT input initialization.
389 // AliDetector::SetTreeAddress() is always called.
391 // If CradlesBranch_Bit is set (see AliPHOS::MakeBranch) than fTreePHOS is
394 AliDetector::SetTreeAddress();
396 if( !TestBit(CradlesBranch_Bit) )
399 if( NULL==(fTreePHOS=(TTree*)gDirectory->Get((char*)(fTreeName.Data())) ) )
401 Error("SetTreeAddress","Can not find Tree \"%s\"\n",fTreeName.Data());
405 TBranch *branch = fTreePHOS->GetBranch(fBranchNameOfCradles.Data());
408 Error("SetTreeAddress","Can not find branch %s in TTree:%s",fBranchNameOfCradles.Data(),fTreeName.Data());
412 branch->SetAddress(&fCradles);
415 //______________________________________________________________________________
417 AliPHOSCradle *AliPHOS::GetCradleOfTheParticle(const TVector3 &p,const TVector3 &v) const
419 // For a given direction 'p' and source point 'v' returns pointer to AliPHOSCradle
420 // in that direction or NULL if AliPHOSCradle was not found.
422 for( int m=0; m<fCradles->GetEntries(); m++ )
424 AliPHOS *PHOS = (AliPHOS *)this; // Removing 'const'...
425 AliPHOSCradle *cradle = (AliPHOSCradle *)PHOS->fCradles->operator[](m);
428 const float d = cradle->GetRadius();
429 cradle->GetXY(p,v,d,x,y,l);
431 if( l>0 && TMath::Abs(x)<cradle->GetNz ()*cradle->GetCellSideSize()/2
432 && TMath::Abs(y)<cradle->GetNphi()*cradle->GetCellSideSize()/2 )
439 //______________________________________________________________________________
441 void AliPHOS::Reconstruction(Float_t signal_step, UInt_t min_signal_reject)
443 // Call AliPHOSCradle::Reconstruction(Float_t signal_step, UInt_t min_signal_reject)
444 // for all AliPHOSCradles.
446 for( int i=0; i<fCradles->GetEntries(); i++ )
447 GetCradle(i).Reconstruction(signal_step,min_signal_reject);
450 //______________________________________________________________________________
452 void AliPHOS::ResetDigits(void)
454 AliDetector::ResetDigits();
456 for( int i=0; i<fCradles->GetEntries(); i++ )
457 ((AliPHOSCradle*)(*fCradles)[i]) -> Clear();
460 //______________________________________________________________________________
462 void AliPHOS::FinishEvent(void)
464 // Called at the end of each 'galice' event.
466 if( NULL!=fTreePHOS )
470 //______________________________________________________________________________
472 void AliPHOS::FinishRun(void)
476 //______________________________________________________________________________
478 void AliPHOS::Print(Option_t *opt)
480 // Print PHOS information.
481 // For each AliPHOSCradle the function AliPHOSCradle::Print(opt) is called.
483 AliPHOS &PHOS = *(AliPHOS *)this; // Removing 'const'...
485 for( int i=0; i<fCradles->GetEntries(); i++ )
487 printf("PHOS cradle %d from %d\n",i+1, fCradles->GetEntries());
488 PHOS.GetCradle(i).Print(opt);
489 printf( "---------------------------------------------------\n");
493 //______________________________________________________________________________
494 void AliPHOS::SetFlags(Float_t p1,Float_t p2,Float_t p3,Float_t p4,
495 Float_t p5,Float_t p6,Float_t p7,Float_t p8,Float_t p9)
508 //______________________________________________________________________________
509 void AliPHOS::SetCell(Float_t p1,Float_t p2,Float_t p3,Float_t p4,
510 Float_t p5,Float_t p6,Float_t p7,Float_t p8,Float_t p9)
523 //______________________________________________________________________________
524 void AliPHOS::SetRadius(Float_t radius)
529 //______________________________________________________________________________
530 void AliPHOS::SetCradleSize(Int_t nz, Int_t nphi, Int_t ncradles)
534 PHOSsize[2]=ncradles;
537 //______________________________________________________________________________
538 void AliPHOS::SetCradleA(Float_t angle)
543 //______________________________________________________________________________
544 void AliPHOS::SetExtra(Float_t p1,Float_t p2,Float_t p3,Float_t p4,
545 Float_t p5,Float_t p6,Float_t p7,Float_t p8,Float_t p9)
558 //______________________________________________________________________________
559 void AliPHOS::SetTextolitWall(Float_t dx, Float_t dy, Float_t dz)
566 //______________________________________________________________________________
567 void AliPHOS::SetInnerAir(Float_t dx, Float_t dy, Float_t dz)
574 //______________________________________________________________________________
575 void AliPHOS::SetFoam(Float_t dx, Float_t dy, Float_t dz, Float_t dr)
583 ClassImp(AliPHOSCradle)
585 //______________________________________________________________________________
587 AliPHOSCradle::AliPHOSCradle(void) {}
589 //______________________________________________________________________________
591 AliPHOSCradle::AliPHOSCradle( int Geometry ,
592 float CrystalSideSize ,
593 float CrystalLength ,
594 float WrapThickness ,
602 fGeometry (Geometry),
604 // fChargedTracksInPIN (),
605 fCrystalSideSize (CrystalSideSize),
606 fCrystalLength (CrystalLength),
607 fWrapThickness (WrapThickness),
608 fAirThickness (AirThickness),
609 fPIN_SideSize (PIN_SideSize),
610 fPIN_Length (PIN_Length),
616 fCellEnergy = TH2F("CellE","Energy deposition in a cells",fNz,0,fNz,fNphi,0,fNphi);
617 fCellEnergy .SetDirectory(0);
618 fChargedTracksInPIN = TH2S("PINCtracks","Amount of charged tracks in PIN",fNz,0,fNz,fNphi,0,fNphi);
619 fChargedTracksInPIN .SetDirectory(0);
622 //______________________________________________________________________________
624 AliPHOSCradle::~AliPHOSCradle(void) // 28.12.1998
626 fGammasReconstructed.Delete();
627 fParticles .Delete();
630 //______________________________________________________________________________
632 void AliPHOSCradle::Clear(Option_t *)
634 // Clear digit. information.
636 fCellEnergy .Reset();
637 fChargedTracksInPIN .Reset();
638 GetParticles() .Delete();
639 GetParticles() .Compress();
640 GetGammasReconstructed() .Delete();
641 GetGammasReconstructed() .Compress();
645 //______________________________________________________________________________
647 void AliPHOSCradle::GetXY(const TVector3 &p,const TVector3 &v,float R,float &x,float &y,float &l) const
649 // This function calculates hit position (x,y) in the CRADLE cells plain from particle in
650 // the direction given by 'p' (not required to be normalized) and start point
651 // given by 3-vector 'v'. So the particle trajectory is t(l) = v + p*l
652 // were 'l' is a number (distance from 'v' to CRADLE cells plain) and 't' is resulting
653 // three-vector of trajectory point.
655 // After the call to this function user should test that l>=0 (the particle HITED the
656 // plain) and (x,y) are in the region of CRADLE:
659 // AliPHOSCradle cradle(......);
660 // TVector3 p(....), v(....);
662 // cradle.GetXY(p,v,x,y,l);
663 // if( l<0 || TMath::Abs(x)>cradle.GetNz() *cradle.GetCellSideSize()/2
664 // || TMath::Abs(y)>cradle.GetNphi()*cradle.GetCellSideSize()/2 )
665 // cout << "Outside the CRADLE.\n";
667 // We have to create three vectors:
668 // s - central point on the PHOS surface
669 // n1 - first vector in CRADLE plain
670 // n2 - second vector in CRADLE plain
671 // This three vectors are orthonormalized.
673 double phi = fPhi/180*TMath::Pi();
674 TVector3 n1( 0.0 , 0.0 , 1.0 ), // Z direction (X)
675 n2( -sin(phi) , cos(phi) , 0 ), // around beam (Y)
676 s ( R*cos(phi) , R*sin(phi) , 0 ); // central point
678 const double l1_min = 1e-2;
680 p_n1 = p*n1, // * - scalar product.
687 if ( TMath::Abs(l1=p.X()-n1.X()*p_n1-n2.X()*p_n2)>l1_min )
688 { l = (-v.X()+s.X()+n1.X()*(v_n1-s_n1)+n2.X()*(v_n2-s_n2))/l1; }
689 else if ( TMath::Abs(l1=p.Y()-n1.Y()*p_n1-n2.Y()*p_n2)>l1_min )
690 { l = (-v.Y()+s.Y()+n1.Y()*(v_n1-s_n1)+n2.Y()*(v_n2-s_n2))/l1; }
691 else if ( TMath::Abs(l1=p.Z()-n1.Z()*p_n1-n2.Z()*p_n2)>l1_min )
692 { l = (-v.Z()+s.Z()+n1.Z()*(v_n1-s_n1)+n2.Z()*(v_n2-s_n2))/l1; }
694 // double lx = (-v.X()+s.X()+n1.X()*(v.dot(n1)-s.dot(n1))+n2.X()*(v.dot(n2)-s.dot(n2)))/
695 // (p.X()-n1.X()*p.dot(n1)-n2.X()*p.dot(n2)),
696 // ly = (-v.Y()+s.Y()+n1.Y()*(v.dot(n1)-s.dot(n1))+n2.Y()*(v.dot(n2)-s.dot(n2)))/
697 // (p.Y()-n1.Y()*p.dot(n1)-n2.Y()*p.dot(n2)),
698 // lz = (-v.Z()+s.Z()+n1.Z()*(v.dot(n1)-s.dot(n1))+n2.Z()*(v.dot(n2)-s.dot(n2)))/
699 // (p.Z()-n1.Z()*p.dot(n1)-n2.Z()*p.dot(n2));
700 // cout.form("x: %g %g %g %g\n",lx,-v.X()+s.X()+n1.X()*(v.dot(n1)-s.dot(n1))+n2.X()*(v.dot(n2)-s.dot(n2)),p.X()-n1.X()*p.dot(n1)-n2.X()*p.dot(n2));
701 // cout.form("y: %g %g %g %g\n",lx,-v.Y()+s.Y()+n1.Y()*(v.dot(n1)-s.dot(n1))+n2.Y()*(v.dot(n2)-s.dot(n2)),p.Y()-n1.Y()*p.dot(n1)-n2.Y()*p.dot(n2));
702 // cout.form("z: %g %g %g %g\n",lx,-v.Z()+s.Z()+n1.Z()*(v.dot(n1)-s.dot(n1))+n2.Z()*(v.dot(n2)-s.dot(n2)),p.Z()-n1.Z()*p.dot(n1)-n2.Z()*p.dot(n2));
703 // cout.form("lx,ly,lz = %g,%g,%g\n",lx,ly,lz);
705 x = p_n1*l + v_n1 - s_n1;
706 y = p_n2*l + v_n2 - s_n2;
709 //______________________________________________________________________________
711 void AliPHOSCradle::Print(Option_t *opt)
713 // Print AliPHOSCradle information.
715 // options: 'd' - print energy deposition for EVERY cell
716 // 'p' - print particles list that hit the cradle
717 // 'r' - print list of reconstructed particles
719 AliPHOSCradle *cr = (AliPHOSCradle *)this; // Removing 'const'...
721 printf("AliPHOSCradle: Nz=%d Nphi=%d, fPhi=%f, E=%g\n",fNz,fNphi,fPhi,
722 cr->fCellEnergy.GetSumOfWeights());
724 if( NULL!=strchr(opt,'d') )
726 printf("\n\nCells Energy (in MeV):\n\n |");
727 for( int x=0; x<fNz; x++ )
731 for( int y=fNphi-1; y>=0; y-- )
734 for( int x=0; x<fNz; x++ )
735 printf("%6d",(int)(cr->fCellEnergy.GetBinContent(cr->fCellEnergy.GetBin(x,y))*1000));
741 if( NULL!=strchr(opt,'p') )
743 printf("This cradle was hit by %d particles\n",
744 ((AliPHOSCradle*)this)->GetParticles().GetEntries());
745 TObjArray &p=((AliPHOSCradle*)this)->GetParticles();
746 for( int i=0; i<p.GetEntries(); i++ )
747 ((AliPHOSgamma*)(p[i]))->Print();
750 if( NULL!=strchr(opt,'p') )
752 printf("Amount of reconstructed gammas is %d\n",
753 ((AliPHOSCradle*)this)->GetGammasReconstructed().GetEntries());
755 TObjArray &p=((AliPHOSCradle*)this)->GetGammasReconstructed();
756 for( int i=0; i<p.GetEntries(); i++ )
757 ((AliPHOSgamma*)(p[i]))->Print();
761 //______________________________________________________________________________
763 void AliPHOSCradle::Distortion(const TH2F *Noise, const TH2F *Stochastic, const TH2F *Calibration)
765 // This function changes histogram of cell energies fCellEnergy on the base of input
766 // histograms Noise, Stochastic, Calibration. The histograms must have
769 //////////////////////////////////
770 // Testing the histograms size. //
771 //////////////////////////////////
773 if( fNz!=fCellEnergy.GetNbinsX() || fNphi!=fCellEnergy.GetNbinsY() )
775 printf ("Bad size of CellEnergy! Must be: Nz x Nphi = %d x %d\n"
776 "but size of CellEnergy is: %d x %d\n",
777 fNz,fNphi,fCellEnergy.GetNbinsX(),fCellEnergy.GetNbinsY());
781 if( fNz!=fChargedTracksInPIN.GetNbinsX() || fNphi!=fChargedTracksInPIN.GetNbinsY() )
783 printf ("Bad size of ChargedTracksInPIN! Must be: Nz x Nphi = %d x %d\n"
784 "but size of ChargedTracksInPIN is: %d x %d\n",
785 fNz,fNphi,fChargedTracksInPIN.GetNbinsX(),fChargedTracksInPIN.GetNbinsY());
789 if( NULL!=Noise && (fNz!=Noise->GetNbinsX() || fNphi!=Noise->GetNbinsX()) )
791 printf ("Bad size of Noise! Must be: Nz x Nphi = %d x %d\n"
792 "but size of Noise is: %d x %d\n",
793 fNz,fNphi,fChargedTracksInPIN.GetNbinsX(),fChargedTracksInPIN.GetNbinsY());
797 if( NULL!=Stochastic && (fNz!=Stochastic->GetNbinsX() || fNphi!=Stochastic->GetNbinsX()) )
799 printf ("Bad size of Stochastic! Must be: Nz x Nphi = %d x %d\n"
800 "but size of Stochastic is: %d x %d\n",
801 fNz,fNphi,fChargedTracksInPIN.GetNbinsX(),fChargedTracksInPIN.GetNbinsY());
805 if( NULL!=Calibration && (fNz!=Calibration->GetNbinsX() || fNphi!=Calibration->GetNbinsX()) )
807 printf ("Bad size of Calibration! Must be: Nz x Nphi = %d x %d\n"
808 "but size of Calibration is: %d x %d\n",
809 fNz,fNphi,fChargedTracksInPIN.GetNbinsX(),fChargedTracksInPIN.GetNbinsY());
817 for( int y=0; y<fNphi; y++ )
818 for( int x=0; x<fNz; x++ )
820 const int n = fCellEnergy.GetBin(x,y); // Bin number
823 Float_t E_old=fCellEnergy.GetBinContent(n), E_new=E_old;
825 if( NULL!=Stochastic )
826 E_new = r.Gaus(E_old,sqrt(E_old)*GetDistortedValue(Stochastic,n));
828 if( NULL!=Calibration )
829 E_new *= GetDistortedValue(Calibration,n);
832 E_new += GetDistortedValue(Noise,n);
834 fCellEnergy.SetBinContent(n,E_new);
838 ////////////////////////////////////////////////////////////////////////////////
840 TH2F* AliPHOSCradle::CreateHistForDistortion(const char *name, const char *title,
842 Float_t MU_mu, Float_t MU_sigma,
843 Float_t SIGMA_mu, Float_t SIGMA_sigma)
845 // Create (new TH2F(...)) histogram with information (for every bin) that will
846 // be used for VALUE creation.
847 // Two values will be created for each bin:
848 // MU = TRandom::Gaus(MU_mu,MU_sigma)
850 // SIGMA = TRandom::Gaus(SIGMA_mu,SIGMA_sigma)
851 // The VALUE in a particluar bin will be equal
852 // VALUE = TRandom::Gaus(MU,SIGMA)
854 // Do not forget to delete the histogram at the end of the work.
856 TH2F *h = new TH2F( name,title, Nx,1,Nx, Ny,1,Ny );
859 Error("CreateHistForDistortion","Can not create the histogram");
864 for( int y=0; y<Ny; y++ )
865 for( int x=0; x<Nx; x++ )
867 const int n = h->GetBin(x,y);
868 h->SetBinContent(n,r.Gaus( MU_mu, MU_sigma));
869 h->SetBinError (n,r.Gaus(SIGMA_mu,SIGMA_sigma));
875 ////////////////////////////////////////////////////////////////////////////////
877 Float_t AliPHOSCradle::GetDistortedValue(const TH2F *h, UInt_t n)
879 return r.Gaus(((TH2F*)h)->GetBinContent(n),n);
882 ////////////////////////////////////////////////////////////////////////////////
883 //______________________________________________________________________________
886 #define common_for_event_storing COMMON_FOR_EVENT_STORING
888 #define common_for_event_storing common_for_event_storing_
893 enum { crystals_matrix_amount_max=4, crystals_in_matrix_amount_max=40000 };
895 // Event-independent information
896 UShort_t crystals_matrix_amount_PHOS,
898 amount_of_crystals_on_Z,
899 amount_of_crystals_on_PHI;
903 matrix_coordinate_Z [crystals_matrix_amount_max],
904 matrix_coordinate_PHI [crystals_matrix_amount_max];
906 UShort_t crystals_amount_with_amplitudes [crystals_matrix_amount_max],
907 crystals_amplitudes_Iad [crystals_matrix_amount_max]
908 [crystals_in_matrix_amount_max][2];
909 } common_for_event_storing;
911 // integer*4 crystals_amount_max,crystals_in_matrix_amount_max,
912 // + crystals_matrix_amount_max
913 // parameter (crystals_matrix_amount_max=4)
914 // parameter (crystals_in_matrix_amount_max=40000)
915 // parameter (crystals_amount_max =crystals_matrix_amount_max*
916 // + crystals_in_matrix_amount_max)
918 // * All units are in GeV, cm, radian
919 // real crystal_amplitudes_unit, radius_unit,
920 // + crystal_size_unit, crystal_length_unit,
921 // + matrix_coordinate_Z_unit, matrix_coordinate_PHI_unit
922 // integer crystal_amplitudes_in_units_min
923 // parameter (crystal_amplitudes_in_units_min = 1)
924 // parameter (crystal_amplitudes_unit = 0.001 ) ! 1.0 MeV
925 // parameter (radius_unit = 0.1 ) ! 0.1 cm
926 // parameter (crystal_size_unit = 0.01 ) ! 0.01 cm
927 // parameter (crystal_length_unit = 0.01 ) ! 0.01 cm
928 // parameter (matrix_coordinate_Z_unit = 0.1 ) ! 0.1 cm
929 // parameter (matrix_coordinate_PHI_unit = 1e-4 ) ! 1e-4 radian
931 // integer*2 crystals_matrix_amount_PHOS, crystal_matrix_type,
932 // + amount_of_crystals_on_Z, amount_of_crystals_on_PHI,
933 // + crystals_amount_with_amplitudes, crystals_amplitudes_Iad
934 // integer*4 event_number
936 // real radius, crystal_size, crystal_length,
937 // + matrix_coordinate_Z, matrix_coordinate_PHI
939 // real crystals_amplitudes, crystals_energy_total
940 // integer event_file_unit_number
942 // common /common_for_event_storing/
943 // + ! Event-independent information
944 // + crystals_matrix_amount_PHOS,
945 // + crystal_matrix_type,
946 // + amount_of_crystals_on_Z,
947 // + amount_of_crystals_on_PHI,
951 // + matrix_coordinate_Z (crystals_matrix_amount_max),
952 // + matrix_coordinate_PHI (crystals_matrix_amount_max),
954 // + ! Event-dependent information
956 // + crystals_amount_with_amplitudes
957 // + (crystals_matrix_amount_max),
958 // + crystals_amplitudes_Iad (2,crystals_in_matrix_amount_max,
959 // + crystals_matrix_amount_max),
961 // + ! These information don't store in data file
962 // + crystals_amplitudes (crystals_amount_max),
963 // + crystals_energy_total,
964 // + event_file_unit_number
967 // parameter (NGp=1000,nsps=10,nvertmax=1000)
968 // COMMON /GAMMA/KG,MW(ngp),ID(ngp),JD(ngp),E(ngp),E4(ngp),
969 // , XW(ngp),YW(ngp),ES(nsps,ngp),ET(nsps,ngp),ISsd(ngp),
970 // , IGDEV(ngp),ZGDEV(ngp),sigexy(3,ngp),Emimx(2,nsps,ngp),
971 // , kgfix,igfix(ngp),cgfix(3,ngp),sgfix(3,ngp),hiw(ngp),
972 // , wsw(nsps,ngp),h1w(ngp),h0w(ngp),raxay(5,ngp),
973 // , sigmaes0(nsps,ngp),dispeces(nsps,ngp),
978 #define rcgamma RCGAMMA
980 #define rcgamma rcgamma_
985 enum {NGP=1000, nsps=10, nvertmax=1000};
986 int recons_gammas_amount, mw[NGP],ID[NGP],JD[NGP];
987 float E[NGP], E4[NGP], XW[NGP], YW[NGP], ES[NGP][nsps],ET[NGP][nsps],ISsd[NGP],
988 igdev[NGP],Zgdev[NGP];
989 // sigexy(3,ngp),Emimx(2,nsps,ngp),
990 // , kgfix,igfix(ngp),cgfix(3,ngp),sgfix(3,ngp),hiw(ngp),
991 // , wsw(nsps,ngp),h1w(ngp),h0w(ngp),raxay(5,ngp),
992 // , sigmaes0(nsps,ngp),dispeces(nsps,ngp),
997 #define reconsfirst RECONSFIRST
998 #define type_of_call _stdcall
1000 #define reconsfirst reconsfirst_
1001 #define type_of_call
1004 extern "C" void type_of_call reconsfirst(const float &,const float &);
1006 void AliPHOSCradle::Reconstruction(Float_t signal_step, UInt_t min_signal_reject)
1008 // Call of PHOS reconstruction program.
1009 // signal_step=0.001 GeV (1MeV)
1010 // min_signal_reject = 15 or 30 MeV
1013 common_for_event_storing.event_number = 0; // We do not know event number?
1014 common_for_event_storing.crystals_matrix_amount_PHOS = 1;
1015 common_for_event_storing.crystal_matrix_type = 1; // 1 - rectangular
1016 common_for_event_storing.amount_of_crystals_on_Z = fNz;
1017 common_for_event_storing.amount_of_crystals_on_PHI = fNphi;
1019 common_for_event_storing.radius = fRadius;
1020 common_for_event_storing.crystal_size = GetCellSideSize();
1021 common_for_event_storing.crystal_length = fCrystalLength;
1023 common_for_event_storing.matrix_coordinate_Z [0] = 0;
1024 common_for_event_storing.matrix_coordinate_PHI [0] = fPhi;
1026 #define k common_for_event_storing.crystals_amount_with_amplitudes[0]
1029 for( int y=0; y<fNphi; y++ )
1030 for( int x=0; x<fNz; x++ )
1032 UInt_t n = fCellEnergy.GetBin(x,y);
1033 UInt_t signal = (int) (fCellEnergy.GetBinContent(n)/signal_step);
1034 if( signal>=min_signal_reject )
1036 common_for_event_storing.crystals_amplitudes_Iad[0][k][0] = signal;
1037 common_for_event_storing.crystals_amplitudes_Iad[0][k][1] = x + y*fNz;
1043 GetGammasReconstructed().Delete();
1044 GetGammasReconstructed().Compress();
1046 const float stochastic_term = 0.03, // per cents over sqrt(E); E in GeV
1047 electronic_noise = 0.01; // GeV
1048 reconsfirst(stochastic_term,electronic_noise); // Call of reconstruction program.
1050 for( int i=0; i<rcgamma.recons_gammas_amount; i++ )
1052 // new (GetGammasReconstructed().UncheckedAt(i) ) AliPHOSgamma;
1053 // AliPHOSgamma &g = *(AliPHOSgamma*)(GetGammasReconstructed().UncheckedAt(i));
1055 AliPHOSgamma *gggg = new AliPHOSgamma;
1058 Error("Reconstruction","Can not create AliPHOSgamma");
1062 GetGammasReconstructed().Add(gggg);
1063 AliPHOSgamma &g=*gggg;
1065 Float_t thetta, alpha, betta, R=fRadius+rcgamma.Zgdev[i]/10;
1067 g.fX = rcgamma.YW[i]/10;
1068 g.fY = rcgamma.XW[i]/10;
1069 g.fE = rcgamma.E [i];
1071 thetta = atan(g.fX/R);
1073 alpha = atan(g.fY/R);
1074 betta = fPhi/180*TMath::Pi() + alpha;
1076 g.fPx = g.fE * cos(thetta) * cos(betta);
1077 g.fPy = g.fE * cos(thetta) * sin(betta);
1078 g.fPz = g.fE * sin(thetta);
1082 //______________________________________________________________________________
1083 //______________________________________________________________________________
1084 //______________________________________________________________________________
1085 //______________________________________________________________________________
1086 //______________________________________________________________________________
1088 ClassImp(AliPHOSgamma)
1090 //______________________________________________________________________________
1092 void AliPHOSgamma::Print(Option_t *)
1094 float mass = fE*fE - fPx*fPx - fPy*fPy - fPz*fPz;
1099 mass = -sqrt(-mass);
1101 printf("XY=(%+7.2f,%+7.2f) (%+7.2f,%+7.2f,%+7.2f;%7.2f) mass=%8.4f Ipart=%2d\n",
1102 fX,fY,fPx,fPy,fPz,fE,mass,fIpart);
1105 //______________________________________________________________________________
1107 AliPHOSgamma &AliPHOSgamma::operator=(const AliPHOSgamma &g)
1120 //______________________________________________________________________________
1121 //______________________________________________________________________________
1122 //______________________________________________________________________________
1123 //______________________________________________________________________________
1124 //______________________________________________________________________________
1126 ClassImp(AliPHOShit)
1128 //______________________________________________________________________________
1130 AliPHOShit::AliPHOShit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
1131 AliHit(shunt, track)
1134 for (i=0;i<5;i++) fVolume[i] = vol[i];
1141 //______________________________________________________________________________