FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndPhoGunShort Class Reference

#include <PndPhoGunShort.h>

Inheritance diagram for PndPhoGunShort:

Public Member Functions

 PndPhoGunShort ()
 
 PndPhoGunShort (Int_t verbose)
 
virtual ~PndPhoGunShort ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *option)
 
virtual void Finish ()
 
void SetOutputFile (TString outName="lutnode.root")
 
void SetEVDepth (Float_t EVdepth=30.)
 
void SetNmcp (Float_t nmcp=5.)
 
void SetShiftedPix (Bool_t spix=kFALSE)
 

Private Member Functions

void ProcessPhotonHit ()
 
void InitLut ()
 
void SetDefaultParameters ()
 
Double_t InBarCoordSyst (TVector3, TVector3 *, TVector3 *, TVector3 *, TVector3 *)
 
Double_t FindReflectionType (Double_t, Double_t, Double_t, TString)
 

Private Attributes

Int_t fDetectorID
 
TClonesArray * fBarPointArray
 
TClonesArray * fPDPointArray
 
TClonesArray * fPDHitArray
 
TClonesArray * fMCArray
 
TClonesArray * fEVPointArray
 
TClonesArray * fDigiArray
 
TList * fHistoList
 
TFile * froot
 
PndDrcDigifDigi
 
PndGeoDrcfGeo
 Basic geometry data of barrel DRC. More...
 
PndGeoHandlingfGeoH
 
Int_t fVerbose
 
Int_t nevents
 
Int_t ambiguity
 
TString fOutputName
 
TString fTableName
 
Double_t ftilt
 
Int_t fBarId
 
Double_t fR
 
Double_t fRBottom
 
Double_t fHThick
 
Double_t fBboxNum
 
Double_t fPipehAngle
 
Double_t fBarBoxGap
 
Double_t fLength
 
Double_t fDphi
 
Double_t fEVlen
 
Double_t fEVdz
 
Double_t fpi
 
Double_t fEVdrop
 
Double_t fPixelSize
 
Float_t fNmcp
 
Bool_t fShiftPix
 
Double_t fNoDD
 
Double_t fNoU0
 
Double_t fNoU1
 
Double_t fNoU2
 
Double_t fNoU3
 
Double_t fNoB
 
Double_t fNoUB
 
Double_t fNoBU0
 
Double_t fNoBU1
 
Double_t fNoBU2
 
Double_t fNoUU0
 
Double_t fNoUU1
 
Double_t fNoUU2
 
Double_t fNoUU3
 
Double_t fNoUUU0
 
Double_t fNoUUU1
 
Double_t fNoUUU2
 
Double_t fNoUUU3
 
Double_t fNoUUU4
 
Double_t fNoBUU0
 
Double_t fNoBUU1
 
Double_t fNoBUU2
 
Double_t fNoBUB
 
Double_t fNoUBU
 
Double_t fNoTotal
 
Double_t fNweirdPhotons
 
Double_t flambdah
 
Double_t fPixIndex
 
Double_t ftime
 
Double_t fPhiRot
 
TVector3 fPphoInit
 
TVector3 fStartVertex
 
TVector3 fPDSec
 
TVector3 fEVSec
 
TVector3 fEvSec
 
Double_t fPDPhi
 
Double_t fEVPhi
 
Double_t fPhiRotEV
 
Double_t Ang_pipe
 
Double_t Rout1
 
Double_t Rin1
 
Double_t Rin2
 
Double_t PlanB [9]
 
Double_t PlanU [9]
 
Double_t PlanR [6]
 
Double_t determint1
 
Double_t determint2
 
Double_t determint3
 
Double_t determint4
 
TString ReflectionType
 
TString ReflName
 
TArrayD fmatrixdata
 
Double_t fZin
 
Double_t fLowZ
 
Double_t fkxBar
 
Double_t fkyBar
 
Double_t fkzBar
 
TVector3 fPphoB
 
TVector3 fBBver1
 
TVector3 fBBver2
 
TVector3 fBBver3
 
TVector3 fBBver4
 
PndDrcEVPointEVpt
 
PndDrcEVPointEVt
 
PndDrcPDPointPpt
 
PndDrcPDHitpdhit
 
PndMCTracktr
 
TClonesArray * fLut [5]
 
TFile * fFile
 
TTree * fTree
 

Detailed Description

PhoGun.h

Class for Analysising DRC Cherenkov Photons fired in the DIRC barrel directly

Definition at line 47 of file PndPhoGunShort.h.

Constructor & Destructor Documentation

PndPhoGunShort::PndPhoGunShort ( )

Default constructor

Definition at line 40 of file PndPhoGunShort.cxx.

References fGeo, and fGeoH.

41 :FairTask("PndPhoGunShort")
42 {
43  fGeo = new PndGeoDrc();
44  fGeoH=NULL;
45 
46 }
PndGeoDrc * fGeo
Basic geometry data of barrel DRC.
PndGeoHandling * fGeoH
PndPhoGunShort::PndPhoGunShort ( Int_t  verbose)

Constructor with verbosity

Definition at line 49 of file PndPhoGunShort.cxx.

References fGeo, fGeoH, fVerbose, and verbose.

50  :FairTask("PndPhoGunShort")
51 {
52  fVerbose = verbose;
53  fGeo = new PndGeoDrc();
54  fGeoH=NULL;
55 }
#define verbose
PndGeoDrc * fGeo
Basic geometry data of barrel DRC.
PndGeoHandling * fGeoH
PndPhoGunShort::~PndPhoGunShort ( )
virtual

Destructor

Definition at line 57 of file PndPhoGunShort.cxx.

References fGeo, fGeoH, and fHistoList.

58 {
59  if (fGeo) delete fGeo;
60  if (fGeoH) delete fGeoH;
61  fHistoList->Delete();
62  delete fHistoList;
63 }
PndGeoDrc * fGeo
Basic geometry data of barrel DRC.
TList * fHistoList
PndGeoHandling * fGeoH

Member Function Documentation

void PndPhoGunShort::Exec ( Option_t *  option)
virtual

Executed task

Definition at line 185 of file PndPhoGunShort.cxx.

References fDetectorID, fGeoH, fVerbose, nevents, ProcessPhotonHit(), and PndGeoHandling::SetVerbose().

186 {
187  nevents++;
188 
190 
191  fDetectorID = 0;
192  /*if(nevents%1000==0)*/ cout<<"Event # "<<nevents<<endl;
194 }
PndGeoHandling * fGeoH
void SetVerbose(Int_t v)
Double_t PndPhoGunShort::FindReflectionType ( Double_t  xev,
Double_t  yev,
Double_t  zev,
TString  ReflType 
)
private

Definition at line 343 of file PndPhoGunShort.cxx.

References Ang_pipe, cos(), determint1, determint2, fDphi, fEVlen, fmatrixdata, fPhiRotEV, fpi, fZin, PlanB, PlanU, ReflectionType, Rin1, Rin2, Rout1, and sin().

Referenced by ProcessPhotonHit().

343  {
344 
345  //cout<<"fPhiRotEV = "<<fPhiRotEV<< ", fDphi = "<<fDphi<<", Ang_pipe = "<<Ang_pipe<<endl;
346 
347  PlanB[0]=(Rin1/cos((Ang_pipe/16.)*fpi/180.))*cos((fPhiRotEV-fDphi/2.)*fpi/180.);
348  PlanB[1]=(Rin1/cos((Ang_pipe/16.)*fpi/180.))*sin((fPhiRotEV-fDphi/2.)*fpi/180.);
349  PlanB[2]=fZin;//-119.6015;//fGeo->barBoxZUp()-0.30075;
350  PlanB[3]=(Rin1/cos((Ang_pipe/16.)*fpi/180.))*cos((fPhiRotEV-fDphi/2.)*fpi/180.);
351  PlanB[4]=(Rin1/cos((Ang_pipe/16.)*fpi/180.))*sin((fPhiRotEV-fDphi/2.)*fpi/180.);
352  PlanB[5]=fZin-fEVlen;//-149.6015;//fGeo->barBoxZUp()-fGeo->EVlen()-0.30075;
353  PlanB[6]=(Rin1/cos((Ang_pipe/16.)*fpi/180.))*cos((fPhiRotEV+fDphi/2.)*fpi/180.);
354  PlanB[7]=(Rin1/cos((Ang_pipe/16.)*fpi/180.))*sin((fPhiRotEV+fDphi/2.)*fpi/180.);
355  PlanB[8]=fZin-fEVlen;//-149.6015;//fGeo->barBoxZUp()-fGeo->EVlen()-0.30075;
356 
357  PlanU[0]=(Rin2/cos((Ang_pipe/16.)*fpi/180.))*cos((fPhiRotEV-fDphi/2.)*fpi/180.);
358  PlanU[1]=(Rin2/cos((Ang_pipe/16.)*fpi/180.))*sin((fPhiRotEV-fDphi/2.)*fpi/180.);
359  PlanU[2]=fZin;//-119.6015;//fGeo->barBoxZUp()-0.30075;
360  PlanU[3]=(Rout1/cos((Ang_pipe/16.)*fpi/180.))*cos((fPhiRotEV-fDphi/2.)*fpi/180.);
361  PlanU[4]=(Rout1/cos((Ang_pipe/16.)*fpi/180.))*sin((fPhiRotEV-fDphi/2.)*fpi/180.);
362  PlanU[5]=fZin-fEVlen;//-149.6015;//fGeo->barBoxZUp()-fGeo->EVlen()-0.30075;
363  PlanU[6]=(Rout1/cos((Ang_pipe/16.)*fpi/180.))*cos((fPhiRotEV+fDphi/2.)*fpi/180.);
364  PlanU[7]=(Rout1/cos((Ang_pipe/16.)*fpi/180.))*sin((fPhiRotEV+fDphi/2.)*fpi/180.);
365  PlanU[8]=fZin-fEVlen;//-149.6015;//fGeo->barBoxZUp()-fGeo->EVlen()-0.30075;
366 
367 //cout<<Rin1<<" "<<Rin2<<" "<<Rout1<<" "<<Rout1/cos((Ang_pipe/16.)*fpi/180.)<<endl;
368  ReflectionType="";
369  TMatrixD matrix1;
370  fmatrixdata[0]=xev-PlanB[0];
371  fmatrixdata[1]=yev-PlanB[1];
372  fmatrixdata[2]=zev-PlanB[2];
373  fmatrixdata[3]=PlanB[3]-PlanB[0];
374  fmatrixdata[4]=PlanB[4]-PlanB[1];
375  fmatrixdata[5]=PlanB[5]-PlanB[2];
376  fmatrixdata[6]=PlanB[6]-PlanB[0];
377  fmatrixdata[7]=PlanB[7]-PlanB[1];
378  fmatrixdata[8]=PlanB[8]-PlanB[2];
379  matrix1.Use(3,3,fmatrixdata.GetArray());
380  determint1 = matrix1.Determinant();
381 
382  if(determint1 > -0.02 && determint1 < 0.02)ReflectionType="B";
383  fmatrixdata.Reset(0);
384  fmatrixdata[0]=xev-PlanU[0];
385  fmatrixdata[1]=yev-PlanU[1];
386  fmatrixdata[2]=zev-PlanU[2];
387  fmatrixdata[3]=PlanU[3]-PlanU[0];
388  fmatrixdata[4]=PlanU[4]-PlanU[1];
389  fmatrixdata[5]=PlanU[5]-PlanU[2];
390  fmatrixdata[6]=PlanU[6]-PlanU[0];
391  fmatrixdata[7]=PlanU[7]-PlanU[1];
392  fmatrixdata[8]=PlanU[8]-PlanU[2];
393  matrix1.Use(3,3,fmatrixdata.GetArray());
394  determint2 = matrix1.Determinant();
395  if(determint2 > -0.02 && determint2 < 0.02)ReflectionType="U";
396 
397 
398  //cout<<"determinants: "<<determint1<<", "<<determint2<<endl;
399  //cout<<"rin2 = "<<Rin1<<", rin2 = "<<Rin2<<", rout = "<<Rout1<<endl;
400 }
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
Double_t PlanB[9]
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
Double_t determint2
TString ReflectionType
Double_t determint1
Double_t PlanU[9]
Double_t fPhiRotEV
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:52
void PndPhoGunShort::Finish ( )
virtual

Finish task

Definition at line 456 of file PndPhoGunShort.cxx.

References fFile, fLut, and fTree.

457 {
458  fTree->Fill();
459  fTree->Write();
460  fFile->Write();
461 
462  for(Int_t iLut=0; iLut<5; iLut++){
463  fLut[iLut]->Clear();
464  }
465  cout << "-I- PndDrcLutFill: Finish" << endl;
466 }
TClonesArray * fLut[5]
Double_t PndPhoGunShort::InBarCoordSyst ( TVector3  start,
TVector3 *  v1,
TVector3 *  v2,
TVector3 *  v3,
TVector3 *  v4 
)
private

Definition at line 403 of file PndPhoGunShort.cxx.

References Double_t, fDphi, fHThick, fLength, fpi, fPipehAngle, and fR.

Referenced by ProcessPhotonHit().

403  {
404 
405  // this function is used with fStartVertex to find the bar from which the photon originated
406 
407  Double_t startPhi = start.Phi()/fpi*180.; // [degrees]
408  //cout<<"-I- InBarCoordinateSystem: start phi = "<<startPhi<<endl;
409  //cout<<"-I- InBarCoordinateSystem: dphi = "<<Dphi<<endl;
410  Double_t PhiRot = 0.; //[degrees]
411  if(startPhi < 0.){startPhi = 360. + startPhi;}
412  if(startPhi > 0. && startPhi < 90.){
413  PhiRot = TMath::Floor(startPhi/fDphi) *fDphi + fDphi/2.;
414  }
415  if(startPhi > 90. && startPhi < 270.){
416  PhiRot = 90. + fPipehAngle + TMath::Floor((startPhi-90.-fPipehAngle)/fDphi) *fDphi + fDphi/2.;
417  }
418  if(startPhi > 270. && startPhi < 360.){
419  PhiRot = 270. + fPipehAngle + TMath::Floor((startPhi-270.-fPipehAngle)/fDphi) *fDphi + fDphi/2.;
420  }
421  //cout<<"-I- InBarCoordinateSystem: PhiRot = "<<PhiRot<<endl;
422 
423  // create initial barbox:
424 
425  TVector3 ver1, ver2, ver3, ver4;
426  ver1.SetXYZ(fR-fHThick, fLength/2. ,0.);
427  ver2.SetXYZ(fR+fHThick, fLength/2. ,0.);
428  ver3.SetXYZ(fR+fHThick, -fLength/2. ,0.);
429  ver4.SetXYZ(fR-fHThick, -fLength/2. ,0.);
430 
431  ver1.RotateZ(PhiRot/180.*fpi);
432  ver2.RotateZ(PhiRot/180.*fpi);
433  ver3.RotateZ(PhiRot/180.*fpi);
434  ver4.RotateZ(PhiRot/180.*fpi);
435 
436  *v1 = ver1;
437  *v2 = ver2;
438  *v3 = ver3;
439  *v4 = ver4;
440 
441  return PhiRot/180.*fpi;
442 }
Double_t fPipehAngle
Double_t
TVector3 v1
Definition: bump_analys.C:40
TVector3 v2
Definition: bump_analys.C:40
InitStatus PndPhoGunShort::Init ( )
virtual

Definition at line 66 of file PndPhoGunShort.cxx.

References Ang_pipe, PndGeoDrc::barHalfThick(), PndGeoDrc::BBoxAngle(), PndGeoDrc::BBoxGap(), PndGeoDrc::BBoxNum(), PndGeoDrc::boxGap(), PndGeoDrc::boxThick(), cos(), PndGeoDrc::EVdrop(), PndGeoDrc::EVoffset(), fBarBoxGap, fBarPointArray, fBboxNum, fDigiArray, fDphi, fEVdrop, fEVdz, fEVlen, fEVPointArray, fFile, fGeo, fGeoH, fHThick, fLength, fLut, fMCArray, fNmcp, fNoB, fNoBU0, fNoBU1, fNoBU2, fNoBUB, fNoBUU0, fNoBUU1, fNoBUU2, fNoDD, fNoTotal, fNoU0, fNoU1, fNoU2, fNoU3, fNoUB, fNoUBU, fNoUU0, fNoUU1, fNoUU2, fNoUU3, fNoUUU0, fNoUUU1, fNoUUU2, fNoUUU3, fNoUUU4, fNweirdPhotons, fOutputName, fPDHitArray, fPDPointArray, fpi, fPipehAngle, fPixelSize, fR, fRBottom, fShiftPix, fTree, InitLut(), PndGeoHandling::Instance(), PndGeoDrc::McpActiveArea(), PndGeoDrc::McpGap(), PndGeoDrc::McpSize(), nevents, Pi, PndGeoDrc::PipehAngle(), PndGeoDrc::PixelSize(), PndGeoDrc::radius(), Rin1, Rin2, Rout1, and PndGeoHandling::SetParContainers().

67 {
68  cout << " ---------- INITIALIZATION ------------" << endl;
69  nevents = 0;
70  // Get RootManager
71  FairRootManager* ioman = FairRootManager::Instance();
72  if ( ! ioman ) {
73  cout << "-E- PndPhoGunShort::Init: "
74  << "RootManager not instantiated!" << endl;
75  return kFATAL;
76  }
77 
78  // Get input array
79  fMCArray = (TClonesArray*) ioman->GetObject("MCTrack");
80  if ( ! fMCArray ) {
81  cout << "-W- PndPhoGunShort::Init: "
82  << "No MCTrack array!" << endl;
83  return kERROR;
84  }
85  // Get Photon point array
86  fPDPointArray = (TClonesArray*) ioman->GetObject("DrcPDPoint");
87  if ( ! fPDPointArray ) {
88  cout << "-W- PndPhoGunShort::Init: "
89  << "No DrcPDPoint array!" << endl;
90  return kERROR;
91  }
92  // Get input array
93  fPDHitArray = (TClonesArray*) ioman->GetObject("DrcPDHit");
94  if ( ! fPDHitArray ) {
95  cout << "-W- PndPhoGunShort::Init: "
96  << "No DrcPDHit array!" << endl;
97  return kERROR;
98  }
99 
100  fEVPointArray = (TClonesArray*) ioman->GetObject("DrcEVPoint");
101  if ( ! fEVPointArray ) {
102  cout << "-W- PndPhoGunShort::Init: "
103  << "No DrcEVPoint array!" << endl;
104  return kERROR;
105  }
106  fBarPointArray = (TClonesArray*) ioman->GetObject("DrcBarPoint");
107  if ( ! fBarPointArray ) {
108  cout << "-W- PndDrcLogLikeli::Init: "
109  << "No DrcBarPoint array!" << endl;
110  return kERROR;
111  }
112  fDigiArray = (TClonesArray*) ioman->GetObject("DrcDigi");
113  if ( ! fDigiArray ) {
114  cout << "-W- PndPhoGunShortP::Init: " << "No DrcDigi array!" << endl;
115  return kERROR;
116  }
117 
118  // Get parameters:
119  fpi = TMath::Pi();
120  fR = fGeo->radius();
122  fBboxNum = fGeo->BBoxNum();
124  fBarBoxGap = fGeo->BBoxGap();
125  fLength = (180. - 2.*fPipehAngle - fBarBoxGap/fR*(fBboxNum/2. - 1.)/fpi*180.)/(fBboxNum/2.) * fR/ 180.*fpi;
126  fDphi = fGeo->BBoxAngle();//21.65 degrees = 2.*(180. - 2*fPipehAngle)/ fBboxNum; //[degrees]
128 
129  //Double_t EVlen = gGeoManager->GetVolume("DrcEVSensor")->GetShape()->Dz();
130  //fEVlen = fGeo->EVlen();
131 
132  fEVdrop = fGeo->EVdrop();
133  fEVlen =fEVdz;
134  fRBottom = fR - fHThick - fEVdrop + fR*(1.-cos(fDphi/2./180.*fpi))/cos(fDphi/2./180.*fpi);
135  //cout<<"fR = "<<fR<<", fHThick = "<<fHThick<<", fEVdrop = "<<fEVdrop<<", fRBottom = "<<fRBottom<<", dphi = "<<fDphi<<endl;
136 // cout<<"EV depth = "<<fEVlen<<", n mcp = "<<fNmcp<<endl;
137 
138  Rin1 = fR - fHThick - fEVdrop - fGeo->boxGap() - fGeo->boxThick();
139  Rin2 = fR + fHThick + fGeo->EVoffset()+ fGeo->boxGap() + fGeo->boxThick() ;
141 
142  if(fShiftPix){
143  Rin1 = fR - fHThick - fEVdrop - fGeo->boxGap() - fGeo->boxThick()- fPixelSize/2.;
145  cout<<"shifted pixels"<<endl;
146  }
147 
148  Ang_pipe = 180. - 2.*fPipehAngle;
149  fNoDD=0;
150  fNoU0=0;fNoU1=0;fNoU2=0;fNoU3=0;
151  fNoB=0;
152  fNoBU0=0;fNoBU1=0;fNoBU2=0;
153  fNoUB=0;
154  fNoUU0=0;fNoUU1=0;fNoUU2=0;fNoUU3=0;
156  fNoBUU0=0;fNoBUU1=0;fNoBUU2=0;
157  fNoUBU=0;
158  fNoBUB=0;
159  fNoTotal=0;
160  fNweirdPhotons = 0;
161 
162  fFile = TFile::Open(fOutputName,"RECREATE");
163  fTree = new TTree("dircsim","Look-up table for DIRC");
164  for(Int_t iLut=0; iLut<5; iLut++){
165  fLut[iLut] = new TClonesArray("PndDrcLutNode");
166  fTree->Branch(Form("LUT%d",iLut),&fLut[iLut],256000,0);
167  }
168 
169 
170  InitLut();
171 
172  if ( fGeoH == NULL ) fGeoH = PndGeoHandling::Instance();
174 
175  //fGeoH->GetGeoManager();
176  //const Double_t *truuu = gGeoManager->GetVolume("DrcEVSensor")->GetShape()->GetTransform()->GetTranslation();
177  //fLowZ = truuu[2]-fGeo->EVlen();
178  //cout<<"posZ = "<<fLowZ<<endl;
179 
180  cout << "-I- PndPhoGunShort: Intialization successfull" << endl;
181  return kSUCCESS;
182 }
Double_t fPipehAngle
TClonesArray * fBarPointArray
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
TClonesArray * fMCArray
Double_t fBarBoxGap
Double_t BBoxNum()
Definition: PndGeoDrc.h:136
Double_t fNweirdPhotons
virtual void SetParContainers()
TClonesArray * fPDPointArray
TClonesArray * fLut[5]
Double_t BBoxGap()
Definition: PndGeoDrc.h:130
Double_t EVdrop()
Definition: PndGeoDrc.h:142
Double_t boxGap()
Definition: PndGeoDrc.h:116
PndGeoDrc * fGeo
Basic geometry data of barrel DRC.
Double_t EVoffset()
Definition: PndGeoDrc.h:145
Double_t PipehAngle()
Definition: PndGeoDrc.h:139
TClonesArray * fPDHitArray
PndGeoHandling * fGeoH
Double_t McpActiveArea()
Definition: PndGeoDrc.h:166
Double_t McpSize()
Definition: PndGeoDrc.h:163
Double_t BBoxAngle()
Definition: PndGeoDrc.h:133
TClonesArray * fDigiArray
Double_t boxThick()
Definition: PndGeoDrc.h:120
static PndGeoHandling * Instance()
TClonesArray * fEVPointArray
Double_t barHalfThick()
Definition: PndGeoDrc.h:96
Double_t Pi
Double_t McpGap()
Definition: PndGeoDrc.h:169
Double_t PixelSize()
Definition: PndGeoDrc.h:175
Double_t radius()
Definition: PndGeoDrc.h:92
Double_t fPixelSize
void PndPhoGunShort::InitLut ( )
private

Definition at line 444 of file PndPhoGunShort.cxx.

References fLut, and n.

Referenced by Init().

445 {
446  Int_t Nnodes = 30000;
447  for(Int_t iLut=0; iLut<5; iLut++){
448  TClonesArray &fLuta = *fLut[iLut];
449  for (Long64_t n=0; n<Nnodes; n++) {
450  new((fLuta)[n]) PndDrcLutNode(-1);
451  }
452  }
453 }
int n
TClonesArray * fLut[5]
void PndPhoGunShort::ProcessPhotonHit ( )
private

Definition at line 197 of file PndPhoGunShort.cxx.

References ambiguity, At, Bool_t, EVpt, EVt, fabs(), fBarId, fBarPointArray, fBBver1, fBBver2, fBBver3, fBBver4, fDigi, fDigiArray, fDphi, fEVlen, fEVPhi, fEVPointArray, fEVSec, FindReflectionType(), fkxBar, fkyBar, fkzBar, fLut, fmatrixdata, fMCArray, fNoB, fNoBU0, fNoBU1, fNoBU2, fNoBUB, fNoBUU0, fNoBUU1, fNoBUU2, fNoDD, fNoTotal, fNoU0, fNoU1, fNoU2, fNoU3, fNoUB, fNoUBU, fNoUU0, fNoUU1, fNoUU2, fNoUU3, fNoUUU0, fNoUUU1, fNoUUU2, fNoUUU3, fNoUUU4, fNweirdPhotons, fPDHitArray, fPDPointArray, fPhiRot, fPhiRotEV, fpi, fPixIndex, fPphoB, fPphoInit, fStartVertex, ftime, fZin, PndDrcBarPoint::GetBarId(), PndDrcPDPoint::GetBarPointID(), PndDrcDigi::GetDetectorId(), PndDrcDigi::GetIndex(), PndMCTrack::GetMomentum(), PndDrcPDHit::GetPosition(), PndDrcPDHit::GetRefIndex(), PndDrcDigi::GetSensorId(), PndMCTrack::GetStartVertex(), InBarCoordSyst(), CAMath::Nint(), pdhit, Ppt, ReflectionType, ReflName, and tr.

Referenced by Exec().

198 {
199  Int_t EVEntry;
200  Int_t SelectionName=0;
201  Int_t PhiSec=0;
202  fmatrixdata.Set(9);
203  fmatrixdata.Reset(0);
204 
205  // Loop over PndDrcPDHits
206  for(Int_t k=0; k<fPDHitArray->GetEntriesFast(); k++) {
207 
208  pdhit = (PndDrcPDHit*)fPDHitArray->At(k);
209 
210  Int_t mcPDRef= pdhit->GetRefIndex();
211  fDigi = (PndDrcDigi*) fDigiArray->At(mcPDRef);
212  Int_t pointID= fDigi->GetIndex(0);
213  Ppt = (PndDrcPDPoint*)fPDPointArray->At(pointID);
214 
215  Int_t trID= Ppt->GetTrackID();
216  tr = (PndMCTrack*)fMCArray->At(trID);
217 
218  EVpt = (PndDrcEVPoint*)fEVPointArray->At(mcPDRef);
219 
220  if(trID!=EVpt->GetTrackID())cout<<"different track IDs"<<endl;
221  EVEntry = fEVPointArray->GetEntriesFast();
222  if(trID!=EVpt->GetTrackID())continue;
224 // fBarId = fBarPoint->GetDetectorID();
225  fBarId = fBarPoint->GetBarId();
226 
227  // production point of the photon
229 
230  // find PHIrot to get to the bar coord. syst:
232 
233  // check if the last of the EV points is on the PD plane as the reflection was in the grease layer/window/photocathode, then exclude this point
234  if(((PndDrcEVPoint*)fEVPointArray->At(EVEntry-1))->GetZ() < ((PndDrcEVPoint*)fEVPointArray->At(0))->GetZ()-fEVlen+0.01){EVEntry = EVEntry - 1;}
235  //check if there are EVpoints on the PD plane and these points are not the last in the array, then reject the photon
236  Bool_t wePho = kFALSE;
237  for(Int_t k2 = 1; k2 < fEVPointArray->GetEntriesFast()-1; k2++){
238  if((((PndDrcEVPoint*)fEVPointArray->At(k2))->GetZ()) < ((PndDrcEVPoint*)fEVPointArray->At(0))->GetZ()-fEVlen+0.01){
239  //cout<<"this is a weird photon - has EV points at the PD plane!!!"<<endl;
240  fNweirdPhotons += 1;
241  wePho = kTRUE;
242  }
243  }
244  if(wePho == kTRUE)continue;
245 
246  if(EVEntry>4){
247  fNweirdPhotons += 1;
248  continue;
249  }
250  ReflName="";
251  if(EVEntry==1)ReflName="DD";
252  fZin = ((PndDrcEVPoint*)fEVPointArray->At(0))->GetZ();
253  if(EVEntry>1){
254  for(Int_t k1 = 1; k1 < EVEntry; k1++){
255  Int_t EVS=0;
256  EVt = (PndDrcEVPoint*)fEVPointArray->At(k1);
257  fEVSec.SetXYZ(EVt->GetX(),EVt->GetY(),EVt->GetZ());
259  fPhiRotEV=fPhiRotEV*180./fpi;
260  fEVPhi=fEVSec.Phi()*180./fpi;
261  //cout<<"EV phi = "<<fEVPhi<<", fPhiRot = "<<fPhiRot<<", fDphi = "<<fDphi<<endl;
262  PhiSec = TMath::Nint(fabs(fEVPhi-fPhiRot*180./fpi)/fDphi);
263  if(fEVPhi<0.){
264  fEVPhi=360.+fEVPhi;
265  PhiSec = TMath::Nint(fabs(fEVPhi-fPhiRot*180./fpi-360.)/fDphi);
266  }
267  //cout<<"phi sec = "<<PhiSec<<endl;
268  FindReflectionType (EVt->GetX(),EVt->GetY(),EVt->GetZ(),ReflectionType);
269  ReflName.Append(ReflectionType);
270  }
271  if(EVEntry==2)ReflName.Append("1");
272  if(EVEntry==3)ReflName.Append("2");
273  if(EVEntry==4)ReflName.Append("3");
274  if(EVEntry==5)ReflName.Append("4");
275  if(EVEntry>=6)ReflName.Append("XX");
276  ReflName.Prepend(Form("%d",PhiSec));
277  }
278  // 24 ambiguities in total:
279  if(ReflName == "DD"){fNoDD +=1;SelectionName=1;ambiguity =1;}
280  if(ReflName == "0B1"){fNoB +=1;SelectionName=1;ambiguity =2;}
281  if(ReflName == "0U1"){fNoU0 +=1;SelectionName=1;ambiguity =3;}
282  if(ReflName == "1U1"){fNoU1 +=1;SelectionName=1;ambiguity =4;}
283  if(ReflName == "2U1"){fNoU2 +=1;SelectionName=1;ambiguity =5;}
284  if(ReflName == "3U1"){fNoU3 +=1;SelectionName=1;ambiguity =6;}
285  if(ReflName == "0BU2"){fNoBU0 +=1;SelectionName=1;ambiguity =7;}
286  if(ReflName == "1BU2"){fNoBU1 +=1;SelectionName=1;ambiguity =8;}
287  if(ReflName == "2BU2"){fNoBU2 +=1;SelectionName=1;ambiguity =9;}
288  if(ReflName == "0UB2"){fNoUB +=1; SelectionName=1;ambiguity =10;}
289  if(ReflName == "0UU2"){fNoUU0 +=1; SelectionName=1;ambiguity =11;}
290  if(ReflName == "1UU2"){fNoUU1 +=1;SelectionName=1;ambiguity =12;}
291  if(ReflName == "2UU2"){fNoUU2 +=1;SelectionName=1;ambiguity =13;}
292  if(ReflName == "3UU2"){fNoUU3 +=1; SelectionName=1;ambiguity =14;}
293  if(ReflName == "0UUU3"){fNoUUU0 +=1; SelectionName=1;ambiguity =15;}
294  if(ReflName == "1UUU3"){fNoUUU1 +=1; SelectionName=1;ambiguity =16;}
295  if(ReflName == "2UUU3"){fNoUUU2 +=1; SelectionName=1;ambiguity =17;}
296  if(ReflName == "3UUU3"){fNoUUU3 +=1; SelectionName=1;ambiguity =18;}
297  if(ReflName == "4UUU3"){fNoUUU4 +=1; SelectionName=1;ambiguity =19;}
298  if(ReflName == "0BUU3"){fNoBUU0 +=1; SelectionName=1;ambiguity =20;}
299  if(ReflName == "1BUU3"){fNoBUU1 +=1; SelectionName=1;ambiguity =21;}
300  if(ReflName == "2BUU3"){fNoBUU2 +=1; SelectionName=1;ambiguity =22;}
301  if(ReflName == "BUB3"){fNoBUB +=1; SelectionName=1;ambiguity =23;}
302  if(ReflName == "UBU3"){fNoUBU +=1; SelectionName=1;ambiguity =24;}
303  if(ReflName == ""){ambiguity =-99;}
304 
305  if(SelectionName==0.){continue;}
306 
308 // fPixIndex = pdhit->GetDetectorID();
309  ftime = Ppt->GetTime();//pdhit->GetTime();
310 
311  // Photon initial momentum in the bar coord. syst.
312  fPphoInit = tr->GetMomentum();
313  //fPphoB = (fGeoH->MasterToLocalShortId(fPphoInit, fBarId) - fGeoH->MasterToLocalShortId((0.,0.,0.),fBarId)).Unit();
314  // fPphoB.Print();
315  // Photon initial momentum in the bar coord. syst.
316  fPphoB = (tr->GetMomentum()).Unit();
317  fPphoB.RotateZ(-fPhiRot);
318 // fkxBar = -fPphoB.Y();
319 // fkyBar = fPphoB.X();
320 // fkzBar = fPphoB.Z();
321  fkxBar = fPphoB.X();
322  fkyBar = fPphoB.Y();
323  fkzBar = fPphoB.Z();
324  fPphoB.SetXYZ(fkxBar,fkyBar,fkzBar);
325 
327 // ((PndDrcLutNodeH*)(fLut->At(fPixIndex)))->AddEntry(fPphoB,ambiguity,ftime);
328 // ((PndDrcLutNodeH*)(fLut->At(fPixIndex)))->SetPos(pdhit->GetPosition());
329 
330 
331 // Double_t rrr = sqrt(pow(pdhit->GetX(),2) + pow(pdhit->GetY(),2));
332 
333  //if(ambiguity>2){
334 // cout<<"pix "<<fPixIndex<<", reflection type = "<<ReflName<<", ambiguity = "<<ambiguity<<", kx = "<<fkxBar<<", R = "<<rrr<<endl;
335  //}
336 
337  fNoTotal +=1;
338 
339  }// photon hits
340 }
PndMCTrack * tr
TClonesArray * fBarPointArray
TClonesArray * fMCArray
Double_t fNweirdPhotons
Int_t GetBarId() const
TVector3 fPphoInit
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
TClonesArray * fPDPointArray
TClonesArray * fLut[5]
TString ReflectionType
TClonesArray * fPDHitArray
Double_t fPixIndex
Int_t GetBarPointID() const
Definition: PndDrcPDPoint.h:56
virtual Int_t GetRefIndex()
Definition: PndDrcPDHit.h:57
PndDrcPDPoint * Ppt
Double_t InBarCoordSyst(TVector3, TVector3 *, TVector3 *, TVector3 *, TVector3 *)
PndDrcEVPoint * EVt
TClonesArray * fDigiArray
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
TVector3 GetPosition() const
Definition: PndDrcPDHit.h:58
PndDrcEVPoint * EVpt
TClonesArray * fEVPointArray
Double_t fPhiRotEV
PndDrcPDHit * pdhit
Int_t GetSensorId() const
Definition: PndDrcDigi.h:91
PndDrcDigi * fDigi
TVector3 fStartVertex
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
Int_t GetDetectorId() const
Definition: PndDrcDigi.h:90
cout<<"the Event No is "<< i<< endl;{{if(hit_array->GetEntriesFast()!=mc_array->GetEntriesFast()) continue;PndSdsHit *hit=(PndSdsHit *) hit_array-> At(j)
Definition: anaLmdCluster.C:71
int Nint(float x)
Definition: PndCAMath.h:117
Double_t FindReflectionType(Double_t, Double_t, Double_t, TString)
Int_t GetIndex(int i=0) const
Definition: PndDrcDigi.h:94
void PndPhoGunShort::SetDefaultParameters ( )
private

Set the parameters to the default values.

void PndPhoGunShort::SetEVDepth ( Float_t  EVdepth = 30.)
inline

Definition at line 70 of file PndPhoGunShort.h.

References fEVdz.

70 {fEVdz = EVdepth;}
void PndPhoGunShort::SetNmcp ( Float_t  nmcp = 5.)
inline

Definition at line 71 of file PndPhoGunShort.h.

References fNmcp.

71 {fNmcp = nmcp;}
void PndPhoGunShort::SetOutputFile ( TString  outName = "lutnode.root")
inline

Definition at line 68 of file PndPhoGunShort.h.

References fOutputName.

68 {fOutputName = outName;}
void PndPhoGunShort::SetShiftedPix ( Bool_t  spix = kFALSE)
inline

Definition at line 72 of file PndPhoGunShort.h.

References fShiftPix.

72 {fShiftPix = spix;}

Member Data Documentation

Int_t PndPhoGunShort::ambiguity
private

Definition at line 108 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

Double_t PndPhoGunShort::Ang_pipe
private

Definition at line 153 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), and Init().

Double_t PndPhoGunShort::determint1
private

Definition at line 154 of file PndPhoGunShort.h.

Referenced by FindReflectionType().

Double_t PndPhoGunShort::determint2
private

Definition at line 154 of file PndPhoGunShort.h.

Referenced by FindReflectionType().

Double_t PndPhoGunShort::determint3
private

Definition at line 154 of file PndPhoGunShort.h.

Double_t PndPhoGunShort::determint4
private

Definition at line 154 of file PndPhoGunShort.h.

PndDrcEVPoint* PndPhoGunShort::EVpt
private

Definition at line 169 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

PndDrcEVPoint* PndPhoGunShort::EVt
private

Definition at line 170 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

Double_t PndPhoGunShort::fBarBoxGap
private

Definition at line 121 of file PndPhoGunShort.h.

Referenced by Init().

Int_t PndPhoGunShort::fBarId
private

Definition at line 113 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TClonesArray* PndPhoGunShort::fBarPointArray
private

Definition at line 82 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fBboxNum
private

Definition at line 119 of file PndPhoGunShort.h.

Referenced by Init().

TVector3 PndPhoGunShort::fBBver1
private

Definition at line 164 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TVector3 PndPhoGunShort::fBBver2
private

Definition at line 165 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TVector3 PndPhoGunShort::fBBver3
private

Definition at line 166 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TVector3 PndPhoGunShort::fBBver4
private

Definition at line 167 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

Int_t PndPhoGunShort::fDetectorID
private

Definition at line 80 of file PndPhoGunShort.h.

Referenced by Exec().

PndDrcDigi* PndPhoGunShort::fDigi
private

Definition at line 94 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TClonesArray* PndPhoGunShort::fDigiArray
private

Definition at line 87 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fDphi
private

Definition at line 123 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), InBarCoordSyst(), Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fEVdrop
private

Definition at line 127 of file PndPhoGunShort.h.

Referenced by Init().

Double_t PndPhoGunShort::fEVdz
private

Definition at line 125 of file PndPhoGunShort.h.

Referenced by Init(), and SetEVDepth().

Double_t PndPhoGunShort::fEVlen
private

Definition at line 124 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fEVPhi
private

Definition at line 152 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TClonesArray* PndPhoGunShort::fEVPointArray
private

Definition at line 86 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

TVector3 PndPhoGunShort::fEVSec
private

Definition at line 151 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TVector3 PndPhoGunShort::fEvSec
private

Definition at line 151 of file PndPhoGunShort.h.

TFile* PndPhoGunShort::fFile
private

Definition at line 176 of file PndPhoGunShort.h.

Referenced by Finish(), and Init().

PndGeoDrc* PndPhoGunShort::fGeo
private

Basic geometry data of barrel DRC.

Definition at line 95 of file PndPhoGunShort.h.

Referenced by Init(), PndPhoGunShort(), and ~PndPhoGunShort().

PndGeoHandling* PndPhoGunShort::fGeoH
private

Definition at line 96 of file PndPhoGunShort.h.

Referenced by Exec(), Init(), PndPhoGunShort(), and ~PndPhoGunShort().

TList* PndPhoGunShort::fHistoList
private

Definition at line 91 of file PndPhoGunShort.h.

Referenced by ~PndPhoGunShort().

Double_t PndPhoGunShort::fHThick
private

Definition at line 118 of file PndPhoGunShort.h.

Referenced by InBarCoordSyst(), and Init().

Double_t PndPhoGunShort::fkxBar
private

Definition at line 160 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

Double_t PndPhoGunShort::fkyBar
private

Definition at line 160 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

Double_t PndPhoGunShort::fkzBar
private

Definition at line 160 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

Double_t PndPhoGunShort::flambdah
private

Definition at line 145 of file PndPhoGunShort.h.

Double_t PndPhoGunShort::fLength
private

Definition at line 122 of file PndPhoGunShort.h.

Referenced by InBarCoordSyst(), and Init().

Double_t PndPhoGunShort::fLowZ
private

Definition at line 158 of file PndPhoGunShort.h.

TClonesArray* PndPhoGunShort::fLut[5]
private

Definition at line 175 of file PndPhoGunShort.h.

Referenced by Finish(), Init(), InitLut(), and ProcessPhotonHit().

TArrayD PndPhoGunShort::fmatrixdata
private

Definition at line 156 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), and ProcessPhotonHit().

TClonesArray* PndPhoGunShort::fMCArray
private

Definition at line 85 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Float_t PndPhoGunShort::fNmcp
private

Definition at line 129 of file PndPhoGunShort.h.

Referenced by Init(), and SetNmcp().

Double_t PndPhoGunShort::fNoB
private

Definition at line 134 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoBU0
private

Definition at line 136 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoBU1
private

Definition at line 136 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoBU2
private

Definition at line 136 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoBUB
private

Definition at line 140 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoBUU0
private

Definition at line 139 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoBUU1
private

Definition at line 139 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoBUU2
private

Definition at line 139 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoDD
private

Definition at line 132 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoTotal
private

Definition at line 142 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoU0
private

Definition at line 133 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoU1
private

Definition at line 133 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoU2
private

Definition at line 133 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoU3
private

Definition at line 133 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUB
private

Definition at line 135 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUBU
private

Definition at line 141 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUU0
private

Definition at line 137 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUU1
private

Definition at line 137 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUU2
private

Definition at line 137 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUU3
private

Definition at line 137 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUUU0
private

Definition at line 138 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUUU1
private

Definition at line 138 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUUU2
private

Definition at line 138 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUUU3
private

Definition at line 138 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNoUUU4
private

Definition at line 138 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fNweirdPhotons
private

Definition at line 143 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

TString PndPhoGunShort::fOutputName
private

Definition at line 110 of file PndPhoGunShort.h.

Referenced by Init(), and SetOutputFile().

TClonesArray* PndPhoGunShort::fPDHitArray
private

Definition at line 84 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fPDPhi
private

Definition at line 152 of file PndPhoGunShort.h.

TClonesArray* PndPhoGunShort::fPDPointArray
private

Definition at line 83 of file PndPhoGunShort.h.

Referenced by Init(), and ProcessPhotonHit().

TVector3 PndPhoGunShort::fPDSec
private

Definition at line 151 of file PndPhoGunShort.h.

Double_t PndPhoGunShort::fPhiRot
private

Definition at line 148 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

Double_t PndPhoGunShort::fPhiRotEV
private

Definition at line 152 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fpi
private

Definition at line 126 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), InBarCoordSyst(), Init(), and ProcessPhotonHit().

Double_t PndPhoGunShort::fPipehAngle
private

Definition at line 120 of file PndPhoGunShort.h.

Referenced by InBarCoordSyst(), and Init().

Double_t PndPhoGunShort::fPixelSize
private

Definition at line 128 of file PndPhoGunShort.h.

Referenced by Init().

Double_t PndPhoGunShort::fPixIndex
private

Definition at line 146 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TVector3 PndPhoGunShort::fPphoB
private

Definition at line 161 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TVector3 PndPhoGunShort::fPphoInit
private

Definition at line 149 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

Double_t PndPhoGunShort::fR
private

Definition at line 116 of file PndPhoGunShort.h.

Referenced by InBarCoordSyst(), and Init().

Double_t PndPhoGunShort::fRBottom
private

Definition at line 117 of file PndPhoGunShort.h.

Referenced by Init().

TFile* PndPhoGunShort::froot
private

Definition at line 92 of file PndPhoGunShort.h.

Bool_t PndPhoGunShort::fShiftPix
private

Definition at line 130 of file PndPhoGunShort.h.

Referenced by Init(), and SetShiftedPix().

TVector3 PndPhoGunShort::fStartVertex
private

Definition at line 150 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TString PndPhoGunShort::fTableName
private

Definition at line 111 of file PndPhoGunShort.h.

Double_t PndPhoGunShort::ftilt
private

Definition at line 112 of file PndPhoGunShort.h.

Double_t PndPhoGunShort::ftime
private

Definition at line 147 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TTree* PndPhoGunShort::fTree
private

Definition at line 177 of file PndPhoGunShort.h.

Referenced by Finish(), and Init().

Int_t PndPhoGunShort::fVerbose
private

Verbosity level

Definition at line 105 of file PndPhoGunShort.h.

Referenced by Exec(), and PndPhoGunShort().

Double_t PndPhoGunShort::fZin
private

Definition at line 157 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), and ProcessPhotonHit().

Int_t PndPhoGunShort::nevents
private

Definition at line 107 of file PndPhoGunShort.h.

Referenced by Exec(), and Init().

PndDrcPDHit* PndPhoGunShort::pdhit
private

Definition at line 172 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

Double_t PndPhoGunShort::PlanB[9]
private

Definition at line 153 of file PndPhoGunShort.h.

Referenced by FindReflectionType().

Double_t PndPhoGunShort::PlanR[6]
private

Definition at line 153 of file PndPhoGunShort.h.

Double_t PndPhoGunShort::PlanU[9]
private

Definition at line 153 of file PndPhoGunShort.h.

Referenced by FindReflectionType().

PndDrcPDPoint* PndPhoGunShort::Ppt
private

Definition at line 171 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

TString PndPhoGunShort::ReflectionType
private

Definition at line 155 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), and ProcessPhotonHit().

TString PndPhoGunShort::ReflName
private

Definition at line 155 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().

Double_t PndPhoGunShort::Rin1
private

Definition at line 153 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), and Init().

Double_t PndPhoGunShort::Rin2
private

Definition at line 153 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), and Init().

Double_t PndPhoGunShort::Rout1
private

Definition at line 153 of file PndPhoGunShort.h.

Referenced by FindReflectionType(), and Init().

PndMCTrack* PndPhoGunShort::tr
private

Definition at line 173 of file PndPhoGunShort.h.

Referenced by ProcessPhotonHit().


The documentation for this class was generated from the following files: