Class Reference for E1039 Core & Analysis Software
SQVertexing.cxx
Go to the documentation of this file.
1 #include "SQVertexing.h"
2 
3 #include <phfield/PHFieldConfig_v3.h>
4 #include <phfield/PHFieldUtility.h>
5 #include <phgeom/PHGeomUtility.h>
6 
8 #include <phool/PHNodeIterator.h>
9 #include <phool/PHIODataNode.h>
10 #include <phool/PHRandomSeed.h>
11 #include <phool/getClass.h>
12 #include <phool/recoConsts.h>
13 #include <interface_main/SQEvent.h>
18 #include <GenFit/FieldManager.h>
19 #include <GenFit/MaterialEffects.h>
20 #include <GenFit/TGeoMaterialInterface.h>
21 
22 #include "SRecEvent.h"
23 #include "GFTrack.h"
24 
25 #include <iostream>
26 #include <vector>
27 
28 namespace
29 {
30  //static flag to indicate the initialized has been done
31  static bool inited = false;
32 
33  static double Z_TARGET;
34  static double Z_DUMP;
35  static double Z_UPSTREAM;
36 
37  static double X_BEAM;
38  static double Y_BEAM;
39  static double SIGX_BEAM;
40  static double SIGY_BEAM;
41 
42  //initialize global variables
43  void initGlobalVariables()
44  {
45  if(!inited)
46  {
47  inited = true;
48 
50  Z_TARGET = rc->get_DoubleFlag("Z_TARGET");
51  Z_DUMP = rc->get_DoubleFlag("Z_DUMP");
52  Z_UPSTREAM = rc->get_DoubleFlag("Z_UPSTREAM");
53 
54  X_BEAM = rc->get_DoubleFlag("X_BEAM");
55  Y_BEAM = rc->get_DoubleFlag("Y_BEAM");
56  SIGX_BEAM = rc->get_DoubleFlag("SIGX_BEAM");
57  SIGY_BEAM = rc->get_DoubleFlag("SIGY_BEAM");
58  }
59  }
60 };
61 
62 SQVertexing::SQVertexing(const std::string& name, int sign1, int sign2):
63  SubsysReco(name),
64  legacyContainer_in(false),
65  legacyContainer_out(false),
66  enableSingleRetracking(false),
67  gfield(nullptr),
68  geom_file_name(""),
69  recEvent(nullptr),
70  recTrackVec(nullptr),
71  recDimuonVec(nullptr),
72  charge1(sign1),
73  charge2(sign2)
74 {}
75 
77 {}
78 
80 {
81  initGlobalVariables();
83 }
84 
86 {
87  int ret = GetNodes(topNode);
88  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
89 
90  ret = MakeNodes(topNode);
91  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
92 
93  ret = InitField(topNode);
94  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
95 
96  ret = InitGeom(topNode);
97  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
98 
99  //Initialize random seed
101  rndm.SetSeed(rc->get_IntFlag("RUNNUMBER"));
102 
104 }
105 
107 {
108  if(legacyContainer_out) recEvent->clearDimuons();
109 
110  std::vector<int> trackIDs1;
111  std::vector<int> trackIDs2;
112  int nTracks = legacyContainer_in ? recEvent->getNTracks() : recTrackVec->size();
113  if(Verbosity() > 10) std::cout << "SQVertexing::process_event(): nTracks = " << nTracks << std::endl;
114  for(int i = 0; i < nTracks; ++i)
115  {
116  SRecTrack* recTrack = legacyContainer_in ? &(recEvent->getTrack(i)) : dynamic_cast<SRecTrack*>(recTrackVec->at(i));
117  if(!recTrack->isKalmanFitted()) continue;
118 
119  if(enableSingleRetracking)
120  {
121  TVector3 pos, mom;
122  double chi2 = refitTrkToVtx(recTrack, Z_TARGET, &pos, &mom);
123  if(chi2 > 0)
124  {
125  recTrack->setChisqTarget(chi2);
126  recTrack->setTargetPos(pos);
127  recTrack->setTargetMom(mom);
128  }
129  }
130 
131  if(recTrack->getCharge() == charge1) trackIDs1.push_back(i);
132  if(recTrack->getCharge() == charge2) trackIDs2.push_back(i);
133 
134  }
135  if(trackIDs1.empty() || trackIDs2.empty()) return Fun4AllReturnCodes::EVENT_OK;
136  if(Verbosity() > 10) std::cout << " N of trackIDs1 & trackIDs1 = " << trackIDs1.size() << " & " << trackIDs2.size() << std::endl;
137 
138  for(int i = 0; i < trackIDs1.size(); ++i)
139  {
140  for(int j = charge1 == charge2 ? i+1 : 0; j < trackIDs2.size(); ++j)
141  {
142  //A protection, probably not really needed
143  if(trackIDs1[i] == trackIDs2[j]) continue;
144 
145  SRecTrack* recTrack1 = legacyContainer_in ? &(recEvent->getTrack(trackIDs1[i])) : dynamic_cast<SRecTrack*>(recTrackVec->at(trackIDs1[i]));
146  SRecTrack* recTrack2 = legacyContainer_in ? &(recEvent->getTrack(trackIDs2[j])) : dynamic_cast<SRecTrack*>(recTrackVec->at(trackIDs2[j]));
147 
148  SRecDimuon recDimuon;
149  if(processOneDimuon(recTrack1, recTrack2, recDimuon))
150  {
151  //recDimuon.set_rec_dimuon_id(legacyContainer ? recEvent->getNDimuons() : recDimuonVec->size());
152  recDimuon.set_track_id_pos(trackIDs1[i]);
153  recDimuon.set_track_id_neg(trackIDs2[j]);
154 
155  if(legacyContainer_out)
156  recEvent->insertDimuon(recDimuon);
157  else
158  recDimuonVec->push_back(&recDimuon);
159  }
160  }
161  }
162 
164 }
165 
167 {
169 }
170 
171 int SQVertexing::GetNodes(PHCompositeNode* topNode)
172 {
173  if(legacyContainer_in || legacyContainer_out)
174  {
175  recEvent = findNode::getClass<SRecEvent>(topNode, "SRecEvent");
176  if(!recEvent)
177  {
178  std::cerr << Name() << ": failed finding SRecEvent node, abort." << std::endl;
180  }
181  }
182  else
183  {
184  recTrackVec = findNode::getClass<SQTrackVector>(topNode, "SQRecTrackVector");
185  if(!recTrackVec)
186  {
187  std::cerr << Name() << ": failed finding SQRecTrackVector node, abort." << std::endl;
189  }
190  }
191 
193 }
194 
195 int SQVertexing::MakeNodes(PHCompositeNode* topNode)
196 {
197  if(!legacyContainer_out)
198  {
199  PHNodeIterator iter(topNode);
200  PHCompositeNode* dstNode = static_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
201  if(!dstNode)
202  {
203  std::cerr << Name() << ": cannot locate DST node, abort." << std::endl;
205  }
206 
207  recDimuonVec = new SQDimuonVector_v1();
208  TString dimuonVecName = TString::Format("SQRecDimuonVector_%s%s", (charge1 == 1 ? "P" : "M"), (charge2 == 1 ? "P" : "M"));
209  dstNode->addNode(new PHIODataNode<PHObject>(recDimuonVec, dimuonVecName.Data(), "PHObject"));
210  }
212 }
213 
214 double SQVertexing::swimTrackToVertex(SQGenFit::GFTrack& track, double z, TVector3* pos, TVector3* mom)
215 {
216  double chi2 = track.swimToVertex(z, pos, mom);
217  return chi2;
218 }
219 
220 double SQVertexing::refitTrkToVtx(SQGenFit::GFTrack& track, double z, TVector3* pos, TVector3* mom)
221 {
222  if(Verbosity() > 10) std::cout << "SQVertexing::refitTrkToVtx(): z = " << z << std::endl;
223  gfield->setOffset(0.);
224 
225  double z_offset_prev = 1.E9;
226  double z_offset_curr = 0.;
227  double chi2 = 0.;
228  const int n_iter_max = 100;
229  int n_iter = 0;
230  while(fabs(z_offset_curr - z_offset_prev) > 0.1)
231  {
232  chi2 = swimTrackToVertex(track, z, pos, mom);
233 
234  //update scatter plane location
235  z_offset_prev = z_offset_curr;
236  z_offset_curr = calcZsclp(mom->Mag());
237 
238  gfield->setOffset(-z_offset_curr);
239  if(Verbosity() > 20) std::cout << " i_iter = " << n_iter << ": p = " << mom->Mag() << ", dz = " << z_offset_curr << " - " << z_offset_prev << " = " << z_offset_curr - z_offset_prev << std::endl;
240  if (++n_iter == n_iter_max) {
241  std::cout << "!WARNING! SQVertexing::refitTrkToVtx(): Give up the iteration." << std::endl;
242  break;
243  }
244  }
245  if(Verbosity() > 10) {
246  std::cout << " n_iter = " << n_iter << ", pos = (" << pos->X() << ", " << pos->Y() << ", " << pos->Z() << "), mom = (" << mom->X() << ", " << mom->Y() << ", " << mom->Z() << ")" << std::endl;
247  }
248 
249  //re-adjust FMag bend plane here
250  gfield->setOffset(0.);
251 
252  return chi2;
253 }
254 
255 double SQVertexing::refitTrkToVtx(SRecTrack* track, double z, TVector3* pos, TVector3* mom)
256 {
257  SQGenFit::GFTrack gtrk(*track);
258  return refitTrkToVtx(gtrk, z, pos, mom);
259 }
260 
261 double SQVertexing::calcZsclp(double p)
262 {
263  if(p < 5.) p = 5.;
264  else if(p > 120.) p = 120.;
265 
266  return 301.84 - 1.27137*p + 0.0218294*p*p - 0.000170711*p*p*p + 4.94683e-07*p*p*p*p - 271.;
267 }
268 
269 bool SQVertexing::processOneDimuon(SRecTrack* track1, SRecTrack* track2, SRecDimuon& dimuon)
270 {
271  if(Verbosity() > 10) std::cout << " SQVertexing::processOneDimuon(): " << std::endl;
272  if(Verbosity() > 20) {
273  track1->printGF();
274  track2->printGF();
275  }
276 
277  //Pre-calculated variables
278  dimuon.proj_target_pos = track1->getTargetPos();
279  dimuon.proj_dump_pos = track1->getDumpPos();
280  dimuon.proj_target_neg = track2->getTargetPos();
281  dimuon.proj_dump_neg = track2->getDumpPos();
282  dimuon.chisq_target = track1->getChisqTarget() + track2->getChisqTarget();
283  dimuon.chisq_dump = track1->getChisqDump() + track2->getChisqDump();
284  dimuon.chisq_upstream = track1->getChisqUpstream() + track2->getChisqUpstream();
285  dimuon.chisq_single = 0.; //no longer used
286  dimuon.chisq_vx = 0.; //no longer used
287 
288  //Vertexing part
289  SQGenFit::GFTrack gtrk1(*track1);
290  SQGenFit::GFTrack gtrk2(*track2);
291  double z_vtx = findDimuonZVertex(dimuon, gtrk1, gtrk2);
292  dimuon.vtx.SetXYZ(0., 0., z_vtx);
293 
294  //Vertex info based on the dimuon vertex
295  //TODO: consider using addjustable bend-plane for vertex test as well
296  TVector3 pos, mom;
297  dimuon.chisq_kf = swimTrackToVertex(gtrk1, z_vtx, &pos, &mom);
298  dimuon.p_pos.SetVectM(mom, M_MU);
299  dimuon.vtx_pos = pos;
300 
301  dimuon.chisq_kf += swimTrackToVertex(gtrk2, z_vtx, &pos, &mom);
302  dimuon.p_neg.SetVectM(mom, M_MU);
303  dimuon.vtx_neg = pos;
304 
305  //Test target hypothesis
306  refitTrkToVtx(gtrk1, Z_TARGET, &pos, &mom);
307  dimuon.p_pos_target.SetVectM(mom, M_MU);
308  refitTrkToVtx(gtrk2, Z_TARGET, &pos, &mom);
309  dimuon.p_neg_target.SetVectM(mom, M_MU);
310 
311  //Test dump hypothesis
312  //TODO: consider using addjustable bend-plane for dump test as well
313  swimTrackToVertex(gtrk1, Z_DUMP, &pos, &mom);
314  dimuon.p_pos_dump.SetVectM(mom, M_MU);
315  swimTrackToVertex(gtrk2, Z_DUMP, &pos, &mom);
316  dimuon.p_neg_dump.SetVectM(mom, M_MU);
317 
318  //Randomly flip the sign of Px of one of the muons if processing like-sign muons
319  if(charge1 == charge2)
320  {
321  if(rndm.Rndm() < 0.5)
322  {
323  dimuon.p_pos.SetPx(-dimuon.p_pos.Px());
324  dimuon.p_pos_target.SetPx(-dimuon.p_pos_target.Px());
325  dimuon.p_pos_dump.SetPx(-dimuon.p_pos_dump.Px());
326  }
327  else
328  {
329  dimuon.p_neg.SetPx(-dimuon.p_neg.Px());
330  dimuon.p_neg_target.SetPx(-dimuon.p_neg_target.Px());
331  dimuon.p_neg_dump.SetPx(-dimuon.p_neg_dump.Px());
332  }
333  }
334 
335  return true;
336 }
337 
338 double SQVertexing::findDimuonZVertex(SRecDimuon& dimuon, SQGenFit::GFTrack& track1, SQGenFit::GFTrack& track2)
339 {
340  //TODO: consider using addjustable bend-plane for vertex finding as well
341  double stepsize[3] = {25., 5., 1.};
342 
343  double z_min = 300.;
344  double chi2_min = 1.E9;
345  for(int i = 0; i < 3; ++i)
346  {
347  double z_start, z_end;
348  if(i == 0)
349  {
350  z_start = z_min;
351  z_end = Z_UPSTREAM;
352  }
353  else
354  {
355  z_start = z_min + stepsize[i-1];
356  z_end = z_min - stepsize[i-1];
357  }
358 
359  for(double z = z_start; z > z_end; z = z - stepsize[i])
360  {
361  double chi2_1 = swimTrackToVertex(track1, z);
362  double chi2_2 = swimTrackToVertex(track2, z);
363  double chi2 = chi2_1 > 0 && chi2_2 > 0 ? chi2_1 + chi2_2 : 1.E9;
364  if(chi2 < chi2_min)
365  {
366  z_min = z;
367  chi2_min = chi2;
368  }
369  }
370  }
371 
372  return z_min;
373 }
374 
375 int SQVertexing::InitField(PHCompositeNode* topNode)
376 {
377  try
378  {
379  gfield = dynamic_cast<SQGenFit::GFField*>(genfit::FieldManager::getInstance()->getField());
380  }
381  catch(const std::exception& e)
382  {
383  std::cout << "SQVertexing::InitGeom - Caught an exception " << std::endl;
384  gfield = nullptr;
385  }
386 
387  if(gfield == nullptr)
388  {
390 
391  std::unique_ptr<PHFieldConfig> default_field_cfg(new PHFieldConfig_v3(rc->get_CharFlag("fMagFile"), rc->get_CharFlag("kMagFile"), rc->get_DoubleFlag("FMAGSTR"), rc->get_DoubleFlag("KMAGSTR"), 5.));
392  PHField* phfield = PHFieldUtility::GetFieldMapNode(default_field_cfg.get(), topNode, 0);
393 
394  std::cout << "SQVertexing::InitGeom - creating new GenFit field map" << std::endl;
395  gfield = new SQGenFit::GFField(phfield);
396  genfit::FieldManager::getInstance()->init(gfield);
397  }
398  else
399  {
400  std::cout << "SQVertexing::InitGeom - reading existing GenFit field map" << std::endl;
401  }
402 
404 }
405 
406 int SQVertexing::InitGeom(PHCompositeNode* topNode)
407 {
408  if(geom_file_name != "")
409  {
410  std::cout << "SQVertexing::InitGeom - create geom from " << geom_file_name << std::endl;
411 
412  int ret = PHGeomUtility::ImportGeomFile(topNode, geom_file_name);
413  genfit::MaterialEffects::getInstance()->init(new genfit::TGeoMaterialInterface());
414 
415  if(ret != Fun4AllReturnCodes::EVENT_OK) return ret;
416  }
417  else
418  {
419  std::cout << "SQVertexing::InitGeom - rely on existing TGeo geometry" << std::endl;
420  }
421 
423 }
#define M_MU
Definition: GlobalConsts.h:12
virtual const std::string Name() const
Returns the name of this module.
Definition: Fun4AllBase.h:23
virtual int Verbosity() const
Gets the verbosity of this module.
Definition: Fun4AllBase.h:64
PHBoolean addNode(PHNode *)
PHFieldConfig_v3 implements field configuration information for input a field map file.
static PHField * GetFieldMapNode(const PHFieldConfig *default_config=nullptr, PHCompositeNode *topNode=nullptr, const int verbosity=0)
Get transient PHField from DST nodes. If not found, make a new one based on default_config.
transient DST object for field storage and access
Definition: PHField.h:14
virtual double get_DoubleFlag(const std::string &name) const
Definition: PHFlag.cc:49
virtual int get_IntFlag(const std::string &name) const
Definition: PHFlag.cc:117
virtual const std::string get_CharFlag(const std::string &flag) const
Definition: PHFlag.cc:13
static int ImportGeomFile(PHCompositeNode *topNode, const std::string &geometry_file)
TGeo ROOT/GDML/Macro file -> DST node with automatic file type discrimination based on file names.
virtual void push_back(const SQDimuon *dim)=0
void setOffset(double offset=0.)
Definition: GFField.h:21
double swimToVertex(double z, TVector3 *pos=nullptr, TVector3 *mom=nullptr, TMatrixDSym *cov=nullptr)
Definition: GFTrack.cxx:295
virtual const SQTrack * at(const size_t id) const =0
virtual size_t size() const =0
int InitRun(PHCompositeNode *topNode)
Definition: SQVertexing.cxx:85
int process_event(PHCompositeNode *topNode)
SQVertexing(const std::string &name="SQVertexing", int sign1=1, int sign2=-1)
Definition: SQVertexing.cxx:62
int Init(PHCompositeNode *topNode)
Definition: SQVertexing.cxx:79
int End(PHCompositeNode *topNode)
Called at the end of all processing.
Double_t chisq_single
Definition: SRecEvent.h:398
TVector3 vtx
3-vector vertex position
Definition: SRecEvent.h:379
TLorentzVector p_neg_target
Definition: SRecEvent.h:370
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
virtual void set_track_id_pos(const int a)
Definition: SRecEvent.h:321
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
TLorentzVector p_pos_target
Track momentum projections at different location.
Definition: SRecEvent.h:369
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
TLorentzVector p_pos_dump
Definition: SRecEvent.h:371
TLorentzVector p_neg_dump
Definition: SRecEvent.h:372
void insertDimuon(SRecDimuon dimuon)
Insert dimuon.
Definition: SRecEvent.h:470
Int_t getNTracks()
Get tracks.
Definition: SRecEvent.h:455
void clearDimuons()
Clear the dimuon list.
Definition: SRecEvent.cxx:824
SRecTrack & getTrack(Int_t i)
Definition: SRecEvent.h:456
TVector3 getTargetPos()
Definition: SRecEvent.h:188
Int_t getCharge() const
Gets.
Definition: SRecEvent.h:101
Double_t getChisqTarget()
Definition: SRecEvent.h:199
Double_t getChisqDump()
Definition: SRecEvent.h:198
void setTargetPos(TVector3 pos)
Definition: SRecEvent.h:205
void setChisqTarget(Double_t chisq)
Definition: SRecEvent.h:216
Double_t getChisqUpstream()
Definition: SRecEvent.h:200
void printGF(std::ostream &os=std::cout) const
Definition: SRecEvent.cxx:593
void setTargetMom(TVector3 mom)
Definition: SRecEvent.h:211
Bool_t isKalmanFitted()
Fit status.
Definition: SRecEvent.h:147
TVector3 getDumpPos()
Get mom/pos at a given location.
Definition: SRecEvent.h:186
static recoConsts * instance()
Definition: recoConsts.cc:7