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

#include <PndDrcLutFill.h>

Inheritance diagram for PndDrcLutFill:

Public Member Functions

 PndDrcLutFill ()
 
 PndDrcLutFill (Int_t verbose)
 
 PndDrcLutFill (Int_t verbose, TString outfilename)
 
virtual ~PndDrcLutFill ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *option)
 
virtual void Finish ()
 
void InitLut ()
 
void SetOutputFile (TString outfilename="luttab.root")
 

Private Member Functions

void ProcessPhotonHit ()
 
void SetDefaultParameters ()
 

Private Attributes

PndGeoDrcfGeo
 
Int_t fDetectorID
 
Double_t fBboxNum
 
Double_t fPipehAngle
 
Double_t fDphi
 
TClonesArray * fMCArray
 
TClonesArray * fBarPointArray
 
TClonesArray * fPDPointArray
 
TClonesArray * fEVPointArray
 
TClonesArray * fDigiArray
 
TClonesArray * fPDHitArray
 
TClonesArray * fLut [5]
 
TFile * fFile
 
TTree * fTree
 
PndMCTrackfMCTrack
 
PndDrcBarPointfBarPoint
 
PndDrcPDPointfPDPoint
 
PndDrcDigifDigi
 
PndDrcPDHitfPDHit
 
PndDrcEVPointfEVPoint
 
Int_t fVerbose
 
Int_t nevents
 
TString fOutputFile
 

Detailed Description

Definition at line 26 of file PndDrcLutFill.h.

Constructor & Destructor Documentation

PndDrcLutFill::PndDrcLutFill ( )

Definition at line 22 of file PndDrcLutFill.cxx.

References fOutputFile.

23 :FairTask("PndDrcLutFill")
24 {
25  fOutputFile = "luttab.root";
26 }
TString fOutputFile
Definition: PndDrcLutFill.h:85
PndDrcLutFill::PndDrcLutFill ( Int_t  verbose)

Definition at line 29 of file PndDrcLutFill.cxx.

References fOutputFile, fVerbose, and verbose.

30  :FairTask("PndDrcLutFill",verbose)
31 {
32  fVerbose = verbose;
33  fOutputFile = "luttab.root";
34 }
#define verbose
TString fOutputFile
Definition: PndDrcLutFill.h:85
PndDrcLutFill::PndDrcLutFill ( Int_t  verbose,
TString  outfilename 
)

Definition at line 36 of file PndDrcLutFill.cxx.

References fOutputFile, fVerbose, and verbose.

37  :FairTask("PndDrcLutFill",verbose)
38 {
39  fVerbose = verbose;
40  fOutputFile = outfilename;
41 }
#define verbose
TString fOutputFile
Definition: PndDrcLutFill.h:85
PndDrcLutFill::~PndDrcLutFill ( )
virtual

Definition at line 44 of file PndDrcLutFill.cxx.

45 {
46 
47 }

Member Function Documentation

void PndDrcLutFill::Exec ( Option_t *  option)
virtual

Definition at line 133 of file PndDrcLutFill.cxx.

References fDetectorID, nevents, and ProcessPhotonHit().

134 {
135  nevents++;
136  // if(fVerbose>0 && nevents%1000==0)
137  std::cout<<"Event # "<< nevents<<std::endl;
138  fDetectorID = 0;
140 }
void ProcessPhotonHit()
void PndDrcLutFill::Finish ( )
virtual

Definition at line 235 of file PndDrcLutFill.cxx.

References fFile, fLut, and fTree.

236 {
237  fTree->Fill();
238  fTree->Write();
239  fFile->Write();
240 
241  for(Int_t l=0; l<5; l++){
242  fLut[l]->Clear();
243  }
244  cout << "-I- PndDrcLutFill: Finish" << endl;
245 }
TClonesArray * fLut[5]
Definition: PndDrcLutFill.h:67
InitStatus PndDrcLutFill::Init ( )
virtual

Definition at line 50 of file PndDrcLutFill.cxx.

References PndGeoDrc::BBoxNum(), Double_t, fBarPointArray, fBboxNum, fDigiArray, fDphi, fEVPointArray, fFile, fGeo, fLut, fMCArray, fOutputFile, fPDHitArray, fPDPointArray, fPipehAngle, fTree, InitLut(), nevents, and PndGeoDrc::PipehAngle().

51 {
52  cout << " ---------- INITIALIZATION ------------" << endl;
53  nevents = 0;
54  // Get RootManager
55  FairRootManager* ioman = FairRootManager::Instance();
56  if ( ! ioman ) {
57  cout << "-E- PndDrcLutFill::Init: " << "RootManager not instantiated!" << endl;
58  return kFATAL;
59  }
60 
61  // Get input array
62  fMCArray = (TClonesArray*) ioman->GetObject("MCTrack");
63  if ( ! fMCArray ) {
64  cout << "-W- PndDrcLutFill::Init: " << "No MCTrack array!" << endl;
65  return kERROR;
66  }
67 
68  // Get bar points array
69  fBarPointArray = (TClonesArray*) ioman->GetObject("DrcBarPoint");
70  if ( ! fBarPointArray ) {
71  cout << "-W- PndDrcLutReco::Init: " << "No DrcBarPoint array!" << endl;
72  return kERROR;
73  }
74 
75  // Get Photon point array
76  fPDPointArray = (TClonesArray*) ioman->GetObject("DrcPDPoint");
77  if ( ! fPDPointArray ) {
78  cout << "-W- PndDrcLutFill::Init: " << "No DrcPDPoint array!" << endl;
79  return kERROR;
80  }
81 
82  // Get ev points array
83  fEVPointArray = (TClonesArray*) ioman->GetObject("DrcEVPoint");
84  if ( ! fEVPointArray ) {
85  cout << "-W- PndDrcLutFill::Init: " << "No DrcEVPoint array!" << endl;
86  return kERROR;
87  }
88 
89  // Get digi array
90  fDigiArray = (TClonesArray*) ioman->GetObject("DrcDigi");
91  if ( ! fDigiArray ) {
92  cout << "-W- PndDrcLutFill::Init: " << "No DrcDigi array!" << endl;
93  return kERROR;
94  }
95  // Get input array
96  fPDHitArray = (TClonesArray*) ioman->GetObject("DrcPDHit");
97  if ( ! fPDHitArray ) {
98  cout << "-W- PndDrcLutFill::Init: " << "No DrcPDHit array!" << endl;
99  return kERROR;
100  }
101 
102  fFile = TFile::Open(fOutputFile,"RECREATE");
103  fTree = new TTree("dircsim","Look-up table for DIRC");
104  for(Int_t l=0; l<5; l++){
105  fLut[l] = new TClonesArray("PndDrcLutNode");
106  fTree->Branch(Form("LUT%d",l),&fLut[l],256000,0);
107  }
108 
109  InitLut();
110 
111  fGeo = new PndGeoDrc();
112  fBboxNum = fGeo->BBoxNum();
114  fDphi = 2.*(180. - 2*fPipehAngle)/(Double_t)fGeo->BBoxNum();
115 
116  cout << "-I- PndDrcLutFill: Intialization successfull" << endl;
117  return kSUCCESS;
118 
119 }
Double_t fPipehAngle
Definition: PndDrcLutFill.h:59
Double_t BBoxNum()
Definition: PndGeoDrc.h:136
PndGeoDrc * fGeo
Definition: PndDrcLutFill.h:57
TClonesArray * fBarPointArray
Definition: PndDrcLutFill.h:62
TString fOutputFile
Definition: PndDrcLutFill.h:85
TClonesArray * fEVPointArray
Definition: PndDrcLutFill.h:64
Double_t PipehAngle()
Definition: PndGeoDrc.h:139
Double_t
TClonesArray * fMCArray
Definition: PndDrcLutFill.h:61
TClonesArray * fPDPointArray
Definition: PndDrcLutFill.h:63
TClonesArray * fPDHitArray
Definition: PndDrcLutFill.h:66
Double_t fBboxNum
Definition: PndDrcLutFill.h:59
Double_t fDphi
Definition: PndDrcLutFill.h:59
TClonesArray * fDigiArray
Definition: PndDrcLutFill.h:65
TClonesArray * fLut[5]
Definition: PndDrcLutFill.h:67
void PndDrcLutFill::InitLut ( )

Definition at line 121 of file PndDrcLutFill.cxx.

References fLut, and n.

Referenced by Init().

122 {
123  Int_t Nnodes = 30000;
124  for(Int_t l=0; l<5; l++){
125  TClonesArray &fLuta = *fLut[l];
126  for (Long64_t n=0; n<Nnodes; n++) {
127  new((fLuta)[n]) PndDrcLutNode(-1);
128  }
129  }
130 }
int n
TClonesArray * fLut[5]
Definition: PndDrcLutFill.h:67
void PndDrcLutFill::ProcessPhotonHit ( )
private

Definition at line 143 of file PndDrcLutFill.cxx.

References At, Double_t, fBarPoint, fBarPointArray, fEVPoint, fEVPointArray, fLut, fMCArray, fMCTrack, fPDHit, fPDHitArray, fPDPoint, fPDPointArray, PndDrcBarPoint::GetBarId(), PndDrcPDPoint::GetBarPointID(), PndMCTrack::GetMomentum(), PndDrcEVPoint::GetNormal(), PndDrcPDHit::GetPosition(), PndDrcPDHit::GetSensorId(), i, Pi, and vec.

Referenced by Exec().

144 {
145  // LUT was generated for:
146  //Int_t lutboxId=3; //[R.K. 01/2017] unused variable?
147  Double_t lutboxPhi=10.825;
148 
149  //Int_t nofChPho = 0; //[R.K. 01/2017] unused variable?
150  Double_t path;//, barPhi; //[R.K. 01/2017] unused variable?
151  TVector3 dir, dirm, vec, posInBar;
152 
153  // Loop over PndDrcPDHits
154  for(Int_t k=0; k<fPDHitArray->GetEntriesFast(); k++) {
155  fPDHit = (PndDrcPDHit*)fPDHitArray->At(k);
156  Int_t pointID = fPDHit->GetLink(1).GetIndex();
157  Int_t sensorId = fPDHit->GetSensorId();
158  if(pointID==-1) continue;
159 
160  // Int_t digiID= fPDHit->GetRefIndex();
161  // fDigi = (PndDrcDigi*) fDigiArray->At(digiID);
162  // Int_t sensorId = fDigi->GetSensorId();
163  // Int_t pointID= fDigi->GetIndex(0);
164 
165  fPDPoint = (PndDrcPDPoint*)fPDPointArray->At(pointID);
167  Int_t barId = fBarPoint->GetBarId();
168 
169  Int_t trackID = fPDPoint->GetTrackID();
170  Double_t time = fPDPoint->GetTime();
171  Int_t nev=0;
172  path=0;
173  for(int i=0; i<fEVPointArray->GetEntriesFast(); i++){
175  if(trackID == fEVPoint->GetTrackID()){
176  // if(fEVPoint->GetDetectorID()==9375) {
177  // // get the direction of the phothon at the ent of the bar
178  // fEVPoint->Momentum(dirm);
179  // }
180  nev++;
181  vec = fEVPoint->GetNormal();
182  path += (vec.X()+vec.Y()*10 + vec.Z()*100)*1000*nev;
183  }
184  }
185 
186  fMCTrack = (PndMCTrack*)fMCArray->At(trackID);
187  dir = fMCTrack->GetMomentum().Unit();
188  dir.RotateZ(-lutboxPhi/180.*TMath::Pi());
189 
190  // sensorId = (sensorId/100 + lutboxId*17)*100 + sensorId%100;
191 
192  if(sensorId>30000 || sensorId<0) {
193  std::cout<<"WTQ fPDHit->GetDetectorID() "<<sensorId <<std::endl;
194  continue;
195  }
196 
197 
198  // //======================
199  // posInBar = fMCTrack->GetStartVertex();
200  // Double_t phi = posInBar.Phi()/TMath::Pi()*180;
201  // if(phi < 0) phi = 360 + phi;
202  // if(phi >= 0 && phi < 90) barPhi = TMath::Floor(phi/fDphi) *fDphi + fDphi/2.;
203  // if(phi >= 90 && phi < 270) barPhi = 90 + fPipehAngle + TMath::Floor((phi-90-fPipehAngle)/fDphi) *fDphi + fDphi/2.;
204  // if(phi >= 270 && phi < 360) barPhi = 270 + fPipehAngle + TMath::Floor((phi-270-fPipehAngle)/fDphi) *fDphi + fDphi/2.;
205 
206  // dirm.RotateZ(-barPhi/180.*TMath::Pi());
207  // TVector3 fnX1 = TVector3 (1,0,0);
208  // TVector3 fnY1 = TVector3( 0,1,0);
209  // double criticalAngle = asin(1.00028/fGeo->nQuartz());
210  // for(int u=0; u<8; u++){
211  // if(u == 0) dir = dirm;
212  // if(u == 1) dir.SetXYZ( dirm.X(), dirm.Y(),-dirm.Z());
213  // if(u == 2) dir.SetXYZ( dirm.X(),-dirm.Y(), dirm.Z());
214  // if(u == 3) dir.SetXYZ(-dirm.X(), dirm.Y(), dirm.Z());
215  // if(u == 4) dir.SetXYZ(-dirm.X(),-dirm.Y(), dirm.Z());
216  // if(u == 5) dir.SetXYZ(-dirm.X(), dirm.Y(),-dirm.Z());
217  // if(u == 6) dir.SetXYZ( dirm.X(),-dirm.Y(),-dirm.Z());
218  // if(u == 7) dir = -dirm;
219 
220  // if(dir.Angle(fnX1) < criticalAngle || dir.Angle(fnY1) < criticalAngle){
221  // //std::cout<<"Wrong combination "<<std::endl;
222  // continue;
223  // }
224  // ((PndDrcLutNode*)(fLut->At(sensorId)))->AddEntry(dir);
225  // ((PndDrcLutNode*)(fLut->At(sensorId)))->AddPathId(pathid);
226  // }
227  // //======================
228 
229  ((PndDrcLutNode*)(fLut[barId]->At(sensorId)))->AddEntry(sensorId, dir,path,0,time,fPDHit->GetPosition()); //fDigi->GetDetectorId()
230  }
231 }
PndDrcBarPoint * fBarPoint
Definition: PndDrcLutFill.h:73
Int_t i
Definition: run_full.C:25
Int_t GetSensorId()
Definition: PndDrcPDHit.h:59
Int_t GetBarId() const
TClonesArray * fBarPointArray
Definition: PndDrcLutFill.h:62
PndDrcPDHit * fPDHit
Definition: PndDrcLutFill.h:76
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
PndDrcPDPoint * fPDPoint
Definition: PndDrcLutFill.h:74
TClonesArray * fEVPointArray
Definition: PndDrcLutFill.h:64
Int_t GetBarPointID() const
Definition: PndDrcPDPoint.h:56
Double_t
PndMCTrack * fMCTrack
Definition: PndDrcLutFill.h:72
TVector3 GetNormal()
Definition: PndDrcEVPoint.h:49
TVector3 GetPosition() const
Definition: PndDrcPDHit.h:58
TClonesArray * fMCArray
Definition: PndDrcLutFill.h:61
PndDrcEVPoint * fEVPoint
Definition: PndDrcLutFill.h:77
TClonesArray * fPDPointArray
Definition: PndDrcLutFill.h:63
TClonesArray * fPDHitArray
Definition: PndDrcLutFill.h:66
Double_t Pi
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
dble_vec_t vec[12]
Definition: ranlxd.cxx:380
TClonesArray * fLut[5]
Definition: PndDrcLutFill.h:67
void PndDrcLutFill::SetDefaultParameters ( )
private
void PndDrcLutFill::SetOutputFile ( TString  outfilename = "luttab.root")
inline

Definition at line 51 of file PndDrcLutFill.h.

References fOutputFile.

51 {fOutputFile = outfilename;}
TString fOutputFile
Definition: PndDrcLutFill.h:85

Member Data Documentation

PndDrcBarPoint* PndDrcLutFill::fBarPoint
private

Definition at line 73 of file PndDrcLutFill.h.

Referenced by ProcessPhotonHit().

TClonesArray* PndDrcLutFill::fBarPointArray
private

Definition at line 62 of file PndDrcLutFill.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndDrcLutFill::fBboxNum
private

Definition at line 59 of file PndDrcLutFill.h.

Referenced by Init().

Int_t PndDrcLutFill::fDetectorID
private

Definition at line 58 of file PndDrcLutFill.h.

Referenced by Exec().

PndDrcDigi* PndDrcLutFill::fDigi
private

Definition at line 75 of file PndDrcLutFill.h.

TClonesArray* PndDrcLutFill::fDigiArray
private

Definition at line 65 of file PndDrcLutFill.h.

Referenced by Init().

Double_t PndDrcLutFill::fDphi
private

Definition at line 59 of file PndDrcLutFill.h.

Referenced by Init().

PndDrcEVPoint* PndDrcLutFill::fEVPoint
private

Definition at line 77 of file PndDrcLutFill.h.

Referenced by ProcessPhotonHit().

TClonesArray* PndDrcLutFill::fEVPointArray
private

Definition at line 64 of file PndDrcLutFill.h.

Referenced by Init(), and ProcessPhotonHit().

TFile* PndDrcLutFill::fFile
private

Definition at line 69 of file PndDrcLutFill.h.

Referenced by Finish(), and Init().

PndGeoDrc* PndDrcLutFill::fGeo
private

Definition at line 57 of file PndDrcLutFill.h.

Referenced by Init().

TClonesArray* PndDrcLutFill::fLut[5]
private

Definition at line 67 of file PndDrcLutFill.h.

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

TClonesArray* PndDrcLutFill::fMCArray
private

Definition at line 61 of file PndDrcLutFill.h.

Referenced by Init(), and ProcessPhotonHit().

PndMCTrack* PndDrcLutFill::fMCTrack
private

Definition at line 72 of file PndDrcLutFill.h.

Referenced by ProcessPhotonHit().

TString PndDrcLutFill::fOutputFile
private

Definition at line 85 of file PndDrcLutFill.h.

Referenced by Init(), PndDrcLutFill(), and SetOutputFile().

PndDrcPDHit* PndDrcLutFill::fPDHit
private

Definition at line 76 of file PndDrcLutFill.h.

Referenced by ProcessPhotonHit().

TClonesArray* PndDrcLutFill::fPDHitArray
private

Definition at line 66 of file PndDrcLutFill.h.

Referenced by Init(), and ProcessPhotonHit().

PndDrcPDPoint* PndDrcLutFill::fPDPoint
private

Definition at line 74 of file PndDrcLutFill.h.

Referenced by ProcessPhotonHit().

TClonesArray* PndDrcLutFill::fPDPointArray
private

Definition at line 63 of file PndDrcLutFill.h.

Referenced by Init(), and ProcessPhotonHit().

Double_t PndDrcLutFill::fPipehAngle
private

Definition at line 59 of file PndDrcLutFill.h.

Referenced by Init().

TTree* PndDrcLutFill::fTree
private

Definition at line 70 of file PndDrcLutFill.h.

Referenced by Finish(), and Init().

Int_t PndDrcLutFill::fVerbose
private

Definition at line 83 of file PndDrcLutFill.h.

Referenced by PndDrcLutFill().

Int_t PndDrcLutFill::nevents
private

Definition at line 84 of file PndDrcLutFill.h.

Referenced by Exec(), and Init().


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