]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOS.cxx
Use gMC and not pMC everywhere
[u/mrichter/AliRoot.git] / PHOS / AliPHOS.cxx
1 ////////////////////////////////////////////////
2 //  Manager and hits classes for set:PHOS     //
3 ////////////////////////////////////////////////
4  
5 // --- ROOT system ---
6 #include "TH1.h"
7 #include "TRandom.h"
8 #include "TFile.h"
9 #include "TTree.h"
10 #include "TBRIK.h"
11 #include "TNode.h"
12 #include "TMath.h"
13
14 // --- Standard library ---
15 #include <stdio.h>
16 #include <string.h>
17 #include <stdlib.h>
18
19 // --- galice header files ---
20 #include "AliPHOS.h"
21 #include "AliRun.h"
22 #include "AliMC.h" 
23 #include "TGeant3.h"
24
25 //______________________________________________________________________________
26
27
28 ClassImp(AliPHOS)
29
30 //______________________________________________________________________________
31
32 AliPHOS::~AliPHOS(void)
33 {
34   delete fHits;                 // 28.12.1998
35   delete fTreePHOS;             // 28.12.1998
36   fCradles->Delete();
37   delete fCradles;
38 }
39
40 //______________________________________________________________________________
41
42 AliPHOS::AliPHOS() :
43          fDebugLevel            (0),
44          fTreePHOS              (NULL),
45          fBranchNameOfCradles   ("AliPHOSCradles"),
46          fTreeName              ("PHOS")
47 {
48    fIshunt   = 0;
49
50   if( NULL==(fCradles=new TObjArray) )
51   {
52     Error("AliPHOS","Can not create fCradles");
53     exit(1);
54   }
55   DefPars();
56 }
57  
58 //______________________________________________________________________________
59
60 AliPHOS::AliPHOS(const char *name, const char *title)
61        : AliDetector            (name,title),
62          fDebugLevel            (0),
63          fTreePHOS              (NULL),
64          fBranchNameOfCradles   ("AliPHOSCradles"),
65          fTreeName              ("PHOS")
66 {
67 //Begin_Html
68 /*
69 <img src="picts/aliphos.gif">
70 */
71 //End_Html
72  
73    fHits   = new TClonesArray("AliPHOShit",  405);
74  
75    fIshunt     =  0;
76
77    SetMarkerColor(kGreen);
78    SetMarkerStyle(2);
79    SetMarkerSize(0.4);
80
81   if( NULL==(fCradles=new TObjArray) ) {
82      Error("AliPHOS","Can not create fCradles");
83      exit(1);
84   }
85   DefPars();
86 }
87
88 //______________________________________________________________________________
89
90 void AliPHOS::DefPars()
91
92       PHOSflags[0]=0;
93       PHOSflags[1]=1;
94       PHOSflags[2]=0;
95       PHOSflags[3]=0;
96       PHOSflags[4]=0;
97       PHOSflags[5]=0;
98       PHOSflags[6]=0;
99       PHOSflags[7]=0;
100       PHOSflags[8]=0;
101       PHOScell[0]=2.2;
102       PHOScell[1]=18.;
103       PHOScell[2]=0.01;
104       PHOScell[3]=0.01;
105       PHOScell[4]=1.0;
106       PHOScell[5]=0.1;
107       PHOScell[6]=0.;
108       PHOScell[7]=0.;
109       PHOScell[8]=0.;
110       PHOSradius=460.;
111       PHOSsize[0]=104;
112       PHOSsize[1]=88;
113       PHOSsize[2]=4;
114       PHOScradlesA=0.;
115       PHOSCPV[0]=1.;
116       PHOSCPV[1]=2.;
117       PHOSCPV[2]=0.;
118       PHOSCPV[3]=0.;
119       PHOSCPV[4]=0.;
120       PHOSCPV[5]=0.;
121       PHOSCPV[6]=0.;
122       PHOSCPV[7]=0.;
123       PHOSCPV[8]=0.;
124       PHOSextra[0]=0.001;
125       PHOSextra[1]=6.95;
126       PHOSextra[2]=4.;
127       PHOSextra[3]=5.;
128       PHOSextra[4]=2.;
129       PHOSextra[5]=0.06;
130       PHOSextra[6]=10.;
131       PHOSextra[7]=3.;
132       PHOSextra[8]=1.;
133       PHOSTXW[0]=209.;
134       PHOSTXW[1]=71.;
135       PHOSTXW[2]=250.;
136       PHOSAIR[0]=206.;
137       PHOSAIR[1]=66.;
138       PHOSAIR[2]=244.;
139       PHOSFTI[0]=214.6;
140       PHOSFTI[1]=80.;
141       PHOSFTI[2]=260.;
142       PHOSFTI[3]=467.;
143 }
144 //______________________________________________________________________________
145
146 void AliPHOS::AddHit(Int_t track, Int_t *vol, Float_t *hits)
147 {
148   TClonesArray &lhits = *fHits;
149   new(lhits[fNhits++]) AliPHOShit(fIshunt,track,vol,hits);
150 }
151  
152 //___________________________________________
153 void AliPHOS::BuildGeometry()
154 {
155
156   TNode *Node, *Top;
157
158   const int kColorPHOS = kRed;
159   //
160   Top=gAlice->GetGeometry()->GetNode("alice");
161
162
163   // PHOS
164   Float_t pphi=12.9399462;
165   new TRotMatrix("rot988","rot988",90,-3*pphi,90,90-3*pphi,0,0);
166   new TRotMatrix("rot989","rot989",90,-  pphi,90,90-  pphi,0,0);
167   new TRotMatrix("rot990","rot990",90,   pphi,90,90+  pphi,0,0);
168   new TRotMatrix("rot991","rot991",90, 3*pphi,90,90+3*pphi,0,0);
169   new TBRIK("S_PHOS","PHOS box","void",107.3,40,130);
170   Top->cd();
171   Node = new TNode("PHOS1","PHOS1","S_PHOS",-317.824921,-395.014343,0,"rot988");
172   Node->SetLineColor(kColorPHOS);
173   fNodes->Add(Node);
174   Top->cd();
175   Node = new TNode("PHOS2","PHOS2","S_PHOS",-113.532333,-494.124908,0,"rot989");
176   fNodes->Add(Node);
177   Node->SetLineColor(kColorPHOS);
178   Top->cd();
179   Node = new TNode("PHOS3","PHOS3","S_PHOS", 113.532333,-494.124908,0,"rot990");
180   Node->SetLineColor(kColorPHOS);
181   fNodes->Add(Node);
182   Top->cd();
183   Node = new TNode("PHOS4","PHOS4","S_PHOS", 317.824921,-395.014343,0,"rot991");
184   Node->SetLineColor(kColorPHOS);
185   fNodes->Add(Node);
186 }
187  
188 //___________________________________________
189 void AliPHOS::CreateMaterials()
190 {
191 // *** DEFINITION OF AVAILABLE PHOS MATERIALS *** 
192
193 // CALLED BY : PHOS_MEDIA 
194 // ORIGIN    : NICK VAN EIJNDHOVEN 
195
196
197
198     Int_t   ISXFLD = gAlice->Field()->Integ();
199     Float_t SXMGMX = gAlice->Field()->Max();
200     
201 // --- The PbWO4 crystals --- 
202     Float_t ax[3] = { 207.19,183.85,16. };
203     Float_t zx[3] = { 82.,74.,8. };
204     Float_t wx[3] = { 1.,1.,4. };
205     Float_t dx    = 8.28;
206 // --- Stainless Steel --- 
207     Float_t as[5] = { 55.847,12.011,51.9961,58.69,28.0855 };
208     Float_t zs[5] = { 26.,6.,24.,28.,14. };
209     Float_t ws[5] = { .6392,8e-4,.2,.14,.02 };
210     Float_t ds    = 8.;
211 // --- The polysterene scintillator (CH) --- 
212     Float_t ap[2] = { 12.011,1.00794 };
213     Float_t zp[2] = { 6.,1. };
214     Float_t wp[2] = { 1.,1. };
215     Float_t dp    = 1.032;
216 // --- Tyvek (CnH2n) 
217     Float_t at[2] = { 12.011,1.00794 };
218     Float_t zt[2] = { 6.,1. };
219     Float_t wt[2] = { 1.,2. };
220     Float_t dt    = .331;
221 // --- Polystyrene foam --- 
222     Float_t af[2] = { 12.011,1.00794 };
223     Float_t zf[2] = { 6.,1. };
224     Float_t wf[2] = { 1.,1. };
225     Float_t df    = .12;
226 //--- Foam thermo insulation (actual chemical composition unknown yet!) ---
227     Float_t ati[2] = { 12.011,1.00794 };
228     Float_t zti[2] = { 6.,1. };
229     Float_t wti[2] = { 1.,1. };
230     Float_t dti    = .1;
231 // --- Textolit (actual chemical composition unknown yet!) --- 
232     Float_t atx[2] = { 12.011,1.00794 };
233     Float_t ztx[2] = { 6.,1. };
234     Float_t wtx[2] = { 1.,1. };
235     Float_t dtx    = 1.83;
236
237     Int_t *idtmed = fIdtmed->GetArray()-699;
238
239     AliMixture(  0, "PbWO4$",          ax, zx, dx, -3, wx);
240     AliMixture(  1, "Polystyrene$",    ap, zp, dp, -2, wp);
241     AliMaterial( 2, "Al$",             26.98, 13., 2.7, 8.9, 999);
242 // ---                                Absorption length^ is ignored --- 
243     AliMixture(  3, "Tyvek$",           at, zt, dt, -2, wt);
244     AliMixture(  4, "Foam$",            af, zf, df, -2, wf);
245     AliMixture(  5, "Stainless Steel$", as, zs, ds, 5, ws);
246     AliMaterial( 6, "Si$",              28.09, 14., 2.33, 9.36, 42.3);
247     AliMixture(  7, "Thermo Insul.$",   ati, zti, dti, -2, wti);
248     AliMixture(  8, "Textolit$",        atx, ztx, dtx, -2, wtx);
249     AliMaterial(99, "Air$",             14.61, 7.3, .001205, 30420., 67500);
250
251     AliMedium(0, "PHOS Xtal    $", 0, 1, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
252     AliMedium(1, "CPV scint.   $", 1, 1, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
253     AliMedium(2, "Al parts     $", 2, 0, ISXFLD, SXMGMX, 10., .1, .1, .001, .001);
254     AliMedium(3, "Tyvek wrapper$", 3, 0, ISXFLD, SXMGMX, 10., .1, .1, .001, .001);
255     AliMedium(4, "Polyst. foam $", 4, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
256     AliMedium(5, "Steel cover  $", 5, 0, ISXFLD, SXMGMX, 10., .1, .1, 1e-4, 1e-4);
257     AliMedium(6, "Si PIN       $", 6, 0, ISXFLD, SXMGMX, 10., .1, .1, .01, .01);
258     AliMedium(7, "Thermo Insul.$", 7, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
259     AliMedium(8, "Textolit     $", 8, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
260     AliMedium(99, "Air          $",99, 0, ISXFLD, SXMGMX, 10., 1., .1, .1, 10);
261
262 // --- Generate explicitly delta rays in the steel cover --- 
263     gMC->Gstpar(idtmed[704], "LOSS", 3.);
264     gMC->Gstpar(idtmed[704], "DRAY", 1.);
265 // --- and in aluminium parts --- 
266     gMC->Gstpar(idtmed[701], "LOSS", 3.);
267     gMC->Gstpar(idtmed[701], "DRAY", 1.);
268 }
269  
270 //______________________________________________________________________________
271
272 void AliPHOS::AddPHOSCradles()
273 {
274   Int_t i;
275   for(i=0;i<GetCradlesAmount();i++) {
276     
277     int n = fCradles->GetEntries();
278     fCradles->Add(new AliPHOSCradle( IsVersion(),            // geometry.
279                                      GetCrystalSideSize    (),
280                                      GetCrystalLength      (),
281                                      GetWrapThickness      (),
282                                      GetAirThickness       (),
283                                      GetPIN_SideSize       (),
284                                      GetPIN_Length         (),
285                                      GetRadius             (),
286                                      GetCPV_Thickness      (),
287                                      GetCPV_PHOS_Distance  (),
288                                      GetNz                 (),
289                                      GetNphi               (),
290                                      GetCradleAngle        (i)));
291     
292     if( n+1 != fCradles->GetEntries() || NULL == fCradles->At(n) )
293       {
294         cout << "  Can not create or add AliPHOSCradle.\n";
295         exit(1);
296       }
297   }
298 }
299
300 //______________________________________________________________________________
301
302 Int_t AliPHOS::DistancetoPrimitive(Int_t , Int_t )
303 {
304    return 9999;
305 }
306  
307 //___________________________________________
308 void AliPHOS::Init()
309 {
310   Int_t i;
311   //
312   printf("\n");
313   for(i=0;i<35;i++) printf("*");
314   printf(" PHOS_INIT ");
315   for(i=0;i<35;i++) printf("*");
316   printf("\n");
317   //
318   // Here the ABSO initialisation code (if any!)
319   for(i=0;i<80;i++) printf("*");
320   printf("\n");
321 }
322
323 //______________________________________________________________________________
324
325 void AliPHOS::MakeBranch(Option_t *)
326 {
327 // ROOT output initialization to ROOT file.
328 // 
329 // AliDetector::MakeBranch()  is always called.
330 //
331 // There will be also special tree "PHOS" with one branch "AliPHOSCradles"
332 // if it was set next flag in the galice card file:
333 //  * PHOSflags:    YES: X<>0   NO: X=0
334 //  * PHOSflags(1) : -----X.  Create branch for TObjArray of AliPHOSCradle
335 //     Examples:
336 //     PHOSflags      1.
337 //     PHOSflags 636301.
338 // In that case special bit CradlesBranch_Bit will be set for AliPHOS
339
340   AliDetector::MakeBranch();
341   
342   int i;
343   float t = GetPHOS_flag(0)/10;
344   i = (int) t;
345   i = (int) ((t-i)*10);
346   if( !i )
347     return;
348
349   SetBit(CradlesBranch_Bit);
350
351   if( NULL==(fTreePHOS=new TTree(fTreeName.Data(),"PHOS events tree")) )
352   {
353     Error("MakeBranch","Can not create TTree");
354     exit(1);
355   }
356
357   if( NULL==fTreePHOS->GetCurrentFile() )
358   {
359     Error("MakeBranch","There is no opened ROOT file");
360     exit(1);
361   }
362
363   // Create a new branch in the current Root Tree.
364
365   if( NULL==fTreePHOS->Branch(fBranchNameOfCradles.Data(),"TObjArray",&fCradles,4000,0) )
366   {
367     Error("MakeBranch","Can not create branch");
368     exit(1);
369   }
370
371   printf("The branch %s has been created\n",fBranchNameOfCradles.Data());
372 }
373
374 //______________________________________________________________________________
375
376 void AliPHOS::SetTreeAddress(void)
377 {
378 // ROOT input initialization.
379 //
380 // AliDetector::SetTreeAddress()  is always called.
381 //
382 // If CradlesBranch_Bit is set (see AliPHOS::MakeBranch) than fTreePHOS is
383 // initilized.
384
385   AliDetector::SetTreeAddress();
386
387   if( !TestBit(CradlesBranch_Bit) )
388     return;
389
390   if( NULL==(fTreePHOS=(TTree*)gDirectory->Get((char*)(fTreeName.Data()))  ) )
391   {
392     Error("Can not find Tree \"%s\"\n",fTreeName.Data());
393     exit(1);
394   }
395
396   TBranch *branch = fTreePHOS->GetBranch(fBranchNameOfCradles.Data());
397   if( NULL==branch )
398   {
399     Error("SetTreeAddress","Can not find branch %s in TTree:%s",fBranchNameOfCradles.Data(),fTreeName.Data());
400     exit(1);
401   }
402
403   branch->SetAddress(&fCradles);
404 }
405
406 //______________________________________________________________________________
407
408 AliPHOSCradle *AliPHOS::GetCradleOfTheParticle(const TVector3 &p,const TVector3 &v) const
409 {
410 // For a given direction 'p' and source point 'v' returns pointer to AliPHOSCradle
411 // in that direction or NULL if AliPHOSCradle was not found.
412
413   for( int m=0; m<fCradles->GetEntries(); m++ )
414   {
415     AliPHOS *PHOS = (AliPHOS *)this;     // Removing 'const'...
416     AliPHOSCradle *cradle = (AliPHOSCradle *)PHOS->fCradles->operator[](m);
417
418     float x,y,l;
419     const float d = cradle->GetRadius()-cradle->GetCPV_PHOS_Distance()-cradle->GetCPV_Thikness();
420     cradle->GetXY(p,v,d,x,y,l);
421
422     if( l>0 && TMath::Abs(x)<cradle->GetNz  ()*cradle->GetCellSideSize()/2 
423             && TMath::Abs(y)<cradle->GetNphi()*cradle->GetCellSideSize()/2 )
424       return cradle;
425   }
426
427   return NULL;
428 }
429
430 //______________________________________________________________________________
431
432 void AliPHOS::Reconstruction(Float_t signal_step, UInt_t min_signal_reject)
433 {
434 // Call AliPHOSCradle::Reconstruction(Float_t signal_step, UInt_t min_signal_reject)
435 // for all AliPHOSCradles.
436
437   for( int i=0; i<fCradles->GetEntries(); i++ )
438     GetCradle(i).Reconstruction(signal_step,min_signal_reject);
439 }
440
441 //______________________________________________________________________________
442
443 void AliPHOS::ResetDigits(void)
444 {
445   AliDetector::ResetDigits();
446
447   for( int i=0; i<fCradles->GetEntries(); i++ )
448     ((AliPHOSCradle*)(*fCradles)[i]) -> Clear();
449 }
450
451 //______________________________________________________________________________
452
453 void AliPHOS::FinishEvent(void)
454 {
455 // Called at the end of each 'galice' event.
456
457   if( NULL!=fTreePHOS )
458     fTreePHOS->Fill();
459 }
460
461 //______________________________________________________________________________
462
463 void AliPHOS::FinishRun(void)
464 {
465 }
466
467 //______________________________________________________________________________
468
469 void AliPHOS::Print(Option_t *opt)
470 {
471 // Print PHOS information.
472 // For each AliPHOSCradle the function AliPHOSCradle::Print(opt) is called.
473
474   AliPHOS &PHOS = *(AliPHOS *)this;     // Removing 'const'...
475
476   for( int i=0; i<fCradles->GetEntries(); i++ )
477   {
478     printf("PHOS cradle %d from %d\n",i+1, fCradles->GetEntries());
479     PHOS.GetCradle(i).Print(opt);
480     printf( "---------------------------------------------------\n");
481   }
482 }
483
484 //______________________________________________________________________________
485 void AliPHOS::SetFlags(Float_t p1,Float_t p2,Float_t p3,Float_t p4,
486                        Float_t p5,Float_t p6,Float_t p7,Float_t p8,Float_t p9)
487 {
488   PHOSflags[0]=p1;
489   PHOSflags[1]=p2;
490   PHOSflags[2]=p3;
491   PHOSflags[3]=p4;
492   PHOSflags[4]=p5;
493   PHOSflags[5]=p6;
494   PHOSflags[6]=p7;
495   PHOSflags[7]=p8;
496   PHOSflags[8]=p9;
497 }
498
499 //______________________________________________________________________________
500 void AliPHOS::SetCell(Float_t p1,Float_t p2,Float_t p3,Float_t p4,
501                        Float_t p5,Float_t p6,Float_t p7,Float_t p8,Float_t p9)
502 {
503   PHOScell[0]=p1;
504   PHOScell[1]=p2;
505   PHOScell[2]=p3;
506   PHOScell[3]=p4;
507   PHOScell[4]=p5;
508   PHOScell[5]=p6;
509   PHOScell[6]=p7;
510   PHOScell[7]=p8;
511   PHOScell[8]=p9;
512 }
513
514 //______________________________________________________________________________
515 void AliPHOS::SetRadius(Float_t radius)
516 {
517    PHOSradius=radius;
518 }
519
520 //______________________________________________________________________________
521 void AliPHOS::SetCradleSize(Int_t nz, Int_t nphi, Int_t ncradles)
522 {
523    PHOSsize[0]=nz;
524    PHOSsize[1]=nphi;
525    PHOSsize[2]=ncradles;
526 }
527
528 //______________________________________________________________________________
529 void AliPHOS::SetCradleA(Float_t angle)
530 {
531    PHOScradlesA=angle;
532 }
533
534 //______________________________________________________________________________
535 void AliPHOS::SetCPV(Float_t p1,Float_t p2,Float_t p3,Float_t p4,
536                      Float_t p5,Float_t p6,Float_t p7,Float_t p8,Float_t p9)
537 {
538    PHOSCPV[0] = p1;
539    PHOSCPV[1] = p2;
540    PHOSCPV[2] = p3;
541    PHOSCPV[3] = p4;
542    PHOSCPV[4] = p5;
543    PHOSCPV[5] = p6;
544    PHOSCPV[6] = p7;
545    PHOSCPV[7] = p8;
546    PHOSCPV[8] = p9;
547 }
548
549 //______________________________________________________________________________
550 void AliPHOS::SetExtra(Float_t p1,Float_t p2,Float_t p3,Float_t p4,
551                        Float_t p5,Float_t p6,Float_t p7,Float_t p8,Float_t p9)
552 {
553    PHOSextra[0] = p1;
554    PHOSextra[1] = p2;
555    PHOSextra[2] = p3;
556    PHOSextra[3] = p4;
557    PHOSextra[4] = p5;
558    PHOSextra[5] = p6;
559    PHOSextra[6] = p7;
560    PHOSextra[7] = p8;
561    PHOSextra[8] = p9;
562 }
563
564 //______________________________________________________________________________
565 void AliPHOS::SetTextolitWall(Float_t dx, Float_t dy, Float_t dz)
566 {
567    PHOSTXW[0] = dx;
568    PHOSTXW[1] = dy;
569    PHOSTXW[2] = dz;
570 }
571
572 //______________________________________________________________________________
573 void AliPHOS::SetInnerAir(Float_t dx, Float_t dy, Float_t dz)
574 {
575    PHOSAIR[0] = dx;
576    PHOSAIR[1] = dy;
577    PHOSAIR[2] = dz;
578 }
579
580 //______________________________________________________________________________
581 void AliPHOS::SetFoam(Float_t dx, Float_t dy, Float_t dz, Float_t dr)
582 {
583    PHOSFTI[0] = dx;
584    PHOSFTI[1] = dy;
585    PHOSFTI[2] = dz;
586    PHOSFTI[3] = dr;
587 }
588
589 ClassImp(AliPHOSCradle)
590
591 //______________________________________________________________________________
592
593 AliPHOSCradle::AliPHOSCradle(void) {}
594
595 //______________________________________________________________________________
596
597 AliPHOSCradle::AliPHOSCradle( int   Geometry           ,
598                               float CrystalSideSize    ,
599                               float CrystalLength      ,
600                               float WrapThickness      ,
601                               float AirThickness       ,
602                               float PIN_SideSize       ,
603                               float PIN_Length         ,
604                               float Radius             ,
605                               float CPV_Thickness      ,
606                               float CPV_PHOS_Distance  ,
607                               int   Nz                 ,
608                               int   Nphi               ,
609                               float Angle              ) :
610     fGeometry                   (Geometry),
611 //  fCellEnergy                 (),
612 //  fChargedTracksInPIN         (),
613 //  fCPV_hitsX                  (),
614 //  fCPV_hitsY                  (),
615     fCrystalSideSize            (CrystalSideSize),
616     fCrystalLength              (CrystalLength),
617     fWrapThickness              (WrapThickness),
618     fAirThickness               (AirThickness),
619     fPIN_SideSize               (PIN_SideSize),
620     fPIN_Length                 (PIN_Length),
621     fRadius                     (Radius),
622     fCPV_PHOS_Distance          (CPV_PHOS_Distance),
623     fCPV_Thickness              (CPV_Thickness),
624     fNz                         (Nz),
625     fNphi                       (Nphi),
626     fPhi                        (Angle)
627 {
628         fCellEnergy         = TH2F("CellE","Energy deposition in a cells",fNz,0,fNz,fNphi,0,fNphi);
629         fCellEnergy           .SetDirectory(0);
630         fChargedTracksInPIN = TH2S("PINCtracks","Amount of charged tracks in PIN",fNz,0,fNz,fNphi,0,fNphi);
631         fChargedTracksInPIN   .SetDirectory(0);
632 }
633
634 //______________________________________________________________________________
635
636 AliPHOSCradle::~AliPHOSCradle(void)        // 28.12.1998
637 {
638   fGammasReconstructed.Delete();
639   fParticles          .Delete();
640 }
641
642 //______________________________________________________________________________
643
644 void AliPHOSCradle::Clear(Option_t *)
645 {
646 // Clear digit. information.
647
648   fCellEnergy              .Reset();
649   fChargedTracksInPIN      .Reset();
650   GetParticles()           .Delete();
651   GetParticles()           .Compress();
652   GetGammasReconstructed() .Delete();
653   GetGammasReconstructed() .Compress();
654
655   fCPV_hitsX.Set(0);
656   fCPV_hitsY.Set(0);
657 }
658
659 //______________________________________________________________________________
660
661 void AliPHOSCradle::AddCPVHit(float x,float y)
662 {
663 // Add this hit to the hits list in CPV detector.
664
665   TArrayF a(fCPV_hitsX.GetSize()+1);
666   
667   memcpy(a.GetArray(),fCPV_hitsX.GetArray(),sizeof(Float_t)*fCPV_hitsX.GetSize());
668   a[fCPV_hitsX.GetSize()] = x;
669   fCPV_hitsX = a;
670
671   // It must be:   fCPV_hitsX.GetSize() == fCPV_hitsY.GetSize()
672
673   memcpy(a.GetArray(),fCPV_hitsY.GetArray(),sizeof(Float_t)*fCPV_hitsY.GetSize());
674   a[fCPV_hitsY.GetSize()] = y;
675   fCPV_hitsY = a;
676 }
677
678 //______________________________________________________________________________
679
680 void AliPHOSCradle::GetXY(const TVector3 &p,const TVector3 &v,float R,float &x,float &y,float &l) const
681 {
682 // This function calculates hit position (x,y) in the CRADLE cells plain from particle in
683 // the direction given by 'p' (not required to be normalized) and start point
684 // given by 3-vector 'v'. So the particle trajectory is   t(l) = v + p*l
685 // were 'l' is a number (distance from 'v' to CRADLE cells plain) and 't' is resulting
686 // three-vector of trajectory point.
687 // 
688 // After the call to this function user should test that l>=0 (the particle HITED the
689 // plain) and (x,y) are in the region of CRADLE:
690 // 
691 // Example:
692 //   AliPHOSCradle cradle(......);
693 //   TVector3 p(....), v(....);
694 //   Float_t x,y,l;
695 //   cradle.GetXY(p,v,x,y,l);
696 //   if( l<0 || TMath::Abs(x)>cradle.GetNz()  *cradle.GetCellSideSize()/2
697 //           || TMath::Abs(y)>cradle.GetNphi()*cradle.GetCellSideSize()/2 )
698 //     cout << "Outside the CRADLE.\n";
699
700   // We have to create three vectors:
701   //    s  - central point on the PHOS surface
702   //    n1 - first vector in CRADLE plain
703   //    n2 - second vector in CRADLE plain
704   // This three vectors are orthonormalized.
705
706   double phi = fPhi/180*TMath::Pi();
707   TVector3        n1(   0.0      ,   0.0      , 1.0 ),   // Z direction (X)
708                   n2(  -sin(phi) ,   cos(phi) , 0 ),   // around beam (Y)
709                   s ( R*cos(phi) , R*sin(phi) , 0 );   // central point
710
711   const double l1_min = 1e-2;
712   double l1,
713          p_n1 = p*n1,        // * - scalar product.
714          p_n2 = p*n2,
715          v_n1 = v*n1,
716          v_n2 = v*n2,
717          s_n1 = s*n1, // 0
718          s_n2 = s*n2; // 0
719   
720   if      ( TMath::Abs(l1=p.X()-n1.X()*p_n1-n2.X()*p_n2)>l1_min )
721     { l = (-v.X()+s.X()+n1.X()*(v_n1-s_n1)+n2.X()*(v_n2-s_n2))/l1; }
722   else if ( TMath::Abs(l1=p.Y()-n1.Y()*p_n1-n2.Y()*p_n2)>l1_min )
723     { l = (-v.Y()+s.Y()+n1.Y()*(v_n1-s_n1)+n2.Y()*(v_n2-s_n2))/l1; }
724   else if ( TMath::Abs(l1=p.Z()-n1.Z()*p_n1-n2.Z()*p_n2)>l1_min )
725     { l = (-v.Z()+s.Z()+n1.Z()*(v_n1-s_n1)+n2.Z()*(v_n2-s_n2))/l1; }
726
727 //         double lx = (-v.X()+s.X()+n1.X()*(v.dot(n1)-s.dot(n1))+n2.X()*(v.dot(n2)-s.dot(n2)))/
728 //                     (p.X()-n1.X()*p.dot(n1)-n2.X()*p.dot(n2)),
729 //                ly = (-v.Y()+s.Y()+n1.Y()*(v.dot(n1)-s.dot(n1))+n2.Y()*(v.dot(n2)-s.dot(n2)))/
730 //                     (p.Y()-n1.Y()*p.dot(n1)-n2.Y()*p.dot(n2)),
731 //                lz = (-v.Z()+s.Z()+n1.Z()*(v.dot(n1)-s.dot(n1))+n2.Z()*(v.dot(n2)-s.dot(n2)))/
732 //                     (p.Z()-n1.Z()*p.dot(n1)-n2.Z()*p.dot(n2));
733 //         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));
734 //         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));
735 //         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));
736 //         cout.form("lx,ly,lz =   %g,%g,%g\n",lx,ly,lz);
737
738   x = p_n1*l + v_n1 - s_n1;
739   y = p_n2*l + v_n2 - s_n2;
740 }
741
742 //______________________________________________________________________________
743
744 void AliPHOSCradle::Print(Option_t *opt)
745 {
746 // Print AliPHOSCradle information.
747 // 
748 // options:  'd' - print energy deposition for EVERY cell
749 //           'p' - print particles list that hit the cradle
750 //           'r' - print list of reconstructed particles
751
752   AliPHOSCradle *cr = (AliPHOSCradle *)this;     // Removing 'const'...
753
754   printf("AliPHOSCradle:  Nz=%d  Nphi=%d, fPhi=%f, E=%g, CPV hits amount = %d\n",fNz,fNphi,fPhi,
755        cr->fCellEnergy.GetSumOfWeights(),fCPV_hitsX.GetSize());
756
757   if( NULL!=strchr(opt,'d') )
758   {
759     printf("\n\nCells Energy (in MeV):\n\n   |");
760     for( int x=0; x<fNz; x++ )
761       printf(" %4d|",x+1);
762     printf("\n");
763
764     for( int y=fNphi-1; y>=0; y-- )
765     {
766       printf("%3d|",y+1);
767       for( int x=0; x<fNz; x++ )
768         printf("%6d",(int)(cr->fCellEnergy.GetBinContent(cr->fCellEnergy.GetBin(x,y))*1000));
769       printf("\n");
770     }
771     printf("\n");
772   }
773
774   if( NULL!=strchr(opt,'p') )
775   {
776     printf("This cradle was hit by %d particles\n",
777          ((AliPHOSCradle*)this)->GetParticles().GetEntries());
778     TObjArray &p=((AliPHOSCradle*)this)->GetParticles();
779     for( int i=0; i<p.GetEntries(); i++ )
780       ((AliPHOSgamma*)(p[i]))->Print();
781   }
782
783   if( NULL!=strchr(opt,'p') )
784   {
785     printf("Amount of reconstructed gammas is %d\n",
786          ((AliPHOSCradle*)this)->GetGammasReconstructed().GetEntries());
787
788     TObjArray &p=((AliPHOSCradle*)this)->GetGammasReconstructed();
789     for( int i=0; i<p.GetEntries(); i++ )
790       ((AliPHOSgamma*)(p[i]))->Print();
791   }
792 }
793
794 //______________________________________________________________________________
795
796 void AliPHOSCradle::Distortion(const TH2F *Noise, const TH2F *Stochastic, const TH2F *Calibration)
797 {
798 // This function changes histogram of cell energies fCellEnergy on the base of input
799 // histograms Noise, Stochastic, Calibration. The histograms must have
800 // size Nz x Nphi. 
801
802   //////////////////////////////////
803   // Testing the histograms size. //
804   //////////////////////////////////
805   
806   if( fNz!=fCellEnergy.GetNbinsX() || fNphi!=fCellEnergy.GetNbinsY() )
807   {
808     printf      ("Bad size of CellEnergy!   Must be:   Nz x Nphi = %d x %d\n"
809                  "but size of CellEnergy is:  %d x %d\n",
810                  fNz,fNphi,fCellEnergy.GetNbinsX(),fCellEnergy.GetNbinsY());
811     exit(1);
812   }
813
814   if( fNz!=fChargedTracksInPIN.GetNbinsX() || fNphi!=fChargedTracksInPIN.GetNbinsY() )
815   {
816     printf      ("Bad size of ChargedTracksInPIN!   Must be:   Nz x Nphi = %d x %d\n"
817                  "but size of ChargedTracksInPIN is:  %d x %d\n",
818                  fNz,fNphi,fChargedTracksInPIN.GetNbinsX(),fChargedTracksInPIN.GetNbinsY());
819     exit(1);
820   }
821
822   if( NULL!=Noise && (fNz!=Noise->GetNbinsX() || fNphi!=Noise->GetNbinsX()) )
823   {
824     printf      ("Bad size of Noise!   Must be:   Nz x Nphi = %d x %d\n"
825                  "but size of Noise is:  %d x %d\n",
826                  fNz,fNphi,fChargedTracksInPIN.GetNbinsX(),fChargedTracksInPIN.GetNbinsY());
827     exit(1);
828   }
829
830   if( NULL!=Stochastic && (fNz!=Stochastic->GetNbinsX() || fNphi!=Stochastic->GetNbinsX()) )
831   {
832     printf      ("Bad size of Stochastic!   Must be:   Nz x Nphi = %d x %d\n"
833                  "but size of Stochastic is:  %d x %d\n",
834                  fNz,fNphi,fChargedTracksInPIN.GetNbinsX(),fChargedTracksInPIN.GetNbinsY());
835     exit(1);
836   }
837
838   if( NULL!=Calibration && (fNz!=Calibration->GetNbinsX() || fNphi!=Calibration->GetNbinsX()) )
839   {
840     printf      ("Bad size of Calibration!   Must be:   Nz x Nphi = %d x %d\n"
841                  "but size of Calibration is:  %d x %d\n",
842                  fNz,fNphi,fChargedTracksInPIN.GetNbinsX(),fChargedTracksInPIN.GetNbinsY());
843     exit(1);
844   }
845
846   ////////////////////
847   // Do distortion! //
848   ////////////////////
849
850   for( int y=0; y<fNphi; y++ )
851     for( int x=0; x<fNz; x++ )
852     {
853       const int n = fCellEnergy.GetBin(x,y);   // Bin number
854       static TRandom r;
855     
856       Float_t   E_old=fCellEnergy.GetBinContent(n),   E_new=E_old;
857
858       if( NULL!=Stochastic )
859         E_new   = r.Gaus(E_old,sqrt(E_old)*GetDistortedValue(Stochastic,n));
860
861       if( NULL!=Calibration )
862         E_new  *=  GetDistortedValue(Calibration,n);
863
864       if( NULL!=Noise )
865         E_new  +=  GetDistortedValue(Noise,n);
866
867       fCellEnergy.SetBinContent(n,E_new);
868     }
869 }
870
871 ////////////////////////////////////////////////////////////////////////////////
872
873 TH2F* AliPHOSCradle::CreateHistForDistortion(const char *name, const char *title,
874                                              Int_t Nx, Int_t Ny,
875                                              Float_t MU_mu,    Float_t MU_sigma,
876                                              Float_t SIGMA_mu, Float_t SIGMA_sigma)
877 {
878 // Create (new TH2F(...)) histogram with information (for every bin) that will
879 // be used for VALUE creation.
880 // Two values will be created for each bin:
881 // MU    = TRandom::Gaus(MU_mu,MU_sigma)
882 // and
883 // SIGMA = TRandom::Gaus(SIGMA_mu,SIGMA_sigma)
884 // The VALUE in a particluar bin will be equal
885 // VALUE = TRandom::Gaus(MU,SIGMA)
886 // 
887 // Do not forget to delete the histogram at the end of the work.
888
889   TH2F *h = new TH2F( name,title, Nx,1,Nx, Ny,1,Ny );
890   if( h==NULL )
891   {
892     Error("CreateHistForDistortion","Can not create the histogram");
893     exit(1);
894   }
895   h->SetDirectory(0);
896
897   for( int y=0; y<Ny; y++ )
898     for( int x=0; x<Nx; x++ )
899     {
900       const int n = h->GetBin(x,y);
901       h->SetBinContent(n,r.Gaus(   MU_mu,   MU_sigma));
902       h->SetBinError  (n,r.Gaus(SIGMA_mu,SIGMA_sigma));
903     }
904
905   return h;
906 }
907
908 ////////////////////////////////////////////////////////////////////////////////
909
910 Float_t AliPHOSCradle::GetDistortedValue(const TH2F *h, UInt_t n)
911 {
912   return r.Gaus(((TH2F*)h)->GetBinContent(n),n);
913 }
914
915 ////////////////////////////////////////////////////////////////////////////////
916 //______________________________________________________________________________
917
918 #ifdef WIN32
919   #define common_for_event_storing COMMON_FOR_EVENT_STORING
920 #else
921   #define common_for_event_storing common_for_event_storing_
922 #endif
923
924 extern "C" struct
925 {
926   enum { crystals_matrix_amount_max=4, crystals_in_matrix_amount_max=40000 };
927
928   // Event-independent information
929   UShort_t      crystals_matrix_amount_PHOS,
930                 crystal_matrix_type,
931                 amount_of_crystals_on_Z,
932                 amount_of_crystals_on_PHI;
933   Float_t       radius,
934                 crystal_size,
935                 crystal_length,
936                 matrix_coordinate_Z             [crystals_matrix_amount_max],
937                 matrix_coordinate_PHI           [crystals_matrix_amount_max];
938   UInt_t        event_number;
939   UShort_t      crystals_amount_with_amplitudes [crystals_matrix_amount_max],
940                 crystals_amplitudes_Iad         [crystals_matrix_amount_max]
941                                                 [crystals_in_matrix_amount_max][2];
942 } common_for_event_storing;
943
944 //       integer*4 crystals_amount_max,crystals_in_matrix_amount_max,
945 //      +          crystals_matrix_amount_max
946 //       parameter (crystals_matrix_amount_max=4)
947 //       parameter (crystals_in_matrix_amount_max=40000)
948 //       parameter (crystals_amount_max =crystals_matrix_amount_max*
949 //      +                                crystals_in_matrix_amount_max)
950 // 
951 // * All units are in GeV, cm, radian
952 //       real       crystal_amplitudes_unit, radius_unit,
953 //      +           crystal_size_unit, crystal_length_unit,
954 //      +           matrix_coordinate_Z_unit, matrix_coordinate_PHI_unit
955 //       integer    crystal_amplitudes_in_units_min
956 //       parameter (crystal_amplitudes_in_units_min        = 1)
957 //       parameter (crystal_amplitudes_unit                = 0.001 ) ! 1.0  MeV
958 //       parameter (radius_unit                            = 0.1   ) ! 0.1  cm
959 //       parameter (crystal_size_unit                      = 0.01  ) ! 0.01 cm
960 //       parameter (crystal_length_unit                    = 0.01  ) ! 0.01 cm
961 //       parameter (matrix_coordinate_Z_unit               = 0.1   ) ! 0.1  cm
962 //       parameter (matrix_coordinate_PHI_unit             = 1e-4  ) ! 1e-4 radian
963 // 
964 //       integer*2 crystals_matrix_amount_PHOS, crystal_matrix_type,
965 //      +          amount_of_crystals_on_Z, amount_of_crystals_on_PHI,
966 //      +          crystals_amount_with_amplitudes, crystals_amplitudes_Iad
967 //       integer*4 event_number
968 // 
969 //       real      radius, crystal_size, crystal_length,
970 //      +          matrix_coordinate_Z, matrix_coordinate_PHI
971 // 
972 //       real      crystals_amplitudes, crystals_energy_total
973 //       integer   event_file_unit_number
974 // 
975 //       common /common_for_event_storing/
976 //      + ! Event-independent information
977 //      +        crystals_matrix_amount_PHOS,
978 //      +        crystal_matrix_type,
979 //      +        amount_of_crystals_on_Z,
980 //      +        amount_of_crystals_on_PHI,
981 //      +        radius,
982 //      +        crystal_size,
983 //      +        crystal_length,
984 //      +        matrix_coordinate_Z     (crystals_matrix_amount_max),
985 //      +        matrix_coordinate_PHI   (crystals_matrix_amount_max),
986 //      +
987 //      + ! Event-dependent information
988 //      +        event_number,
989 //      +        crystals_amount_with_amplitudes
990 //      +                                (crystals_matrix_amount_max),
991 //      +        crystals_amplitudes_Iad (2,crystals_in_matrix_amount_max,
992 //      +                                 crystals_matrix_amount_max),
993 //      +        
994 //      + ! These information don't store in data file
995 //      +        crystals_amplitudes     (crystals_amount_max),
996 //      +        crystals_energy_total,
997 //      +        event_file_unit_number
998
999
1000 //      parameter (NGp=1000,nsps=10,nvertmax=1000)
1001 //         COMMON /GAMMA/KG,MW(ngp),ID(ngp),JD(ngp),E(ngp),E4(ngp),
1002 //      ,  XW(ngp),YW(ngp),ES(nsps,ngp),ET(nsps,ngp),ISsd(ngp),
1003 //      ,  IGDEV(ngp),ZGDEV(ngp),sigexy(3,ngp),Emimx(2,nsps,ngp),
1004 //      ,  kgfix,igfix(ngp),cgfix(3,ngp),sgfix(3,ngp),hiw(ngp),
1005 //      ,  wsw(nsps,ngp),h1w(ngp),h0w(ngp),raxay(5,ngp),
1006 //      ,  sigmaes0(nsps,ngp),dispeces(nsps,ngp),
1007 //      ,  igamvert(ngp)
1008
1009
1010 #ifdef WIN32
1011 #define rcgamma RCGAMMA
1012 #else
1013 #define rcgamma rcgamma_
1014 #endif
1015
1016 extern "C" struct
1017 {
1018   enum {NGP=1000, nsps=10, nvertmax=1000};
1019   int   recons_gammas_amount, mw[NGP],ID[NGP],JD[NGP];
1020   float E[NGP], E4[NGP], XW[NGP], YW[NGP], ES[NGP][nsps],ET[NGP][nsps],ISsd[NGP],
1021         igdev[NGP],Zgdev[NGP];
1022 //      sigexy(3,ngp),Emimx(2,nsps,ngp),
1023 //   ,  kgfix,igfix(ngp),cgfix(3,ngp),sgfix(3,ngp),hiw(ngp),
1024 //   ,  wsw(nsps,ngp),h1w(ngp),h0w(ngp),raxay(5,ngp),
1025 //   ,  sigmaes0(nsps,ngp),dispeces(nsps,ngp),
1026 //   ,  igamvert(ngp)
1027 } rcgamma;
1028
1029 #ifdef WIN32
1030 #define reconsfirst RECONSFIRST
1031 #define type_of_call _stdcall
1032 #else
1033 #define reconsfirst reconsfirst_
1034 #define type_of_call
1035 #endif
1036
1037 extern "C" void type_of_call reconsfirst(const float &,const float &);
1038
1039 void AliPHOSCradle::Reconstruction(Float_t signal_step, UInt_t min_signal_reject)
1040 {
1041 // Call of PHOS reconstruction program.
1042 // signal_step=0.001  GeV (1MeV)
1043 // min_signal_reject = 15 or 30 MeV
1044
1045
1046   common_for_event_storing.event_number                       = 0;  // We do not know event number?
1047   common_for_event_storing.crystals_matrix_amount_PHOS        = 1;
1048   common_for_event_storing.crystal_matrix_type                = 1; // 1 - rectangular
1049   common_for_event_storing.amount_of_crystals_on_Z            = fNz;
1050   common_for_event_storing.amount_of_crystals_on_PHI          = fNphi;
1051
1052   common_for_event_storing.radius                             = fRadius;
1053   common_for_event_storing.crystal_size                       = GetCellSideSize();
1054   common_for_event_storing.crystal_length                     = fCrystalLength;
1055
1056   common_for_event_storing.matrix_coordinate_Z            [0] = 0;
1057   common_for_event_storing.matrix_coordinate_PHI          [0] = fPhi;
1058
1059   #define  k    common_for_event_storing.crystals_amount_with_amplitudes[0] 
1060   k=0;
1061
1062   for( int y=0; y<fNphi; y++ )
1063     for( int x=0; x<fNz; x++ )
1064     {
1065       UInt_t    n       = fCellEnergy.GetBin(x,y);
1066       UInt_t    signal  = (int) (fCellEnergy.GetBinContent(n)/signal_step);
1067       if( signal>=min_signal_reject )
1068       {
1069         common_for_event_storing.crystals_amplitudes_Iad[0][k][0] = signal;
1070         common_for_event_storing.crystals_amplitudes_Iad[0][k][1] = x + y*fNz;
1071         k++;
1072       }
1073     }
1074   #undef  k
1075
1076   GetGammasReconstructed().Delete();
1077   GetGammasReconstructed().Compress();
1078
1079   const float   stochastic_term   = 0.03,        // per cents over sqrt(E);  E in GeV
1080                 electronic_noise  = 0.01;        // GeV
1081   reconsfirst(stochastic_term,electronic_noise); // Call of reconstruction program.
1082
1083   for( int i=0; i<rcgamma.recons_gammas_amount; i++ )
1084   {
1085 //     new (GetGammasReconstructed().UncheckedAt(i) ) AliPHOSgamma;
1086 //     AliPHOSgamma &g = *(AliPHOSgamma*)(GetGammasReconstructed().UncheckedAt(i));
1087
1088     AliPHOSgamma *gggg = new AliPHOSgamma;
1089     if( NULL==gggg )
1090     {
1091       Error("Reconstruction","Can not create AliPHOSgamma");
1092       exit(1);
1093     }
1094
1095     GetGammasReconstructed().Add(gggg);
1096     AliPHOSgamma &g=*gggg;
1097     
1098     Float_t thetta, alpha, betta, R=fRadius+rcgamma.Zgdev[i]/10;
1099
1100     g.fX      = rcgamma.YW[i]/10;
1101     g.fY      = rcgamma.XW[i]/10;
1102     g.fE      = rcgamma.E [i];
1103
1104     thetta      = atan(g.fX/R);
1105
1106     alpha = atan(g.fY/R);
1107     betta = fPhi/180*TMath::Pi() + alpha;
1108
1109     g.fPx = g.fE * cos(thetta) * cos(betta);
1110     g.fPy = g.fE * cos(thetta) * sin(betta);
1111     g.fPz = g.fE * sin(thetta);
1112   }
1113 }
1114
1115 //______________________________________________________________________________
1116 //______________________________________________________________________________
1117 //______________________________________________________________________________
1118 //______________________________________________________________________________
1119 //______________________________________________________________________________
1120
1121 ClassImp(AliPHOSgamma)
1122
1123 //______________________________________________________________________________
1124
1125 void AliPHOSgamma::Print(Option_t *)
1126 {
1127   float mass = fE*fE - fPx*fPx - fPy*fPy - fPz*fPz;
1128
1129   if( mass>=0 )
1130     mass =  sqrt( mass);
1131   else
1132     mass = -sqrt(-mass);
1133
1134   printf("XY=(%+7.2f,%+7.2f)  (%+7.2f,%+7.2f,%+7.2f;%7.2f)  mass=%8.4f  Ipart=%2d\n",
1135           fX,fY,fPx,fPy,fPz,fE,mass,fIpart);
1136 }
1137
1138 //______________________________________________________________________________
1139
1140 AliPHOSgamma &AliPHOSgamma::operator=(const AliPHOSgamma &g)
1141 {
1142   fX           = g.fX;
1143   fY           = g.fY;
1144   fE           = g.fE;
1145   fPx          = g.fPx;
1146   fPy          = g.fPy;
1147   fPz          = g.fPz;
1148   fIpart       = g.fIpart;
1149
1150   return *this;
1151 }
1152
1153 //______________________________________________________________________________
1154 //______________________________________________________________________________
1155 //______________________________________________________________________________
1156 //______________________________________________________________________________
1157 //______________________________________________________________________________
1158
1159 ClassImp(AliPHOShit)
1160
1161 //______________________________________________________________________________
1162
1163 AliPHOShit::AliPHOShit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
1164 AliHit(shunt, track)
1165 {
1166    Int_t i;
1167    for (i=0;i<5;i++) fVolume[i] = vol[i];
1168    fX       = hits[0];
1169    fY       = hits[1];
1170    fZ       = hits[2];
1171    fELOS    = hits[3];
1172 }
1173  
1174 //______________________________________________________________________________