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

Class for conversion points to hits. More...

#include <PndGemIdealHitProducer.h>

Inheritance diagram for PndGemIdealHitProducer:

Public Member Functions

 PndGemIdealHitProducer ()
 
 PndGemIdealHitProducer (const char *name, Int_t iVerbose)
 
virtual ~PndGemIdealHitProducer ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 

Private Member Functions

virtual void SetParContainers ()
 
virtual void Finish ()
 
 ClassDef (PndGemIdealHitProducer, 1)
 

Private Attributes

PndGemDigiParfDigiPar
 
Int_t fTNofEvents
 
Int_t fTNofPoints
 
Int_t fTNofHits
 
TClonesArray * fPointArray
 Input array of PndGemMCPoints. More...
 
TClonesArray * fHitArray
 Output array of PndGemHits. More...
 

Detailed Description

Class for conversion points to hits.

Author
R.Karabowicz r.kar.nosp@m.abow.nosp@m.icz@g.nosp@m.si.d.nosp@m.e
Date
17/03/2009 is a dummy hit reconstruction task, which converts each point into a hit (1:1 conversion) with a resolution defined by strip pitches and angles taken from the parameter file

Definition at line 30 of file PndGemIdealHitProducer.h.

Constructor & Destructor Documentation

PndGemIdealHitProducer::PndGemIdealHitProducer ( )

Default constructor

Definition at line 50 of file PndGemIdealHitProducer.cxx.

References fDigiPar, fTNofEvents, fTNofHits, and fTNofPoints.

50  :
51  FairTask("Ideal GEM hit Producer") {
52  fDigiPar = NULL;
53 
54  fTNofEvents = 0;
55  fTNofPoints = 0;
56  fTNofHits = 0;
57 }
PndGemIdealHitProducer::PndGemIdealHitProducer ( const char *  name,
Int_t  iVerbose 
)

Constructor

Definition at line 60 of file PndGemIdealHitProducer.cxx.

References fDigiPar, fTNofEvents, fTNofHits, and fTNofPoints.

61  : FairTask(name, iVerbose) {
62  fDigiPar = NULL;
63 
64  fTNofEvents = 0;
65  fTNofPoints = 0;
66  fTNofHits = 0;
67 }
TString name
Int_t iVerbose
PndGemIdealHitProducer::~PndGemIdealHitProducer ( )
virtual

Destructor

Definition at line 71 of file PndGemIdealHitProducer.cxx.

References fDigiPar, and fHitArray.

71  {
72  if(fHitArray){
73  fHitArray->Delete();
74  delete fHitArray;
75  }
76  if ( fDigiPar) delete fDigiPar;
77 }
TClonesArray * fHitArray
Output array of PndGemHits.

Member Function Documentation

PndGemIdealHitProducer::ClassDef ( PndGemIdealHitProducer  ,
 
)
private
void PndGemIdealHitProducer::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 128 of file PndGemIdealHitProducer.cxx.

References CAMath::Cos(), Double_t, fDigiPar, fHitArray, fPointArray, fTNofEvents, fTNofHits, fTNofPoints, PndGemSensor::GetD(), PndGemSensor::GetInnerRadius(), PndGemDigiPar::GetNodeName(), PndGemSensor::GetOuterRadius(), PndGemSensor::GetPitch(), PndGemDigiPar::GetSensor(), PndGemMCPoint::GetSensorId(), PndGemSensor::GetStripAngle(), PndGemSensor::GetZ0(), gGeoManager, kGemHit, Pi, pos, sensor, sigma, CAMath::Sin(), CAMath::Sqrt(), CAMath::Tan(), and TString.

128  {
129 
130  // Reset output array
131  if( !fHitArray ) Fatal("Exec", "No hit array");
132 
133  fTNofEvents++;
134 
136 
137  Int_t nofHits = 0;
138  fHitArray->Delete();
139  Double_t dr, dp;
140  Double_t radius, innerR;
141 
142  for ( Int_t iPoint = 0 ; iPoint < fPointArray->GetEntriesFast() ; iPoint++ ) {
143  PndGemMCPoint* currentPndGemMCPoint = (PndGemMCPoint*)fPointArray->At(iPoint);
144 
145  Double_t posIn[3] = {currentPndGemMCPoint->GetX(),
146  currentPndGemMCPoint->GetY(),
147  currentPndGemMCPoint->GetZ()};
148 
149  Int_t sensorId = currentPndGemMCPoint->GetSensorId();
150 
151  TString nodeName = fDigiPar->GetNodeName(sensorId);
152 
153  gGeoManager->cd(nodeName.Data());
154  TGeoNode* curNode = gGeoManager->GetCurrentNode();
155 
156  sensor = (PndGemSensor*)fDigiPar->GetSensor(sensorId);
157 
158  fTNofPoints++;
159 
160  Double_t locPosIn[4];
161 
162  curNode->MasterToLocal(posIn,locPosIn);
163 
164  radius = TMath::Sqrt(locPosIn[0]*locPosIn[0]+locPosIn[1]*locPosIn[1]);
165  innerR = sensor->GetInnerRadius();
166  if ( sensor->GetStripAngle(0) == 0. && sensor->GetStripAngle(1) == 90. ) {
167  dp = sensor->GetPitch(0)*radius/innerR/TMath::Sqrt(12.);
168  dr = sensor->GetPitch(1)/TMath::Sqrt(12.);
169  }
170 
171  if ( sensor->GetStripAngle(0) == -sensor->GetStripAngle(1) ) {
172  Double_t pitch = ( sensor->GetPitch(0)+sensor->GetPitch(1) ) / 2.;
173  dp = pitch*radius/innerR/TMath::Sqrt(12.);
174  Double_t cosAng = TMath::Cos(TMath::DegToRad()*sensor->GetStripAngle(0));
175  Double_t xh = - innerR*cosAng + TMath::Sqrt(innerR*innerR*cosAng*cosAng-innerR*innerR+radius*radius);
176  Double_t backAng = TMath::ACos((radius*radius+xh*xh-innerR*innerR)/(2.*radius*xh));
177  dr = dp / TMath::Tan(backAng);
178  }
179 
180  Double_t phiAValue = TMath::ATan(locPosIn[0]/locPosIn[1]);
181  if ( locPosIn[1] < 0 ) phiAValue += TMath::Pi();
182  else if ( locPosIn[0] < 0 ) phiAValue += 2.*TMath::Pi();
183 
184  Double_t rSmear, pSmear;
185 
186  gRandom->Rannor(rSmear,pSmear);
187  rSmear = radius + rSmear*dr;
188  pSmear = 0. + pSmear*dp;
189 
190  if ( rSmear < innerR || rSmear > sensor->GetOuterRadius() ) rSmear = radius;
191 
192  locPosIn[0] = rSmear*TMath::Sin(phiAValue) + pSmear*TMath::Cos(phiAValue);
193  locPosIn[1] = rSmear*TMath::Cos(phiAValue) - pSmear*TMath::Sin(phiAValue);
194 
195  TVector3 pos(locPosIn[0],locPosIn[1],sensor->GetZ0());
196 
197  Double_t sigma = dp;
198  if ( dr > sigma ) sigma = dr;
199  TVector3 dpos(sigma, sigma, sensor->GetD());
200 
201  new ((*fHitArray)[nofHits++]) PndGemHit(kGemHit,
202  pos, dpos,
203  -1, -1, dr, dp, iPoint);
204  fTNofHits++;
205  } // end of loop over Points
206 }
TVector3 pos
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
Double_t GetInnerRadius() const
Definition: PndGemSensor.h:102
static T Sin(const T &x)
Definition: PndCAMath.h:42
float Tan(float x)
Definition: PndCAMath.h:165
TGeoVolume * sensor
Double_t sigma[nsteps]
Definition: dedx_bands.C:65
TClonesArray * fHitArray
Output array of PndGemHits.
TGeoManager * gGeoManager
static T Cos(const T &x)
Definition: PndCAMath.h:43
PndGemSensor * GetSensor(Int_t stationNr, Int_t sensorNr)
Int_t GetSensorId() const
Definition: PndGemMCPoint.h:90
Double_t
TString GetNodeName(Int_t sensorId)
TClonesArray * fPointArray
Input array of PndGemMCPoints.
Double_t GetD() const
Definition: PndGemSensor.h:104
Double_t GetZ0() const
Definition: PndGemSensor.h:100
Double_t GetOuterRadius() const
Definition: PndGemSensor.h:103
Double_t GetPitch(Int_t index) const
Definition: PndGemSensor.h:106
Double_t Pi
Double_t GetStripAngle(Int_t index) const
Definition: PndGemSensor.h:105
void PndGemIdealHitProducer::Finish ( )
privatevirtual

Virtual method Finish

Definition at line 209 of file PndGemIdealHitProducer.cxx.

References Double_t, fDigiPar, fHitArray, fTNofEvents, fTNofHits, fTNofPoints, and PndGemDigiPar::GetNSensors().

209  {
210  if ( fHitArray ) fHitArray->Delete();
211 
212  cout << "-------------------- " << fName.Data() << " : Summary ---------------" << endl;
213  cout << " Events: " << setw(10) << fTNofEvents << endl;
214  cout << " Points: " << setw(10) << fTNofPoints << " ( " << (Double_t)fTNofPoints/((Double_t)fTNofEvents) << " per event )" << endl;
215  cout << " Hits: " << setw(10) << fTNofHits << " ( " << (Double_t)fTNofHits /((Double_t)fTNofEvents) << " per event )" << endl;
216  cout << " --> ( " << (Double_t)fTNofHits /((Double_t)fTNofEvents)/((Double_t)fDigiPar->GetNSensors()) << " per sensor )" << endl;
217  cout << " --> ( " << (Double_t)fTNofHits /((Double_t)fTNofPoints ) << " per point )" << endl;
218  cout << "---------------------------------------------------------------------" << endl;
219 }
TGeoVolume * sensor
TClonesArray * fHitArray
Output array of PndGemHits.
Double_t
TClonesArray * point
Definition: anaLmdDigi.C:29
InitStatus PndGemIdealHitProducer::Init ( )
virtual

Virtual method Init

Definition at line 97 of file PndGemIdealHitProducer.cxx.

References fHitArray, and fPointArray.

97  {
98 
99  std::cout << " INITIALIZATION OF Ideal Gem Hit Producer***************" << std::endl;
100 
101  // Get RootManager
102  FairRootManager* ioman = FairRootManager::Instance();
103  if( !ioman ) {
104  std::cout << "-E- PndGemIdealHitProducer::Init: "
105  << "RootManager not instantiated!" << std::endl;
106  return kFATAL;
107  }
108 
109  // Get input array
110  fPointArray = (TClonesArray*) ioman->GetObject("GEMPoint");
111  if( !fPointArray ) {
112  std::cout << "-W- PndGemIdealHitProducer::Init: "
113  << "Array of GEMPoints not found!" << std::endl;
114  return kERROR;
115  }
116 
117  // Create and register output array
118  fHitArray = new TClonesArray("PndGemHit");
119 
120  ioman->Register("GEMHit","GEM ideal hit",fHitArray,kTRUE);
121 
122  std::cout << "-I- PndGemIdealHitProducer: Intialization successfull" << std::endl;
123  return kSUCCESS;
124 
125 }
TClonesArray * fHitArray
Output array of PndGemHits.
TClonesArray * fPointArray
Input array of PndGemMCPoints.
void PndGemIdealHitProducer::SetParContainers ( )
privatevirtual

Private method Smear(...) smears the position vectorGet parameter containers

Definition at line 80 of file PndGemIdealHitProducer.cxx.

References fDigiPar, and run.

80  {
81 
82  // Get run and runtime database
83  FairRunAna* run = FairRunAna::Instance();
84  if ( ! run ) Fatal("SetParContainers", "No analysis run");
85 
86  FairRuntimeDb* db = run->GetRuntimeDb();
87  if ( ! db ) Fatal("SetParContainers", "No runtime database");
88 
89  // Get GEM digitisation parameter container
90  fDigiPar = (PndGemDigiPar*)(db->getContainer("PndGemDetectors"));
91 
92 }
Int_t run
Definition: autocutx.C:47
Digitization Parameter Class for GEM part.
Definition: PndGemDigiPar.h:31

Member Data Documentation

PndGemDigiPar* PndGemIdealHitProducer::fDigiPar
private

Public method AddHit(...) adds another hit to the hit array

Definition at line 56 of file PndGemIdealHitProducer.h.

Referenced by Exec(), Finish(), PndGemIdealHitProducer(), SetParContainers(), and ~PndGemIdealHitProducer().

TClonesArray* PndGemIdealHitProducer::fHitArray
private

Output array of PndGemHits.

Definition at line 64 of file PndGemIdealHitProducer.h.

Referenced by Exec(), Finish(), Init(), and ~PndGemIdealHitProducer().

TClonesArray* PndGemIdealHitProducer::fPointArray
private

Input array of PndGemMCPoints.

Definition at line 63 of file PndGemIdealHitProducer.h.

Referenced by Exec(), and Init().

Int_t PndGemIdealHitProducer::fTNofEvents
private

Definition at line 59 of file PndGemIdealHitProducer.h.

Referenced by Exec(), Finish(), and PndGemIdealHitProducer().

Int_t PndGemIdealHitProducer::fTNofHits
private

Definition at line 61 of file PndGemIdealHitProducer.h.

Referenced by Exec(), Finish(), and PndGemIdealHitProducer().

Int_t PndGemIdealHitProducer::fTNofPoints
private

Definition at line 60 of file PndGemIdealHitProducer.h.

Referenced by Exec(), Finish(), and PndGemIdealHitProducer().


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