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

#include <CbmGeaneTrEmc.h>

Inheritance diagram for FairGeaneTrEmc:

Public Member Functions

 FairGeaneTrEmc ()
 
 ~FairGeaneTrEmc ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 

Private Member Functions

 ClassDef (FairGeaneTrEmc, 1)
 

Private Attributes

TClonesArray * fPointArray1
 
TClonesArray * fPointArray2
 
PndMCTrackfPoint1
 
PndEmcPointfPoint2
 
TTree * t
 
TFile * f
 
TClonesArray * fTrackParIni
 
TClonesArray * fTrackParGeane
 
TClonesArray * fTrackParFinal
 
TGeant3 * gMC3
 
Int_t fEvent
 
FairGeanePro * fPro
 

Detailed Description

Definition at line 19 of file CbmGeaneTrEmc.h.

Constructor & Destructor Documentation

FairGeaneTrEmc::FairGeaneTrEmc ( )

Default constructor

Definition at line 20 of file CbmGeaneTrEmc.cxx.

20  :
21  FairTask("Test") { }
FairGeaneTrEmc::~FairGeaneTrEmc ( )

Destructor

Definition at line 27 of file CbmGeaneTrEmc.cxx.

27 { }

Member Function Documentation

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

Virtual method Exec

Definition at line 85 of file CbmGeaneTrEmc.cxx.

References Bool_t, CAMath::Cos(), Double_t, fParticle, fPhi, fPoint1, fPoint2, fPointArray1, fPointArray2, fPro, fTrackParFinal, fTrackParGeane, fTrackParIni, fX, fY, fZ, PndEmcPoint::GetModule(), PndMCTrack::GetMomentum(), PndMCTrack::GetPdgCode(), PndMCTrack::GetStartVertex(), and CAMath::Sin().

85  {
86 // cout << "FairGeaneTrEmc::Exec" << endl;
87  fTrackParGeane->Delete();
88  fTrackParIni->Delete();
89  fTrackParFinal->Delete();
90 
91  cout<<fPointArray1->GetEntriesFast()<<"/"<<fPointArray2->GetEntriesFast()<<endl;
92 
93  fPoint1 = (PndMCTrack *)fPointArray1->At(0);
94  Int_t trId=0;
95  fPoint2=0;
96  for (Int_t k=0; k<fPointArray2->GetEntriesFast(); k++) {
97  fPoint2 = (PndEmcPoint *)fPointArray2->At(k);
98  // cout << "loop for " << k << " / " << fPoint2->GetTrackID() << endl;
99  //if(fPoint2->GetTrackID()==trId && (fPoint2->GetModule()==1 || fPoint2->GetModule()==2) ) break;
100  if(fPoint2->GetTrackID()==trId) break;
101  fPoint2=0;
102  }
103 
104  if(fPoint2==0) return;
105 
106  if (fPoint2->GetModule()<3)
107  {
108  fPro->PropagateToVolume("Emc12",0,1);
109  }
110  else if (fPoint2->GetModule()==3)
111  {
112  fPro->PropagateToVolume("Emc3",0,1);
113  }
114  else if (fPoint2->GetModule()==4)
115  {
116  fPro->PropagateToVolume("Emc4",0,1);
117  }
118  else
119  {
120  return;
121  }
122 
123  TVector3 StartPos = fPoint1->GetStartVertex();
124 // TVector3 StartPos = TVector3 (0.,0.,0.);
125  TVector3 StartPosErr = TVector3(0,0,0);
126  TVector3 StartMom = fPoint1->GetMomentum();
127 // TVector3 StartMom = TVector3 (1.,0.1,0.1);
128  TVector3 StartMomErr = TVector3(0,0,0);
129 
130  TVector3 EndPos =TVector3 (fPoint2->GetX(),fPoint2->GetY(),fPoint2->GetZ());
131  TVector3 EndPosErr=TVector3(0,0,0);
132  TVector3 EndMom= TVector3 (fPoint2->GetPx(),fPoint2->GetPy(),fPoint2->GetPz());
133  TVector3 EndMomErr=TVector3(0,0,0);
134 
135  Int_t PDGCode=fPoint1->GetPdgCode();
136  TDatabasePDG *fdbPDG= TDatabasePDG::Instance();
137  TParticlePDG *fParticle= fdbPDG->GetParticle(PDGCode);
138  Double_t fCharge= fParticle->Charge();
139 
140  TClonesArray& clref1 = *fTrackParIni;
141  Int_t size1 = clref1.GetEntriesFast();
142  FairTrackParH *fStart= new (clref1[size1]) FairTrackParH(StartPos, StartMom, StartPosErr, StartMomErr, fCharge);
143 
144  TClonesArray& clref = *fTrackParGeane;
145  Int_t size = clref.GetEntriesFast();
146  FairTrackParH *fRes= new(clref[size]) FairTrackParH();
147 
148  TClonesArray& clref2 = *fTrackParFinal;
149  Int_t size2 = clref2.GetEntriesFast();
150  FairTrackParH *fFinal= new(clref2[size2]) FairTrackParH(EndPos, EndMom, EndPosErr, EndMomErr, fCharge);
151 
152  cout << "fStart = " << endl;
153  fStart->Print();
154  Bool_t rc = fPro->Propagate(fStart, fRes,PDGCode);
155 
156  cout << "fFinal = " << endl;
157  fFinal->Print();
158  cout << "fEnd = " << endl;
159  fRes->Print();
160 
161  Double_t fLm,fPhi,cLm,sLm,cphi,sphi,fX_sc,fY_sc,fZ_sc,fX,fY,fZ;
162 
163  fLm = fRes->GetLambda();
164  fPhi= fRes->GetPhi();
165  fX = fFinal->GetX();
166  fY = fFinal->GetY();
167  fZ = fFinal->GetZ();
168  cLm= TMath::Cos(fLm);
169  sLm= TMath::Sin(fLm);
170  cphi= TMath::Cos(fPhi);
171  sphi= TMath::Sin(fPhi);
172 
173  fX_sc = fX*cphi*cLm+ fY*cLm*sphi+fZ*sLm;
174  fY_sc = fY*cphi-fX*sphi;
175  fZ_sc = fZ*cLm-fY*sLm*sphi-fX*sLm*cphi;
176 
177  fFinal->SetX_sc(fX_sc);
178  fFinal->SetY_sc(fY_sc);
179  fFinal->SetZ_sc(fZ_sc);
180 
181  // }
182  // }
183  // if (fPointArray1) fPointArray1->Delete();
184  //if (fPointArray2) fPointArray2->Delete();
185 
186 }
represents a mc hit in an emc crystal
Definition: PndEmcPoint.h:19
Short_t GetModule() const
Definition: PndEmcPoint.h:60
TClonesArray * fTrackParIni
Definition: CbmGeaneTrEmc.h:44
Int_t GetPdgCode() const
Definition: PndMCTrack.h:73
static T Sin(const T &x)
Definition: PndCAMath.h:42
Double_t fX
Definition: PndCaloDraw.cxx:34
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
static T Cos(const T &x)
Definition: PndCAMath.h:43
Double_t fZ
Definition: PndCaloDraw.cxx:34
Int_t * fParticle
Definition: run_full.C:24
Double_t
FairGeanePro * fPro
Definition: CbmGeaneTrEmc.h:51
Double_t const fPhi
TClonesArray * fPointArray2
Definition: CbmGeaneTrEmc.h:36
TClonesArray * fTrackParFinal
Definition: CbmGeaneTrEmc.h:46
TClonesArray * fPointArray1
Definition: CbmGeaneTrEmc.h:35
Double_t fY
Definition: PndCaloDraw.cxx:34
PndMCTrack * fPoint1
Definition: CbmGeaneTrEmc.h:38
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
PndEmcPoint * fPoint2
Definition: CbmGeaneTrEmc.h:39
TClonesArray * fTrackParGeane
Definition: CbmGeaneTrEmc.h:45
InitStatus FairGeaneTrEmc::Init ( )
virtual

Virtual method Init

Definition at line 32 of file CbmGeaneTrEmc.cxx.

References fPointArray1, fPointArray2, fPro, fTrackParFinal, fTrackParGeane, and fTrackParIni.

32  {
33 
34  // Get RootManager
35  FairRootManager* ioman = FairRootManager::Instance();
36  if ( ! ioman ) {
37  cout << "-E- FairGeaneTrEmc::Init: "
38  << "RootManager not instantised!" << endl;
39  return kFATAL;
40  }
41 
42  // Get input array
43  fPointArray1 = (TClonesArray*) ioman->GetObject("MCTrack");
44  if ( ! fPointArray1 ) {
45  cout << "-W- FairGeaneTrEmc::Init: "
46  << "No MCTrack array!" << endl;
47  return kERROR;
48  }
49 
50 
51  fPointArray2 = (TClonesArray*) ioman->GetObject("EmcPoint");
52  if ( ! fPointArray2 ) {
53  cout << "-W- FairGeaneTrEmc::Init: "
54  << "No EmcPoint array!" << endl;
55  return kERROR;
56  }
57 
58  fTrackParGeane = new TClonesArray("FairTrackParH");
59  ioman->Register("GeaneTrackPar","Geane", fTrackParGeane, kTRUE);
60 
61  fTrackParIni = new TClonesArray("FairTrackParH");
62  ioman->Register("GeaneTrackIni","Geane", fTrackParIni, kTRUE);
63 
64  fTrackParFinal = new TClonesArray("FairTrackParH");
65  ioman->Register("GeaneTrackFinal","Geane", fTrackParFinal, kTRUE);
66 
67  // Create and register output array
68 
69  fPro = new FairGeanePro();
70 
71 
72  // if the vis manager is available then initialize it!
73 
74  /* FairTrajFilter *fTrajFilter = FairTrajFilter::Instance();
75  if(fTrajFilter) fTrajFilter->Init("GeaneTrk", "Geane");
76 */
77 
78  return kSUCCESS;
79 
80 }
TClonesArray * fTrackParIni
Definition: CbmGeaneTrEmc.h:44
FairGeanePro * fPro
Definition: CbmGeaneTrEmc.h:51
TClonesArray * fPointArray2
Definition: CbmGeaneTrEmc.h:36
TClonesArray * fTrackParFinal
Definition: CbmGeaneTrEmc.h:46
TClonesArray * fPointArray1
Definition: CbmGeaneTrEmc.h:35
TClonesArray * fTrackParGeane
Definition: CbmGeaneTrEmc.h:45

Member Data Documentation

TFile* FairGeaneTrEmc::f
private

Definition at line 42 of file CbmGeaneTrEmc.h.

Int_t FairGeaneTrEmc::fEvent
private

Definition at line 50 of file CbmGeaneTrEmc.h.

PndMCTrack* FairGeaneTrEmc::fPoint1
private

Definition at line 38 of file CbmGeaneTrEmc.h.

Referenced by Exec().

PndEmcPoint* FairGeaneTrEmc::fPoint2
private

Definition at line 39 of file CbmGeaneTrEmc.h.

Referenced by Exec().

TClonesArray* FairGeaneTrEmc::fPointArray1
private

Input array of Points

Definition at line 35 of file CbmGeaneTrEmc.h.

Referenced by Exec(), and Init().

TClonesArray* FairGeaneTrEmc::fPointArray2
private

Definition at line 36 of file CbmGeaneTrEmc.h.

Referenced by Exec(), and Init().

FairGeanePro* FairGeaneTrEmc::fPro
private

Definition at line 51 of file CbmGeaneTrEmc.h.

Referenced by Exec(), and Init().

TClonesArray* FairGeaneTrEmc::fTrackParFinal
private

Definition at line 46 of file CbmGeaneTrEmc.h.

Referenced by Exec(), and Init().

TClonesArray* FairGeaneTrEmc::fTrackParGeane
private

Definition at line 45 of file CbmGeaneTrEmc.h.

Referenced by Exec(), and Init().

TClonesArray* FairGeaneTrEmc::fTrackParIni
private

Output array of Hits

Definition at line 44 of file CbmGeaneTrEmc.h.

Referenced by Exec(), and Init().

TGeant3* FairGeaneTrEmc::gMC3
private

Definition at line 48 of file CbmGeaneTrEmc.h.

TTree* FairGeaneTrEmc::t
private

Definition at line 41 of file CbmGeaneTrEmc.h.


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