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

#include <PndLmdGeaneTask.h>

Inheritance diagram for PndLmdGeaneTask:

Public Member Functions

 PndLmdGeaneTask ()
 
 PndLmdGeaneTask (Double_t pBeam, TVector3 IP, bool is_prefiltered=true)
 Set up beam momuntum value and position of interaction point. More...
 
 ~PndLmdGeaneTask ()
 
virtual void SetParContainers ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void Finish ()
 
FairTrackParP * PropToPlane (FairTrackParP *fStartPst, double zpos, int dir, bool &isProp)
 

Private Member Functions

void Register ()
 
void Reset ()
 
void ProduceHits ()
 
 ClassDef (PndLmdGeaneTask, 2)
 

Private Attributes

FairGeanePro * fPro
 
PndGeoHandlingfGeoH
 
int fEventNr
 
bool fUseMVDPoint
 
Double_t fPbeam
 
Int_t fPDGid
 
TVector3 vtx
 
FairField * pndField
 
std::map< int, std::vector< int > > fTrackPixHitIdMap
 
std::map< int, std::vector< int > > fTrackStripHitIdMap
 
std::string track_branch_name
 
TClonesArray * fTracks
 
TClonesArray * fTrackParGeane
 
TClonesArray * fTrackParIni
 
TClonesArray * fTrackParFinal
 
TClonesArray * fDetName
 

Detailed Description

Definition at line 35 of file PndLmdGeaneTask.h.

Constructor & Destructor Documentation

PndLmdGeaneTask::PndLmdGeaneTask ( )

Default constructor

Definition at line 41 of file PndLmdGeaneTask.cxx.

References track_branch_name.

42  : FairTask("Geane Task for PANDA Lmd"), fEventNr(0), fUseMVDPoint(false) {
43  track_branch_name = "LMDPndTrackFilt";
44 }
std::string track_branch_name
PndLmdGeaneTask::PndLmdGeaneTask ( Double_t  pBeam,
TVector3  IP,
bool  is_prefiltered = true 
)

Set up beam momuntum value and position of interaction point.

Definition at line 47 of file PndLmdGeaneTask.cxx.

References fPbeam, fPDGid, track_branch_name, and vtx.

49  : FairTask("Geane Task for PANDA Lmd"), fEventNr(0), fUseMVDPoint(false) {
50  fPDGid = -2212;
51  fPbeam = pBeam;
52  cout << "Beam Momentum for particle with PDGid#" << fPDGid << " this run is "
53  << fPbeam << endl;
54  vtx = IP;
55  cout << "Interaction Point:" << endl;
56  vtx.Print();
57 
58  if (is_prefiltered)
59  track_branch_name = "LMDPndTrackFilt";
60  else
61  track_branch_name = "LMDPndTrack";
62 }
std::string track_branch_name
PndLmdGeaneTask::~PndLmdGeaneTask ( )

Destructor

Definition at line 65 of file PndLmdGeaneTask.cxx.

65 {}

Member Function Documentation

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

Virtual method Exec

Propagate to the PCA (a space point) in 7 steps


Definition at line 112 of file PndLmdGeaneTask.cxx.

References Bool_t, fDetName, fEventNr, fPbeam, fPro, fTrackParFinal, fTrackParGeane, fTrackParIni, fTracks, fVerbose, PndTrack::GetParamFirst(), i, PropToPlane(), sqrt(), and vtx.

112  {
113  fTrackParGeane->Delete();
114  fTrackParIni->Delete();
115  fTrackParFinal->Delete();
116  fDetName->Delete();
117 
118  if (fVerbose > 2) {
119  cout << " ---- Info: " << fEventNr << endl;
120  }
121  fEventNr++;
122  // Charge & mass of particle
123  Int_t PDGCode = -2212; // antiproton
124 
125  // go through all tracks
126  int glI = fTracks->GetEntriesFast();
127  // cout<<"glI = "<<glI<<endl;
128  Int_t counterGeaneTrk = 0;
129 
130  for (Int_t i = 0; i < glI; i++) {
131  TVector3 StartPos, StartPosErr, StartMom, StartMomErr, StartO, StartU,
132  StartV;
133 
134  PndTrack* recTrack = (PndTrack*)(fTracks->At(i));
135  FairTrackParP fFittedTrkP = recTrack->GetParamFirst();
136 
137  TVector3 PosRecLMD(fFittedTrkP.GetX(), fFittedTrkP.GetY(),
138  fFittedTrkP.GetZ());
139  if (fFittedTrkP.GetZ() > 1130 || fFittedTrkP.GetZ() < 1000 )
140  continue; // TEST: skip trks from 2nd plane. TODO: check are they fine???
141  TVector3 MomRecLMD(fFittedTrkP.GetPx(), fFittedTrkP.GetPy(),
142  fFittedTrkP.GetPz());
143  MomRecLMD *=
144  fPbeam / MomRecLMD.Mag(); // external assumption about mom magnitude
145  fFittedTrkP.SetPx(MomRecLMD.X());
146  fFittedTrkP.SetPy(MomRecLMD.Y());
147  fFittedTrkP.SetPz(MomRecLMD.Z());
148  double covMARS[6][6];
149  fFittedTrkP.GetMARSCov(covMARS);
150  TVector3 errMomRecLMD(sqrt(covMARS[0][0]), sqrt(covMARS[1][1]),
151  sqrt(covMARS[2][2]));
152  TVector3 errPosRecLMD(sqrt(covMARS[3][3]), sqrt(covMARS[4][4]),
153  sqrt(covMARS[5][5]));
154 
155  StartPos = PosRecLMD;
156  StartPosErr = errPosRecLMD;
157  StartMom = MomRecLMD;
158  StartMomErr = errMomRecLMD;
159 
162  // Comment: seems back propagation in one step (for 11 m) is too much
163  // try smoothing it by small steps, which follow mag.field
164  const int nstep = 7;
165  double zbend[nstep] = {661, 660.5, 660., 659,
166  319, 316, 220}; // entarance and exit mag.field
167  if (fVerbose > 2) {
168  cout << "------------------------------------------" << endl;
169  cout << "StartPos:" << endl;
170  StartPos.Print();
171  cout << "StartPosErr:" << endl;
172  StartPosErr.Print();
173 
174  cout << "StartMom: " << StartMom.Mag() << endl;
175  StartMom.Print();
176  cout << "StartMomErr: " << StartMomErr.Mag() << endl;
177  StartMomErr.Print();
178  cout << "" << endl;
179  }
180  TClonesArray& clref1 = *fTrackParIni;
181  Int_t size1 = clref1.GetEntriesFast();
182 
183  FairTrackParP* fStartPst = new FairTrackParP(fFittedTrkP); // REC
184  for (int js = 1; js < nstep; js++) {
185  bool isProp;
186 
187  FairTrackParP* fResPst =
188  PropToPlane(fStartPst, zbend[js], -1, isProp); // back propagation
189 
190  if (isProp) {
191  delete fStartPst;
192  fStartPst = fResPst;
193  } else
194  break;
195  }
196 
197  // ///and now Propagate to the PCA (a space point) in one step
198  // ---------------------------------
199  int ierr = 0;
200  FairTrackParH* fStart = new (clref1[size1]) FairTrackParH(fStartPst, ierr);
201  TClonesArray& clref = *fTrackParGeane;
202  Int_t size = clref.GetEntriesFast();
204  FairTrackParH* fRes = new (clref[size]) FairTrackParH();
205  fPro->SetPoint(vtx);
206  fPro->PropagateToPCA(1, -1); // back-propagate to point
207  Bool_t isProp = fPro->Propagate(fStart, fRes, PDGCode);
209 
210  // // // /////back propagation to line/plane TEST----
211  // Bool_t isProp;
212  // FairTrackParP* fResPst = PropToXZPlane(fStartPst,0,-1,isProp);//back
213  // propagation to plane TEST
214  // FairTrackParH *fRes = new(clref[size]) FairTrackParH(fResPst,ierr);
215 
216  // // FairTrackParH* fRes = PropToLine(fStart,0,-1,isProp);//back
217  // propagation to line TEST
218 
219  // // // // /////END back propagation to line/plane TEST
220  // ///----------------------------------------------------------------------
221  delete fStartPst;
222 
224 
225  //------- TEST of calculation errors -------
226  TVector3 gPos(fRes->GetX(), fRes->GetY(), fRes->GetZ());
227  TVector3 gMom(fRes->GetPx(), fRes->GetPy(), fRes->GetPz());
228 
229  TVector3 gErrPos(fRes->GetDX(), fRes->GetDY(), fRes->GetDZ());
230  TVector3 gErrMom(fRes->GetDPx(), fRes->GetDPy(), fRes->GetDPz());
231  // cout<<" "<<endl;
232  if (fVerbose > 2) {
233  // cout<<"================= %%%% ===================="<<endl;
234 
235  cout << "gPos:" << endl;
236  gPos.Print();
237  // cout<<"difference between final and initial point =
238  // "<<(-StartPos.Mag()+gPos.Mag())<<endl;
239  cout << "gErrPos:" << endl;
240  gErrPos.Print();
241  // cout<<"difference between final and initial point error =
242  // "<<(-StartPosErr.Mag()+gErrPos.Mag())<<endl;
243 
244  cout << "gMom: " << gMom.Mag() << endl;
245  gMom.Print();
246  // cout<<"difference between initial and final momentum =
247  // "<<(StartMom.Mag()-gMom.Mag())<<endl;
248  cout << "gErrMom: " << gErrMom.Mag() << endl;
249  // cout<<"difference between initial and final momentum error=
250  // "<<(StartMomErr.Mag()-gErrMom.Mag())<<endl;
251  gErrMom.Print();
252 
253  cout << "================= %%%% ====================" << endl;
254  }
255  //------------------------------------------
256 
257  if (isProp == kTRUE) {
258  new ((*fTrackParFinal)[counterGeaneTrk])
259  FairTrackParH(*(fRes)); // save Track
260  // new((*fTrackParFinal)[counterGeaneTrk]) FairTrackParP(*(fRes));
262  counterGeaneTrk++;
263  if (fVerbose > 2) cout << "***** isProp TRUE *****" << endl;
264  } else {
265  if (fVerbose > 2) {
266  cout << "!!! Back-propagation with GEANE didn't return result !!!"
267  << endl;
268  cout << "StartPos:" << endl;
269  StartPos.Print();
270  cout << "StartPosErr:" << endl;
271  StartPosErr.Print();
272 
273  cout << "StartMom: " << StartMom.Mag() << endl;
274  StartMom.Print();
275  cout << "StartMomErr: " << StartMomErr.Mag() << endl;
276  StartMomErr.Print();
277  }
278  new ((*fTrackParFinal)[counterGeaneTrk]) FairTrackParH(); // save NULL
279  // new((*fTrackParFinal)[counterGeaneTrk]) FairTrackParP(); //save
280  //NULL
281  counterGeaneTrk++;
282  }
283  }
284  // fMCTracks->Delete();
285  // fMCHits->Delete();
286  if (fVerbose > 2) cout << "PndLmdGeaneTask::Exec END!" << endl;
287 }
int fVerbose
Definition: poormantracks.C:24
FairGeanePro * fPro
Int_t i
Definition: run_full.C:25
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
TClonesArray * fTrackParGeane
TClonesArray * fTracks
TClonesArray * fDetName
TClonesArray * fTrackParIni
TClonesArray * fTrackParFinal
FairTrackParP * PropToPlane(FairTrackParP *fStartPst, double zpos, int dir, bool &isProp)
FairTrackParP GetParamFirst()
Definition: PndTrack.h:49
void PndLmdGeaneTask::Finish ( )
virtual

Definition at line 289 of file PndLmdGeaneTask.cxx.

289  {
290  // tprop->Write();
291 }
InitStatus PndLmdGeaneTask::Init ( )
virtual

Definition at line 68 of file PndLmdGeaneTask.cxx.

References fDetName, fGeoH, fPro, fTrackParFinal, fTrackParGeane, fTrackParIni, fTracks, PndGeoHandling::Instance(), pndField, and track_branch_name.

Referenced by reco_LMD(), runLumi5Geane(), and runLumiPixel5BackProp().

68  {
69  // Get RootManager
70  FairRootManager* ioman = FairRootManager::Instance();
71  if (!ioman) {
72  std::cout << "-E- PndLmdGeaneTask::Init: "
73  << "RootManager not instantiated!" << std::endl;
74  return kFATAL;
75  }
76 
77  // fTracks = (TClonesArray*) ioman->GetObject("LMDTrack");
78  fTracks = (TClonesArray*)ioman->GetObject(track_branch_name.c_str());
79  if (!fTracks) {
80  std::cout << "-W- PndLmdGeaneTask::Init: "
81  << "No Track"
82  << " array!" << std::endl;
83  return kERROR;
84  }
85 
86  fTrackParGeane = new TClonesArray("FairTrackParH");
87  ioman->Register("GeaneTrackPar", "Geane", fTrackParGeane, kTRUE);
88 
89  fTrackParIni = new TClonesArray("FairTrackParH");
90  ioman->Register("GeaneTrackIni", "Geane", fTrackParIni, kTRUE);
91 
92  fTrackParFinal = new TClonesArray("FairTrackParH");
93  ioman->Register("GeaneTrackFinal", "Geane", fTrackParFinal, kTRUE);
94 
95  fDetName = new TClonesArray("TObjString");
96  ioman->Register("DetName", "Geane", fDetName, kTRUE);
97 
98  fPro = new FairGeanePro();
100  // FairRun* fRun = FairRun::Instance(); //[R.K. 01/2017] unused variable?
101  // FairRuntimeDb* rtdb = fRun->GetRuntimeDb(); //[R.K. 01/2017] unused
102  // variable
103 
104  pndField = FairRunAna::Instance()->GetField();
105 
106  return kSUCCESS;
107 }
PndGeoHandling * fGeoH
FairGeanePro * fPro
static PndGeoHandling * Instance()
TClonesArray * fTrackParGeane
TClonesArray * fTracks
TClonesArray * fDetName
TClonesArray * fTrackParIni
TClonesArray * fTrackParFinal
FairField * pndField
std::string track_branch_name
void PndLmdGeaneTask::ProduceHits ( )
private
FairTrackParP * PndLmdGeaneTask::PropToPlane ( FairTrackParP *  fStartPst,
double  zpos,
int  dir,
bool &  isProp 
)

Definition at line 293 of file PndLmdGeaneTask.cxx.

References fPbeam, fPDGid, and fPro.

Referenced by Exec().

295  {
296  TVector3 stStartPos(fStartPst->GetX(), fStartPst->GetY(), fStartPst->GetZ());
297  if (zpos > 1e3) { // external assumption about mom magnitude [use it while in
298  // BOX and we know mom should be diff]
299  TVector3 MomStartPos(fStartPst->GetPx(), fStartPst->GetPy(),
300  fStartPst->GetPz());
301  MomStartPos *= fPbeam / MomStartPos.Mag();
302  fStartPst->SetPx(MomStartPos.X());
303  fStartPst->SetPy(MomStartPos.Y());
304  fStartPst->SetPz(MomStartPos.Z()); // correct mom.magnitude
305  }
306  // propagate plane-to-plane
307  TVector3 ist(fStartPst->GetIVer());
308  TVector3 jst(fStartPst->GetJVer());
309  TVector3 kst(fStartPst->GetKVer());
310  TVector3 oc(0, 0, zpos);
311  TVector3 dj(1., 0., 0.);
312  TVector3 dk(0., 1., 0.);
313  dj.SetMag(1);
314  dk.SetMag(1);
315  fPro->PropagateFromPlane(jst, kst); // 1st detector plane
316  fPro->PropagateToPlane(oc, dj, dk); // virtual plane at fixed z
317  if (dir < 0) fPro->setBackProp();
318  FairTrackParP* fResPst = new FairTrackParP();
319  isProp = fPro->Propagate(fStartPst, fResPst, fPDGid);
320  return fResPst;
321 }
FairGeanePro * fPro
const Double_t zpos
void PndLmdGeaneTask::Register ( )
private
void PndLmdGeaneTask::Reset ( )
private
void PndLmdGeaneTask::SetParContainers ( )
virtual

Virtual method Init

Definition at line 109 of file PndLmdGeaneTask.cxx.

109 {}

Member Data Documentation

TClonesArray* PndLmdGeaneTask::fDetName
private

Definition at line 78 of file PndLmdGeaneTask.h.

Referenced by Exec(), and Init().

int PndLmdGeaneTask::fEventNr
private

Definition at line 60 of file PndLmdGeaneTask.h.

Referenced by Exec().

PndGeoHandling* PndLmdGeaneTask::fGeoH
private

Definition at line 59 of file PndLmdGeaneTask.h.

Referenced by Init().

Double_t PndLmdGeaneTask::fPbeam
private

Definition at line 62 of file PndLmdGeaneTask.h.

Referenced by Exec(), PndLmdGeaneTask(), and PropToPlane().

Int_t PndLmdGeaneTask::fPDGid
private

Definition at line 63 of file PndLmdGeaneTask.h.

Referenced by PndLmdGeaneTask(), and PropToPlane().

FairGeanePro* PndLmdGeaneTask::fPro
private

Definition at line 58 of file PndLmdGeaneTask.h.

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

TClonesArray* PndLmdGeaneTask::fTrackParFinal
private

Definition at line 77 of file PndLmdGeaneTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndLmdGeaneTask::fTrackParGeane
private

Definition at line 75 of file PndLmdGeaneTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndLmdGeaneTask::fTrackParIni
private

Definition at line 76 of file PndLmdGeaneTask.h.

Referenced by Exec(), and Init().

std::map<int, std::vector<int> > PndLmdGeaneTask::fTrackPixHitIdMap
private

Definition at line 66 of file PndLmdGeaneTask.h.

TClonesArray* PndLmdGeaneTask::fTracks
private

Definition at line 72 of file PndLmdGeaneTask.h.

Referenced by Exec(), and Init().

std::map<int, std::vector<int> > PndLmdGeaneTask::fTrackStripHitIdMap
private

Definition at line 67 of file PndLmdGeaneTask.h.

bool PndLmdGeaneTask::fUseMVDPoint
private

Definition at line 61 of file PndLmdGeaneTask.h.

FairField* PndLmdGeaneTask::pndField
private

Definition at line 65 of file PndLmdGeaneTask.h.

Referenced by Init().

std::string PndLmdGeaneTask::track_branch_name
private

Definition at line 69 of file PndLmdGeaneTask.h.

Referenced by Init(), and PndLmdGeaneTask().

TVector3 PndLmdGeaneTask::vtx
private

Definition at line 64 of file PndLmdGeaneTask.h.

Referenced by Exec(), and PndLmdGeaneTask().


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