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

#include <PndFtsHitProducerMcPointCoordinates.h>

Inheritance diagram for PndFtsHitProducerMcPointCoordinates:

Public Member Functions

 PndFtsHitProducerMcPointCoordinates ()
 
 ~PndFtsHitProducerMcPointCoordinates ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
PndFtsHitAddHit (Int_t detID, Int_t tubeID, Int_t chamberID, Int_t layerID, Int_t skew, Int_t iPoint, TVector3 &pos, TVector3 &dpos, Double_t p, Double_t rsim, Double_t closestDistanceError, Double_t depcharge)
 
PndFtsHitInfoAddHitInfo (Int_t fileNumber, Int_t eventNumber, Int_t trackID, Int_t pointID, Int_t nMerged, Bool_t isFake)
 
void FoldZPosWithResolution (Double_t &zpos, Double_t &zposError, TVector3 localInPos, TVector3 localOutPos)
 
Double_t GetError (Double_t)
 
void SetPersistence (Bool_t persistence)
 

Private Member Functions

void SetParContainers ()
 
 PndFtsHitProducerMcPointCoordinates (const PndFtsHitProducerMcPointCoordinates &L)
 
PndFtsHitProducerMcPointCoordinatesoperator= (const PndFtsHitProducerMcPointCoordinates &)
 
 ClassDef (PndFtsHitProducerMcPointCoordinates, 1)
 

Private Attributes

TClonesArray * fPointArray
 
TClonesArray * fHitArray
 
TObjArray * fVolumeArray
 
TClonesArray * fHitInfoArray
 
Int_t fevtn
 
PndGeoFtsParfFtsParameters
 
Bool_t fPersistence
 
TClonesArray * fTubeArray
 

Detailed Description

Definition at line 17 of file PndFtsHitProducerMcPointCoordinates.h.

Constructor & Destructor Documentation

PndFtsHitProducerMcPointCoordinates::PndFtsHitProducerMcPointCoordinates ( )
PndFtsHitProducerMcPointCoordinates::~PndFtsHitProducerMcPointCoordinates ( )

Destructor

Definition at line 55 of file PndFtsHitProducerMcPointCoordinates.cxx.

55 { }
PndFtsHitProducerMcPointCoordinates::PndFtsHitProducerMcPointCoordinates ( const PndFtsHitProducerMcPointCoordinates L)
private

Member Function Documentation

PndFtsHit * PndFtsHitProducerMcPointCoordinates::AddHit ( Int_t  detID,
Int_t  tubeID,
Int_t  chamberID,
Int_t  layerID,
Int_t  skew,
Int_t  iPoint,
TVector3 &  pos,
TVector3 &  dpos,
Double_t  p,
Double_t  rsim,
Double_t  closestDistanceError,
Double_t  depcharge 
)

Definition at line 269 of file PndFtsHitProducerMcPointCoordinates.cxx.

References fHitArray.

Referenced by Exec().

269  {
270 
271 
272  // see PndFtsHit for hit description
273  TClonesArray& clref = *fHitArray;
274  Int_t size = clref.GetEntriesFast();
275 
276  PndFtsHit *hitnew = new(clref[size]) PndFtsHit(detID, tubeID, chamberID, layerID, skew, iPoint, pos, dpos, p, rsim, closestDistanceError, depcharge);
277  return hitnew;
278 
279 }
TVector3 pos
Double_t p
Definition: anasim.C:58
double skew
PndFtsHitInfo * PndFtsHitProducerMcPointCoordinates::AddHitInfo ( Int_t  fileNumber,
Int_t  eventNumber,
Int_t  trackID,
Int_t  pointID,
Int_t  nMerged,
Bool_t  isFake 
)

Definition at line 283 of file PndFtsHitProducerMcPointCoordinates.cxx.

References fHitInfoArray.

Referenced by Exec().

283  {
284  // see PndFtsHitInfo for hit description
285 
286 
287  TClonesArray& clref = *fHitInfoArray;
288  Int_t size = clref.GetEntriesFast();
289  return new(clref[size]) PndFtsHitInfo(fileNumber, eventNumber, trackID, pointID, nMerged, isFake);
290 }
PndFtsHitProducerMcPointCoordinates::ClassDef ( PndFtsHitProducerMcPointCoordinates  ,
 
)
private
void PndFtsHitProducerMcPointCoordinates::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 116 of file PndFtsHitProducerMcPointCoordinates.cxx.

References AddHit(), AddHitInfo(), Double_t, PndFtsSingleStraw::FastPartToADC(), PndFtsSingleStraw::FastRec(), fevtn, fHitArray, fHitInfoArray, fPointArray, fVerbose, PndFtsPoint::GetChamberID(), GetError(), PndFtsPoint::GetLayerID(), PndFtsPoint::GetMass(), PndFtsPoint::GetPxOut(), PndFtsPoint::GetPyOut(), PndFtsPoint::GetPzOut(), PndFtsPoint::GetTubeID(), PndFtsPoint::GetXInLocal(), PndFtsPoint::GetXOutLocal(), PndFtsPoint::GetYInLocal(), PndFtsPoint::GetYOutLocal(), PndFtsPoint::GetZInLocal(), PndFtsPoint::GetZOutLocal(), point, pos, PndFtsSingleStraw::PutWireXYZ(), skew, PndFtsSingleStraw::TConst(), PndFtsSingleStraw::TInit(), and PndFtsSingleStraw::TrueDist().

116  {
117 
118  //std::cout<<"PndFtsHitProducer Exec ########"<<std::endl;
119  if(fVerbose && fevtn%50==0) cout << "Event Number "<<fevtn<<endl;
120  else if(fVerbose >= 3) cout << "Event Number "<<fevtn<<endl;
121 
122  fevtn++;
123 
124  // Reset output array
125  if ( ! fHitArray ) Fatal("Exec", "No HitArray");
126 
127  fHitArray->Clear();
128  fHitInfoArray->Clear();
129 
130  Int_t detID = 0; // detectorID
131  TVector3 pos, dpos; // position and error vectors
132 
133  // Declare some variables
134  PndFtsPoint* point = NULL;
135 
136  // Loop over FtsPoints
137  Int_t nPoints = fPointArray->GetEntriesFast();
138 
139  // cout << "------------ " << nPoints << endl;
140 
141  for (Int_t iPoint = 0; iPoint < nPoints; iPoint++) {
142  point = (PndFtsPoint*) fPointArray->At(iPoint);
143  if (point == NULL) continue;
144 
145  detID = point->GetDetectorID();
146 
147  // tubeID CHECK added
148  Int_t skew = 0;
149  Int_t tubeID = point->GetTubeID();
150  Int_t chamberID=point->GetChamberID();
151  Int_t layerID=point->GetLayerID();
152  //PndFtsTube *tube = (PndFtsTube*) fTubeArray->At(tubeID); //[R.K. 01/2017] unused variable?
153 
154  //if skewed tube: skew==1
155  if(layerID>=3 && layerID<=6){skew=1;} //skewed tudes fts1
156  if(layerID>=11 && layerID<=14){skew=1;} //skewed tudes fts2
157  if(layerID>=19 && layerID<=22){skew=1;} //skewed tudes fts3
158  if(layerID>=27 && layerID<=30){skew=1;} //skewed tudes fts4
159  if(layerID>=35 && layerID<=38){skew=1;} //skewed tudes fts5
160  if(layerID>=43 && layerID<=46){skew=1;} //skewed tudes fts6
161 
162 
163  double InOut[6];
164  memset(InOut, 0, sizeof(InOut));
165 
166  InOut[0] = point->GetXInLocal();
167  InOut[1] = point->GetYInLocal();
168  InOut[2] = point->GetZInLocal();
169  InOut[3] = point->GetXOutLocal();
170  InOut[4] = point->GetYOutLocal();
171  InOut[5] = point->GetZOutLocal();
172 
173  // single straw tube simulation -----------------------
174  PndFtsSingleStraw fts;
175 
176  //setting the single straw tube simulation constants
177  // 3 options currently available:
178  // TConst(tube radius (cm), gas pressure (bar), Ar%, CO2%)
179  // fts.TConst(0.4, 1, 0.9, 0.1);
180  //fts.TConst(0.5, 1, 0.9, 0.1);//1 bar
181  fts.TConst(0.5, 2, 0.8, 0.2); //2 bar
182 
183  // wire positioning->controllare bene
184  fts.PutWireXYZ(0., 0., -75., 0., 0., 75.);
185 
186  // get particle momentum
187  TVector3 momentum(point->GetPxOut(),point->GetPyOut(),point->GetPzOut()); // GeV/c
188 
189  Double_t GeV=1.;
190  // position in cm (already in cm); momentum in GeV (already in GeV); mass in GeV (already in GeV)
191 
192  // drift time calculation
193 
194  Double_t pulset =-1;
195  //pulset = fts.PartToTime(point->GetMass()/GeV, momentum.Mag()/GeV, InOut);
196 
197  // constant initialization
198  fts.TInit(point->GetMass()/GeV, momentum.Mag()/GeV, InOut);
199 
200  // true radius (cm)
201  Double_t true_rad = fts.TrueDist(InOut);
202 
203  // simulated radius (cm)
204  //Double_t radius = fts.TimnsToDiscm(pulset);
205  //if(radius < 0.) radius = 0.; // CHECK
206  //if(radius <0. ||radius==0.) radius =-999;
207 
208  // fast simulation
209  Double_t radius= fts.FastRec(true_rad,1) ; //,0) standard curve ,1) Juelich exp curve
210  //Juelich is at 2 bar pressure
211  // dE calculation
212  // double depCharge = fts.PartToADC();
213 
214  // dE calculation ------- check
215  // charge calculation
216  Double_t depcharge = fts.FastPartToADC(); // CHECK arbitrary units!
217  // dE/dx calculation postponed
218  //Double_t dedx = -999; //[R.K. 01/2017] unused variable?
219 
220  // stt2: detID, pos, dpos, index come from --------------
221  // stt2 (FairHit):
222  Double_t closestDistanceError = GetError(radius);//calculates the error according to Juelich experimental curves
223  //cout<<"radius "<<radius<<" error "<<closestDistanceError<<endl;
224  //closestDistanceError = 0.0150; //150 microns check this point!
225  //closestDistanceError =TMath::Sqrt(2.)*radius/TMath::Sqrt(12);
226 
227 
228  TVector3 position(point->GetX(), point->GetY(), point->GetZ()); // use this for hits having same coordinates as MC points
229  //TVector3 position = tube->GetPosition(); // use this for realistic hit production
230  pos.SetXYZ(position.X(), position.Y(), position.Z()); // <--- stt1
231 
232  // dpos.SetXYZ(innerStrawDiameter / 2., innerStrawDiameter / 2., GetLongitudinalResolution(position.Z()));
233  dpos.SetXYZ(0.5, 0.5, 3.); // per adesso (stessi che in Ideal:
234  // innerStrawDiameter/2 = 0.5,
235  // longitudinalResolution = 3.)
236 
237  // create hit
238  AddHit(detID, tubeID, chamberID, layerID, skew, iPoint, pos, dpos, pulset, radius, closestDistanceError, depcharge);
239 
240  AddHitInfo(0, 0, point->GetTrackID(), iPoint, 0, kFALSE);
241 
242  }// Loop over MCPoints
243 
244 
245  // Event summary
246  //cout << "-I- PndSttHitProducerRealFast: " << nPoints << " FtsPoints, "
247  //<< nPoints << " Hits created." << endl;
248 
249 }
Double_t GetXOutLocal() const
Definition: PndFtsPoint.h:59
TVector3 pos
Int_t GetChamberID()
Definition: PndFtsPoint.h:89
Double_t GetPzOut() const
Definition: PndFtsPoint.h:54
int fVerbose
Definition: poormantracks.C:24
Double_t GetZOutLocal() const
Definition: PndFtsPoint.h:61
Double_t GetPxOut() const
Definition: PndFtsPoint.h:52
void TInit(Double_t Mass, Double_t Momentum, Double_t InOut[])
Double_t GetYOutLocal() const
Definition: PndFtsPoint.h:60
Double_t GetYInLocal() const
Definition: PndFtsPoint.h:63
Double_t GetZInLocal() const
Definition: PndFtsPoint.h:64
double skew
PndFtsHit * AddHit(Int_t detID, Int_t tubeID, Int_t chamberID, Int_t layerID, Int_t skew, Int_t iPoint, TVector3 &pos, TVector3 &dpos, Double_t p, Double_t rsim, Double_t closestDistanceError, Double_t depcharge)
void TConst(Double_t Radius, Double_t pSTP, Double_t ArP, Double_t CO2P)
Int_t GetTubeID()
Definition: PndFtsPoint.h:86
Double_t TrueDist(Double_t Point[])
Double_t
void PutWireXYZ(Double_t w1, Double_t w2, Double_t w3, Double_t w4, Double_t w5, Double_t w6)
PndFtsHitInfo * AddHitInfo(Int_t fileNumber, Int_t eventNumber, Int_t trackID, Int_t pointID, Int_t nMerged, Bool_t isFake)
Double_t FastRec(Double_t TrueDcm, Int_t Flag)
Double_t GetPyOut() const
Definition: PndFtsPoint.h:53
Double_t GetXInLocal() const
Definition: PndFtsPoint.h:62
Int_t GetLayerID()
Definition: PndFtsPoint.h:90
PndSdsMCPoint * point
Definition: anaLmdCluster.C:72
Double_t GetMass() const
Definition: PndFtsPoint.h:70
void PndFtsHitProducerMcPointCoordinates::FoldZPosWithResolution ( Double_t zpos,
Double_t zposError,
TVector3  localInPos,
TVector3  localOutPos 
)

Definition at line 251 of file PndFtsHitProducerMcPointCoordinates.cxx.

253 {
254 
255 
256  //Double_t //[R.K. 01/2017] unused variable
257  //zPosInStrawFrame = (localOutPos.Z() - localInPos.Z()) / 2.; //[R.K. 01/2017] unused variable?
258 
259  // zposError = gRandom->Gaus(0., GetLongitudinalResolution(zPosInStrawFrame));
260  zposError = gRandom->Gaus(0., 3.); // per adesso (stesso che in Ideal:
261  // longitudinalResolution = 3.)
262 
263  zpos += zposError;
264 }
const Double_t zpos
Double_t PndFtsHitProducerMcPointCoordinates::GetError ( Double_t  TrueDcm)

Definition at line 293 of file PndFtsHitProducerMcPointCoordinates.cxx.

References Double_t.

Referenced by Exec().

293  {
294 
295 
296  // data from julich
297  Double_t resmic=-1;
298  if(TrueDcm < 0.48){
299  resmic = 20. +1.48048e+02
300  -3.35951e+02*TrueDcm
301  -1.87575e+03*pow(TrueDcm,2)
302  +1.92910e+04*pow(TrueDcm,3)
303  -6.90036e+04*pow(TrueDcm,4)
304  +1.07960e+05*pow(TrueDcm,5)
305  -5.90064e+04*pow(TrueDcm,6) ;
306  }
307  else resmic=65.;
308 
309  return resmic*0.0001;
310 }
Double_t
InitStatus PndFtsHitProducerMcPointCoordinates::Init ( )
virtual

Virtual method Init

Definition at line 61 of file PndFtsHitProducerMcPointCoordinates.cxx.

References fevtn, fFtsParameters, fHitArray, fHitInfoArray, PndFtsMapCreator::FillTubeArray(), fPersistence, fPointArray, fTubeArray, fVolumeArray, and gGeoManager.

61  {
62  fevtn=0;
63 
64  std::cout<<"#########################################################"<<std::endl;
65  std::cout<<"PndFtsHitProducerMcPointCoordinates: Init()#######"<<std::endl;
66  std::cout<<"FTS hits are produced with MC point coordinates for hit positions instead of center of straws!\n";
67  std::cout<<"Use this class only for debugging or parameter optimisation purposes!\n\n";
68  std::cout<<"#########################################################"<<std::endl;
69  // Get RootManager
70  FairRootManager* ioman = FairRootManager::Instance();
71  if ( ! ioman ) {
72  cout << "-E- PndFtsHitProducerMcPointCoordinates-wintz::Init: "
73  << "RootManager not instantiated!" << endl;
74  return kFATAL;
75  }
76 
77  // Get input array
78  fPointArray = (TClonesArray*) ioman->GetObject("FTSPoint");
79  if ( ! fPointArray ) {
80  cout << "-W- PndFtsHitProducerMcPointCoordinates::Init: "
81  << "No FTSPoint array!" << endl;
82  return kERROR;
83  }
84 
85  // Create and register output array
86  fHitArray = new TClonesArray("PndFtsHit");
87  ioman->Register("FTSHit","FTS",fHitArray, fPersistence);
88 
89  // Create and register output array
90  fHitInfoArray = new TClonesArray("PndFtsHitInfo");
91  ioman->Register("FTSHitInfo", "FTS", fHitInfoArray, kFALSE);
92 
93  fVolumeArray = gGeoManager->GetListOfVolumes();
94 
95  cout << "-I- PndFtsHitProducerMcPointCoordinates: INITIALIZATION SUCCESSFUL" << endl;
96 
97 
98  //CHECK added
100  fTubeArray = mapper->FillTubeArray();
101  delete mapper;
102 
103  return kSUCCESS;
104 
105 }
TGeoManager * gGeoManager
TClonesArray * FillTubeArray()
this function will be used in PndFtsHitProducesRealFast
PndFtsHitProducerMcPointCoordinates& PndFtsHitProducerMcPointCoordinates::operator= ( const PndFtsHitProducerMcPointCoordinates )
inlineprivate

Definition at line 76 of file PndFtsHitProducerMcPointCoordinates.h.

76 {return *this;}
void PndFtsHitProducerMcPointCoordinates::SetParContainers ( )
private

Definition at line 108 of file PndFtsHitProducerMcPointCoordinates.cxx.

References fFtsParameters, and rtdb.

108  {
109 
110  FairRuntimeDb* rtdb = FairRunAna::Instance()->GetRuntimeDb();
111  fFtsParameters = (PndGeoFtsPar*) rtdb->getContainer("PndGeoFtsPar");
112 }
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
void PndFtsHitProducerMcPointCoordinates::SetPersistence ( Bool_t  persistence)
inline

set persistence flag

Definition at line 47 of file PndFtsHitProducerMcPointCoordinates.h.

References fPersistence.

Member Data Documentation

Int_t PndFtsHitProducerMcPointCoordinates::fevtn
private

Definition at line 64 of file PndFtsHitProducerMcPointCoordinates.h.

Referenced by Exec(), and Init().

PndGeoFtsPar* PndFtsHitProducerMcPointCoordinates::fFtsParameters
private

Definition at line 66 of file PndFtsHitProducerMcPointCoordinates.h.

Referenced by Init(), and SetParContainers().

TClonesArray* PndFtsHitProducerMcPointCoordinates::fHitArray
private

Output array of PndFtsHits

Definition at line 58 of file PndFtsHitProducerMcPointCoordinates.h.

Referenced by AddHit(), Exec(), and Init().

TClonesArray* PndFtsHitProducerMcPointCoordinates::fHitInfoArray
private

Output array of PndFtsHitInfo

Definition at line 63 of file PndFtsHitProducerMcPointCoordinates.h.

Referenced by AddHitInfo(), Exec(), and Init().

Bool_t PndFtsHitProducerMcPointCoordinates::fPersistence
private

object persistence

Definition at line 69 of file PndFtsHitProducerMcPointCoordinates.h.

Referenced by Init(), and SetPersistence().

TClonesArray* PndFtsHitProducerMcPointCoordinates::fPointArray
private

Input array of PndFtsPoints

Definition at line 55 of file PndFtsHitProducerMcPointCoordinates.h.

Referenced by Exec(), and Init().

TClonesArray* PndFtsHitProducerMcPointCoordinates::fTubeArray
private

Definition at line 72 of file PndFtsHitProducerMcPointCoordinates.h.

Referenced by Init().

TObjArray* PndFtsHitProducerMcPointCoordinates::fVolumeArray
private

Definition at line 60 of file PndFtsHitProducerMcPointCoordinates.h.

Referenced by Init().


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