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

#include <PndGemSmearingTask.h>

Inheritance diagram for PndGemSmearingTask:

Public Member Functions

 PndGemSmearingTask ()
 
 PndGemSmearingTask (Double_t sx, Double_t sy, Double_t sz)
 
virtual ~PndGemSmearingTask ()
 
virtual void SetParContainers ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 

Private Member Functions

void InitTransMat ()
 
void smear (TVector3 &pos)
 
void smearLocal (TVector3 &pos)
 
void CalcGFDetPlane (TVector3 &o, TVector3 &u, TVector3 &v)
 
void Register ()
 
void Reset ()
 
void ProduceHits ()
 
 ClassDef (PndGemSmearingTask, 1)
 

Private Attributes

PndGemDigiParfDigiPar
 
TString fBranchName
 
TClonesArray * fPointArray
 
TClonesArray * fMctruthArray
 
TClonesArray * fHitOutputArray
 
Double_t fSigmaX
 
Double_t fSigmaY
 
Double_t fSigmaZ
 
PndGemMCPointfCurrentPndGemMCPoint
 
TGeoHMatrix * fCurrentTransMat
 
TMatrixT< Double_tfHitCovMatrix
 

Detailed Description

Definition at line 33 of file PndGemSmearingTask.h.

Constructor & Destructor Documentation

PndGemSmearingTask::PndGemSmearingTask ( )

Default constructor

Definition at line 29 of file PndGemSmearingTask.cxx.

References fBranchName, fDigiPar, fSigmaX, fSigmaY, and fSigmaZ.

29  :
30  FairTask("Ideal reconstruction task for PANDA PndGem"),
31  fHitCovMatrix(3,3)
32 {
33  fDigiPar = NULL;
34  fSigmaX=0.;
35  fSigmaY=0.;
36  fSigmaZ=0.;
37  fBranchName = "GEMPoint";
38 }
TMatrixT< Double_t > fHitCovMatrix
PndGemDigiPar * fDigiPar
PndGemSmearingTask::PndGemSmearingTask ( Double_t  sx,
Double_t  sy,
Double_t  sz 
)

Definition at line 42 of file PndGemSmearingTask.cxx.

References fBranchName, fDigiPar, fSigmaX, fSigmaY, and fSigmaZ.

42  :
43  FairTask("Ideal reconstruction task for PANDA PndGem"),
44  fHitCovMatrix(3,3)
45 {
46  fDigiPar = NULL;
47  fSigmaX=sx;
48  fSigmaY=sy;
49  fSigmaZ=sz;
50  fBranchName = "GEMPoint";
51 }
TMatrixT< Double_t > fHitCovMatrix
PndGemDigiPar * fDigiPar
PndGemSmearingTask::~PndGemSmearingTask ( )
virtual

Destructor

Definition at line 56 of file PndGemSmearingTask.cxx.

References fDigiPar.

57 {
58  if ( fDigiPar) delete fDigiPar;
59 }
PndGemDigiPar * fDigiPar

Member Function Documentation

void PndGemSmearingTask::CalcGFDetPlane ( TVector3 &  o,
TVector3 &  u,
TVector3 &  v 
)
private

Definition at line 219 of file PndGemSmearingTask.cxx.

References Double_t, fCurrentPndGemMCPoint, fCurrentTransMat, fVerbose, PndGemMCPoint::GetSensorId(), and v.

Referenced by Exec().

220 {
221  Double_t O[3], U[3], V[3], o[3], u[3], v[3];
222  O[0]=oVect.x(); O[1]=oVect.y(); O[2]=oVect.z();
223  U[0]=uVect.x(); U[1]=uVect.y(); U[2]=uVect.z();
224  V[0]=vVect.x(); V[1]=vVect.y(); V[2]=vVect.z();
225 
226  if (fVerbose > 1) {
227  std::cout<<"PndGemSmearingTask::CalcGFDetPlane from Detector sensorId "
228  <<fCurrentPndGemMCPoint->GetSensorId()<<std::endl;
229  }
230  //make transformation
231  fCurrentTransMat->LocalToMaster(O,o);
232  fCurrentTransMat->LocalToMaster(U,u);
233  fCurrentTransMat->LocalToMaster(V,v);
234  oVect.SetXYZ(o[0],o[1],o[2]);
235  uVect.SetXYZ(u[0],u[1],u[2]);
236  vVect.SetXYZ(v[0],v[1],v[2]);
237 }
int fVerbose
Definition: poormantracks.C:24
TGeoHMatrix * fCurrentTransMat
__m128 v
Definition: P4_F32vec4.h:4
Int_t GetSensorId() const
Definition: PndGemMCPoint.h:90
Double_t
PndGemMCPoint * fCurrentPndGemMCPoint
PndGemSmearingTask::ClassDef ( PndGemSmearingTask  ,
 
)
private
void PndGemSmearingTask::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 107 of file PndGemSmearingTask.cxx.

References CalcGFDetPlane(), fCurrentPndGemMCPoint, fHitOutputArray, fPointArray, fSigmaX, fSigmaY, fSigmaZ, fVerbose, InitTransMat(), pos, smearLocal(), and v.

108 {
109  // Fills PndGemHits with the MC Truth
110  // TODO filling of RecoHits, together with the sensor plane
111 
112 
113  if ( ! fHitOutputArray ) Fatal("Exec", "No fHitOutputArray");
114  fHitOutputArray->Clear();
115 
116  std::map<Int_t, PndGemHit*> clusterMap;
117 
118  Int_t nPndGemHits=fPointArray->GetEntriesFast();
119  for(Int_t iMvdPoint=0;iMvdPoint<nPndGemHits;++iMvdPoint)
120  {
122  Int_t trackid=fCurrentPndGemMCPoint->GetTrackID();
123  Int_t size = fHitOutputArray->GetEntriesFast();
124  InitTransMat();
125 
126  // cut on secondaries (deltas) etc
127  if(trackid<0)continue;
128 
129  //set the plane definition inside the local frame
130  //sensor origin in the middle, u^ and v^ are xy plane
131  TVector3 o(0.,0.,0.),u(1.,0.,0.),v(0.,1.,0.);
132  CalcGFDetPlane(o,u,v);
133  TVector3 pos;
134  fCurrentPndGemMCPoint->Position(pos);
135  smearLocal(pos);
136  TVector3 dposLocal(fSigmaX,fSigmaY,fSigmaZ);
137 
138  // TODO here we shall distinguish between strip and pixel sensors
139  // TODO How to handle the covariance matrix? OR do we really use local point
140  // errors. this would avoid two conversations, myabe overload the FairHit
141  // functions for the global error points.
142 
143  // Now the 3D Info is smared inside the FairHit part of PndGemHit
144  new ((*fHitOutputArray)[size]) PndGemHit(fCurrentPndGemMCPoint->GetDetectorID(),
145  pos,dposLocal,iMvdPoint,fCurrentPndGemMCPoint->GetEnergyLoss(),1);
146 
147  }//end for PndMvdiMvdPoint
148 
149  if (fVerbose > 0) {
150  std::cout<<fHitOutputArray->GetEntriesFast() <<" Hits created out of "
151  <<fPointArray->GetEntriesFast() <<" Points"<<std::endl;
152  }
153 
154 }
TVector3 pos
int fVerbose
Definition: poormantracks.C:24
TClonesArray * fHitOutputArray
__m128 v
Definition: P4_F32vec4.h:4
TClonesArray * fPointArray
void CalcGFDetPlane(TVector3 &o, TVector3 &u, TVector3 &v)
PndGemMCPoint * fCurrentPndGemMCPoint
void smearLocal(TVector3 &pos)
InitStatus PndGemSmearingTask::Init ( )
virtual

Definition at line 62 of file PndGemSmearingTask.cxx.

References fBranchName, fHitOutputArray, fMctruthArray, fPointArray, and gGeoManager.

Referenced by runGemSmearing().

63 {
64  // Get RootManager
65  FairRootManager* ioman = FairRootManager::Instance();
66  if ( ! ioman ) {
67  std::cout << "-E- PndGemSmearingTask::Init: "
68  << "RootManager not instantiated!" << std::endl;
69  return kFATAL; }
70 
71  // Get input array
72  fPointArray = (TClonesArray*) ioman->GetObject(fBranchName);
73  if ( ! fPointArray ) {
74  std::cout << "-W- PndGemSmearingTask::Init: "<< "No "<<fBranchName
75  <<" array!" << std::endl;
76  return kERROR; }
77 
78  // Get MCTruth collection
79  fMctruthArray=(TClonesArray*) ioman->GetObject("MCTrack");
80  if(fMctruthArray==0) {
81  std::cout << "-W- PndGemSmearingTask::Init: No McTruth array!" << std::endl;
82  return kERROR; }
83 
84  // Create and register output array
85  fHitOutputArray = new TClonesArray("PndGemHit");
86  ioman->Register("GEMHit", "PndMvd ideal Hits",
87  fHitOutputArray, kTRUE);
88 
89  std::cout << "-I- gGeoManager = "<<gGeoManager << std::endl;
90 
91 
92  return kSUCCESS;
93 }
TClonesArray * fHitOutputArray
TGeoManager * gGeoManager
TClonesArray * fPointArray
TClonesArray * fMctruthArray
void PndGemSmearingTask::InitTransMat ( )
private

smearing and geometry access

Definition at line 156 of file PndGemSmearingTask.cxx.

References fCurrentPndGemMCPoint, fCurrentTransMat, fDigiPar, fVerbose, PndGemDigiPar::GetNodeName(), PndGemMCPoint::GetSensorId(), gGeoManager, and TString.

Referenced by Exec().

157 {
158 // std::cout<<"InitTransMat() with "<<fCurrentPndGemMCPoint->GetDetName()<<std::endl;
159  Int_t sensorId = fCurrentPndGemMCPoint->GetSensorId();
160  TString nodeName = fDigiPar->GetNodeName(sensorId);
161  gGeoManager->cd(nodeName.Data());
162  fCurrentTransMat = gGeoManager->GetCurrentMatrix();
163  if (fVerbose > 1) {
164  fCurrentTransMat->Print("");
165  }
166 }
int fVerbose
Definition: poormantracks.C:24
TGeoHMatrix * fCurrentTransMat
TGeoManager * gGeoManager
Int_t GetSensorId() const
Definition: PndGemMCPoint.h:90
TString GetNodeName(Int_t sensorId)
PndGemMCPoint * fCurrentPndGemMCPoint
PndGemDigiPar * fDigiPar
void PndGemSmearingTask::ProduceHits ( )
private
void PndGemSmearingTask::Register ( )
private
void PndGemSmearingTask::Reset ( )
private
void PndGemSmearingTask::SetParContainers ( )
virtual

Virtual method Init

Definition at line 95 of file PndGemSmearingTask.cxx.

References fDigiPar, and rtdb.

96 {
97  // Get Base Container
98  FairRun* ana = FairRun::Instance();
99  FairRuntimeDb* rtdb=ana->GetRuntimeDb();
100 
101  // Get GEM digitisation parameter container
102  fDigiPar = (PndGemDigiPar*)(rtdb->getContainer("PndGemDetectors"));
103 }
Digitization Parameter Class for GEM part.
Definition: PndGemDigiPar.h:31
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
PndGemDigiPar * fDigiPar
void PndGemSmearingTask::smear ( TVector3 &  pos)
private

smear a 3d vector

Definition at line 169 of file PndGemSmearingTask.cxx.

References Double_t, fSigmaX, fSigmaY, fSigmaZ, fVerbose, x, y, and z.

Referenced by smearLocal().

170 {
172 
173  Double_t sigx=gRandom->Gaus(0,fSigmaX);
174  Double_t sigy=gRandom->Gaus(0,fSigmaY);
175  Double_t sigz=gRandom->Gaus(0,fSigmaZ);
176 
177  Double_t x = pos.x() + sigx;
178  Double_t y = pos.y() + sigy;
179  Double_t z = pos.z() + sigz;
180 
181  if (fVerbose > 1) {
182  std::cout<<"PndGemSmearingTask::smear Point (x,y,z)=("
183  <<pos.x()<<","<<pos.z()<<","<<pos.z()<<") by ("
184  <<fSigmaX<<","<<fSigmaY<<","<<fSigmaZ<<") to ";
185  }
186  pos.SetXYZ(x,y,z);
187  if (fVerbose > 1) {
188  std::cout<<"("<<pos.x()<<","<<pos.y()<<","<<pos.z()<<")"<<std::endl;
189  }
190  return;
191 }
TVector3 pos
int fVerbose
Definition: poormantracks.C:24
Double_t
Double_t z
Double_t x
Double_t y
void PndGemSmearingTask::smearLocal ( TVector3 &  pos)
private

smear a 3d vector in the local sensor plane

Definition at line 194 of file PndGemSmearingTask.cxx.

References Double_t, fCurrentTransMat, fVerbose, and smear().

Referenced by Exec().

195 {
197  if (fVerbose > 1) {
198  std::cout<<"PndGemSmearingTask::smearLocal"<<std::endl;
199  }
200  Double_t posLab[3], posSens[3];
201 
202  posLab[0]=pos.x(); posLab[1]=pos.y(); posLab[2]=pos.z();
203  fCurrentTransMat->MasterToLocal(posLab,posSens);
204 
205  pos.SetXYZ(posSens[0],posSens[1],posSens[2]);
206 
207  smear(pos); // apply a gaussian
208 
209  posSens[0]=pos.x(); posSens[1]=pos.y(); posSens[2]=pos.z();
210  fCurrentTransMat->LocalToMaster(posSens,posLab);
211  pos.SetXYZ(posLab[0],posLab[1],posLab[2]);
212 
213 
214  return;
215 }
TVector3 pos
int fVerbose
Definition: poormantracks.C:24
TGeoHMatrix * fCurrentTransMat
void smear(TVector3 &pos)
Double_t

Member Data Documentation

TString PndGemSmearingTask::fBranchName
private

Definition at line 60 of file PndGemSmearingTask.h.

Referenced by Init(), and PndGemSmearingTask().

PndGemMCPoint* PndGemSmearingTask::fCurrentPndGemMCPoint
private

Definition at line 72 of file PndGemSmearingTask.h.

Referenced by CalcGFDetPlane(), Exec(), and InitTransMat().

TGeoHMatrix* PndGemSmearingTask::fCurrentTransMat
private

Definition at line 73 of file PndGemSmearingTask.h.

Referenced by CalcGFDetPlane(), InitTransMat(), and smearLocal().

PndGemDigiPar* PndGemSmearingTask::fDigiPar
private
TMatrixT<Double_t> PndGemSmearingTask::fHitCovMatrix
private

Definition at line 74 of file PndGemSmearingTask.h.

TClonesArray* PndGemSmearingTask::fHitOutputArray
private

Output array of Hits

Definition at line 67 of file PndGemSmearingTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndGemSmearingTask::fMctruthArray
private

Definition at line 64 of file PndGemSmearingTask.h.

Referenced by Init().

TClonesArray* PndGemSmearingTask::fPointArray
private

Input array of Points

Definition at line 63 of file PndGemSmearingTask.h.

Referenced by Exec(), and Init().

Double_t PndGemSmearingTask::fSigmaX
private

Properties

Definition at line 69 of file PndGemSmearingTask.h.

Referenced by Exec(), PndGemSmearingTask(), and smear().

Double_t PndGemSmearingTask::fSigmaY
private

Definition at line 70 of file PndGemSmearingTask.h.

Referenced by Exec(), PndGemSmearingTask(), and smear().

Double_t PndGemSmearingTask::fSigmaZ
private

Definition at line 71 of file PndGemSmearingTask.h.

Referenced by Exec(), PndGemSmearingTask(), and smear().


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