Class Reference for E1039 Core & Analysis Software
SQTruthVertexing.cxx
Go to the documentation of this file.
1 #include "SQTruthVertexing.h"
2 
4 #include <phool/PHNodeIterator.h>
5 #include <phool/PHIODataNode.h>
6 #include <phool/PHRandomSeed.h>
7 #include <phool/getClass.h>
12 
13 #include "SRecEvent.h"
14 #include "GFTrack.h"
15 
16 #include <iostream>
17 #include <vector>
18 
19 SQTruthVertexing::SQTruthVertexing(const std::string& name):
20  SubsysReco(name),
21  legacyContainer(false),
22  vtxSmearing(-1.),
23  recEvent(nullptr),
24  recTrackVec(nullptr),
25  truthTrackVec(nullptr),
26  truthDimuonVec(nullptr),
27  recDimuonVec(nullptr)
28 {}
29 
31 {}
32 
34 {
36 }
37 
39 {
40  int ret = GetNodes(topNode);
41  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
42 
43  ret = MakeNodes(topNode);
44  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
45 
46  if(vtxSmearing) rndm.SetSeed(PHRandomSeed());
47 
49 }
50 
52 {
53  std::map<int, SRecTrack*> posTracks;
54  std::map<int, SRecTrack*> negTracks;
55  int nTracks = truthTrackVec->size();
56  int nRecTracks = legacyContainer ? recEvent->getNTracks() : recTrackVec->size();
57  for(int i = 0; i < nTracks; ++i)
58  {
59  SQTrack* trk = truthTrackVec->at(i);
60  int recTrackIdx = trk->get_rec_track_id();
61  if(recTrackIdx < 0 || recTrackIdx >= nRecTracks) continue;
62 
63  SRecTrack* recTrack = legacyContainer ? &(recEvent->getTrack(recTrackIdx)) : dynamic_cast<SRecTrack*>(recTrackVec->at(recTrackIdx));
64  double z_vtx = trk->get_pos_vtx().Z() + (vtxSmearing>0. ? rndm.Gaus(0., vtxSmearing) : 0.);
65 
66  if(swimTrackToVertex(recTrack, z_vtx) < 0.) continue;
67  if(recTrack->getCharge() > 0)
68  {
69  posTracks[trk->get_track_id()] = recTrack;
70  }
71  else
72  {
73  negTracks[trk->get_track_id()] = recTrack;
74  }
75  }
76  if(posTracks.empty() || negTracks.empty()) return Fun4AllReturnCodes::EVENT_OK;
77 
78  //if true dimuons exist, only pair the tracks that are associated with the true dimuon, otherwise produce all possible combinations
79  if(truthDimuonVec)
80  {
81  int nTrueDimuons = truthDimuonVec->size();
82  for(int i = 0; i < nTrueDimuons; ++i)
83  {
84  SQDimuon* trueDimuon = truthDimuonVec->at(i);
85  int pid = trueDimuon->get_track_id_pos();
86  int mid = trueDimuon->get_track_id_neg();
87  if(posTracks.find(pid) == posTracks.end() || negTracks.find(mid) == negTracks.end()) continue;
88 
89  SRecDimuon recDimuon;
90  double z_vtx = trueDimuon->get_pos().Z() + (vtxSmearing>0. ? rndm.Gaus(0., vtxSmearing) : 0.);
91  if(!buildRecDimuon(z_vtx, posTracks[pid], negTracks[mid], &recDimuon)) continue;
92  recDimuon.trackID_pos = pid;
93  recDimuon.trackID_neg = mid;
94 
95  trueDimuon->set_rec_dimuon_id(legacyContainer ? recEvent->getNDimuons() : recDimuonVec->size());
96  if(!legacyContainer)
97  recDimuonVec->push_back(&recDimuon);
98  else
99  recEvent->insertDimuon(recDimuon);
100  }
101  }
102  // else
103  // {
104  // for(auto ptrk = posTracks.begin(); ptrk != posTracks.end(); ++ptrk)
105  // {
106  // double z_pos = ptrk->second->getVertexPos().Z();
107  // for(auto mtrk = negTracks.begin(); mtrk != negTracks.end(); ++mtrk)
108  // {
109  // double z_neg = mtrk->second->getVertexPos().Z();
110  // if(fabs(z_pos - z_neg) > 100.) continue;
111 
112  // double z_vtx = 0.5*(z_pos + z_neg) + (vtxSmearing>0. ? rndm.Gaus(0., vtxSmearing) : 0.);
113  // SRecDimuon recDimuon;
114  // if(!buildRecDimuon(z_vtx, ptrk->second, mtrk->second, &recDimuon)) continue;
115  // recDimuon.trackID_pos = ptrk->first;
116  // recDimuon.trackID_neg = mtrk->first;
117 
118  // if(!legacyContainer)
119  // recDimuonVec->push_back(&recDimuon);
120  // else
121  // recEvent->insertDimuon(recDimuon);
122  // }
123  // }
124  // }
125 
127 }
128 
130 {
132 }
133 
134 int SQTruthVertexing::GetNodes(PHCompositeNode* topNode)
135 {
136  truthTrackVec = findNode::getClass<SQTrackVector>(topNode, "SQTruthTrackVector");
137  if(!truthTrackVec)
138  {
139  std::cerr << Name() << ": failed finding truth track info, abort." << std::endl;
141  }
142 
143  truthDimuonVec = findNode::getClass<SQDimuonVector>(topNode, "SQTruthDimuonVector");
144  if(!truthDimuonVec)
145  {
146  std::cout << Name() << ": failed finding truth dimuon info, rec dimuon will be for reference only. " << std::endl;
147  }
148 
149  if(legacyContainer)
150  {
151  recEvent = findNode::getClass<SRecEvent>(topNode, "SRecEvent");
152  }
153  else
154  {
155  recTrackVec = findNode::getClass<SQTrackVector>(topNode, "SQRecTrackVector");
156  }
157 
158  if((!recEvent) && (!recTrackVec))
159  {
160  std::cerr << Name() << ": failed finding rec track info, abort." << std::endl;
162  }
163 
165 }
166 
167 int SQTruthVertexing::MakeNodes(PHCompositeNode* topNode)
168 {
169  if(!legacyContainer)
170  {
171  PHNodeIterator iter(topNode);
172  PHCompositeNode* dstNode = static_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
173  if(!dstNode)
174  {
175  std::cerr << Name() << ": cannot locate DST node, abort." << std::endl;
177  }
178 
179  recDimuonVec = new SQDimuonVector_v1();
180  dstNode->addNode(new PHIODataNode<PHObject>(recDimuonVec, "SQRecDimuonVector", "PHObject"));
181  }
183 }
184 
185 double SQTruthVertexing::swimTrackToVertex(SRecTrack* track, double z, TVector3* pos, TVector3* mom)
186 {
187  SQGenFit::GFTrack gftrk(*track);
188 
189  TVector3 p, m;
190  double chi2 = gftrk.swimToVertex(z, &p, &m);
191  if(chi2 < 0.) return chi2;
192 
193  if(pos == nullptr)
194  {
195  track->setChisqVertex(chi2);
196  track->setVertexPos(p);
197  track->setVertexMom(m);
198  }
199  else
200  {
201  pos->SetXYZ(p.X(), p.Y(), p.Z());
202  mom->SetXYZ(m.X(), m.Y(), m.Z());
203  }
204 
205  return chi2;
206 }
207 
208 bool SQTruthVertexing::buildRecDimuon(double z_vtx, SRecTrack* posTrack, SRecTrack* negTrack, SRecDimuon* dimuon)
209 {
210  TVector3 p_mom, p_pos = posTrack->getVertexPos();
211  double p_chi2;
212  if(fabs(p_pos.Z() - z_vtx) > 1.)
213  {
214  p_chi2 = swimTrackToVertex(posTrack, z_vtx, &p_pos, &p_mom);
215  if(p_chi2 < 0.) return false;
216  }
217  else
218  {
219  p_mom = posTrack->getVertexMom();
220  p_chi2 = posTrack->getChisqVertex();
221  }
222 
223  TVector3 m_mom, m_pos = negTrack->getVertexPos();
224  double m_chi2;
225  if(fabs(m_pos.Z() - z_vtx) > 1.)
226  {
227  m_chi2 = swimTrackToVertex(negTrack, z_vtx, &p_pos, &p_mom);
228  if(m_chi2 < 0.) return false;
229  }
230  else
231  {
232  m_mom = negTrack->getVertexMom();
233  m_chi2 = negTrack->getChisqVertex();
234  }
235 
236  dimuon->p_pos.SetVectM(p_mom, M_MU);
237  dimuon->p_neg.SetVectM(m_mom, M_MU);
238  dimuon->chisq_kf = p_chi2 + m_chi2;
239  dimuon->vtx.SetXYZ(0., 0., z_vtx);
240  dimuon->vtx_pos = p_pos;
241  dimuon->vtx_neg = m_pos;
242  dimuon->proj_target_pos = posTrack->getTargetPos();
243  dimuon->proj_dump_pos = posTrack->getDumpPos();
244  dimuon->proj_target_neg = negTrack->getTargetPos();
245  dimuon->proj_dump_neg = negTrack->getDumpPos();
246  dimuon->chisq_target = posTrack->getChisqTarget() + negTrack->getChisqTarget();
247  dimuon->chisq_dump = posTrack->getChisqDump() + negTrack->getChisqDump();
248  dimuon->chisq_upstream = posTrack->getChisqUpstream() + negTrack->getChisqUpstream();
249  dimuon->chisq_single = 0.;
250  dimuon->chisq_vx = 0.;
251  dimuon->calcVariables();
252 
253  return true;
254 }
#define M_MU
Definition: GlobalConsts.h:12
virtual const std::string Name() const
Returns the name of this module.
Definition: Fun4AllBase.h:23
PHBoolean addNode(PHNode *)
virtual const SQDimuon * at(const size_t id) const =0
virtual size_t size() const =0
virtual void push_back(const SQDimuon *dim)=0
An SQ interface class to hold one true or reconstructed dimuon.
Definition: SQDimuon.h:8
virtual void set_rec_dimuon_id(const int a)=0
virtual int get_track_id_neg() const =0
Return the track ID of the negative track.
virtual int get_track_id_pos() const =0
Return the track ID of the positive track.
virtual TVector3 get_pos() const =0
Return the dimuon position at vertex.
virtual const SQTrack * at(const size_t id) const =0
virtual size_t size() const =0
An SQ interface class to hold one true or reconstructed track.
Definition: SQTrack.h:8
virtual int get_track_id() const =0
Return the track ID, which is unique per event(?).
virtual TVector3 get_pos_vtx() const =0
Return the track position at vertex.
virtual int get_rec_track_id() const =0
Return the track ID of associated reconstructed track. Valid only if this object holds truth track in...
int Init(PHCompositeNode *topNode)
int End(PHCompositeNode *topNode)
Called at the end of all processing.
int process_event(PHCompositeNode *topNode)
SQTruthVertexing(const std::string &name="SQTruthVertexing")
int InitRun(PHCompositeNode *topNode)
Double_t chisq_single
Definition: SRecEvent.h:398
TVector3 vtx
3-vector vertex position
Definition: SRecEvent.h:379
void calcVariables(int choice=0)
Calculate the kinematic vairables, 0 = vertex, 1 = target, 2 = dump.
Definition: SRecEvent.cxx:633
Double_t chisq_dump
Definition: SRecEvent.h:406
Double_t chisq_upstream
Definition: SRecEvent.h:407
TVector3 proj_target_neg
Definition: SRecEvent.h:386
Double_t chisq_vx
Definition: SRecEvent.h:402
Int_t trackID_pos
Index of muon track used in host SRecEvent.
Definition: SRecEvent.h:361
Double_t chisq_target
Chisq of three test position.
Definition: SRecEvent.h:405
TVector3 proj_dump_pos
Definition: SRecEvent.h:385
TVector3 proj_dump_neg
Definition: SRecEvent.h:387
Int_t trackID_neg
Definition: SRecEvent.h:362
TLorentzVector p_neg
Definition: SRecEvent.h:366
Double_t chisq_kf
Vertex fit chisqs.
Definition: SRecEvent.h:401
TLorentzVector p_pos
4-momentum of the muon tracks after vertex fit
Definition: SRecEvent.h:365
TVector3 vtx_neg
Definition: SRecEvent.h:381
TVector3 vtx_pos
Definition: SRecEvent.h:380
TVector3 proj_target_pos
Track projections at different location.
Definition: SRecEvent.h:384
void insertDimuon(SRecDimuon dimuon)
Insert dimuon.
Definition: SRecEvent.h:470
Int_t getNTracks()
Get tracks.
Definition: SRecEvent.h:455
SRecTrack & getTrack(Int_t i)
Definition: SRecEvent.h:456
Int_t getNDimuons()
Get dimuons.
Definition: SRecEvent.h:462
TVector3 getTargetPos()
Definition: SRecEvent.h:188
Int_t getCharge() const
Gets.
Definition: SRecEvent.h:101
Double_t getChisqTarget()
Definition: SRecEvent.h:199
TVector3 getVertexPos()
Definition: SRecEvent.h:196
void setVertexPos(TVector3 pos)
Definition: SRecEvent.h:208
Double_t getChisqDump()
Definition: SRecEvent.h:198
Double_t getChisqUpstream()
Definition: SRecEvent.h:200
void setVertexMom(TVector3 mom)
Definition: SRecEvent.h:214
TVector3 getVertexMom()
Definition: SRecEvent.h:197
TVector3 getDumpPos()
Get mom/pos at a given location.
Definition: SRecEvent.h:186
Double_t getChisqVertex()
Definition: SRecEvent.h:183
void setChisqVertex(Double_t chisq)
Definition: SRecEvent.h:218