Class Reference for E1039 Core & Analysis Software
AnaSimDst.cc
Go to the documentation of this file.
1 #include <iomanip>
2 #include <TFile.h>
3 #include <TTree.h>
8 #include <ktracker/SRecEvent.h>
10 #include <phool/getClass.h>
11 #include <UtilAna/UtilDimuon.h>
12 #include "AnaSimDst.h"
13 using namespace std;
14 
16 {
17  ;
18 }
19 
21 {
23 }
24 
26 {
27  mi_evt = findNode::getClass<SQEvent >(topNode, "SQEvent");
28  mi_evt_true = findNode::getClass<SQMCEvent >(topNode, "SQMCEvent");
29  mi_vec_trk = findNode::getClass<SQTrackVector >(topNode, "SQTruthTrackVector");
30  mi_vec_dim = findNode::getClass<SQDimuonVector>(topNode, "SQTruthDimuonVector");
31  mi_srec = findNode::getClass<SRecEvent >(topNode, "SRecEvent");
32  if (!mi_evt || !mi_evt_true || !mi_vec_trk || !mi_vec_dim) return Fun4AllReturnCodes::ABORTEVENT;
33  if (!mi_srec) {
34  cout << "The SRecEvent node cannot be found in DST and thus won't be analyzed." << endl;
35  }
36 
37  mo_file = new TFile("sim_tree.root", "RECREATE");
38  mo_tree = new TTree("tree", "Created by AnaSimDst");
39  mo_tree->Branch("evt" , &mo_evt);
40  mo_tree->Branch("trk_true", &mo_trk_true);
41  mo_tree->Branch("trk_reco", &mo_trk_reco);
42  mo_tree->Branch("dim_true", &mo_dim_true);
43  mo_tree->Branch("dim_reco", &mo_dim_reco);
44 
46 }
47 
49 {
50  static unsigned int n_evt = 0;
51  if (++n_evt % 100000 == 0) cout << n_evt << endl;
52  else if (n_evt % 10000 == 0) cout << " . " << flush;
53 
57  mo_evt.proc_id = mi_evt_true->get_process_id();
58  for (int ii = 0; ii < 4; ii++) {
59  mo_evt.par_id [ii] = mi_evt_true->get_particle_id (ii);
60  mo_evt.par_mom[ii] = mi_evt_true->get_particle_momentum(ii);
61  }
62  mo_evt.weight = mi_evt_true->get_weight();
63  mo_evt.trig_bits = mi_evt->get_trigger();
64  mo_evt.n_dim_true = mi_vec_dim->size();
65 
66  if (mi_srec) {
67  mo_evt.rec_stat = mi_srec->getRecStatus();
68  mo_evt.n_dim_reco = mi_srec->getNDimuons();
69  } else {
70  mo_evt.rec_stat = 0;
71  mo_evt.n_dim_reco = 0;
72  }
73 
77  IdMap_t id_trk_t2r;
78  if (mi_srec) FindTrackRelation(id_trk_t2r);
79  mo_trk_true.clear();
80  mo_trk_reco.clear();
81  for (unsigned int ii = 0; ii < mi_vec_trk->size(); ii++) {
82  SQTrack* trk = mi_vec_trk->at(ii);
83  TrackData td;
84  td.charge = trk->get_charge();
85  td.pos_vtx = trk->get_pos_vtx();
86  td.mom_vtx = trk->get_mom_vtx();
87  mo_trk_true.push_back(td);
88 
89  if (mi_srec) {
90  TrackData tdr;
91  if (id_trk_t2r[ii] >= 0) {
92  SRecTrack* trk_reco = &mi_srec->getTrack(id_trk_t2r[ii]);
93  tdr.charge = trk_reco->getCharge();
94  tdr.pos_vtx = trk_reco->getVertex();
95  tdr.mom_vtx = trk_reco->getMomentumVertex();
96  }
97  mo_trk_reco.push_back(tdr);
98  }
99  }
100 
104  IdMap_t id_dim_t2r;
105  if (mi_srec) FindDimuonRelation(id_dim_t2r);
106  mo_dim_true.clear();
107  mo_dim_reco.clear();
108  for (unsigned int ii = 0; ii < mi_vec_dim->size(); ii++) {
109  SQDimuon* dim = mi_vec_dim->at(ii);
110  DimuonData dd;
111  dd.pdg_id = dim->get_pdg_id();
112  dd.pos = dim->get_pos();
113  dd.mom = dim->get_mom();
114  dd.mom_pos = dim->get_mom_pos();
115  dd.mom_neg = dim->get_mom_neg();
116  UtilDimuon::CalcVar(dim, dd.mass, dd.pT, dd.x1, dd.x2, dd.xF, dd.costh, dd.phi);
117  mo_dim_true.push_back(dd);
118 
119  if (mi_srec) {
120  DimuonData ddr;
121  if (id_dim_t2r[ii] >= 0) {
122  SRecDimuon dim_reco = mi_srec->getDimuon(id_dim_t2r[ii]);
123  ddr.pos = dim_reco.vtx;
124  ddr.mom = dim_reco.p_pos + dim_reco.p_neg;
125  ddr.mom_pos = dim_reco.p_pos;
126  ddr.mom_neg = dim_reco.p_neg;
127  ddr.x1 = dim_reco.x1;
128  ddr.x2 = dim_reco.x2;
129  }
130  mo_dim_reco.push_back(ddr);
131  }
132  }
133 
134  mo_tree->Fill();
136 }
137 
139 {
140  mo_file->cd();
141  mo_file->Write();
142  mo_file->Close();
144 }
145 
146 void AnaSimDst::FindTrackRelation(IdMap_t& id_map)
147 {
148  id_map.clear();
149  for (unsigned int i_true = 0; i_true < mi_vec_trk->size(); i_true++) {
150  SQTrack* trk_true = mi_vec_trk->at(i_true);
151  int ch_true = trk_true->get_charge();
152  double mom_true = trk_true->get_mom_vtx().Mag();
153 
154  int i_reco_best = -1;
155  double mom_diff_best;
156  for (int i_reco = 0; i_reco < mi_srec->getNTracks(); i_reco++) {
157  SRecTrack* trk_reco = &mi_srec->getTrack(i_reco);
158  if (trk_reco->getCharge() != ch_true) continue;
159  double mom_diff = fabs(trk_reco->getMomentumVertex().Mag() - mom_true);
160  if (i_reco_best < 0 || mom_diff < mom_diff_best) {
161  i_reco_best = i_reco;
162  mom_diff_best = mom_diff;
163  }
164  }
165  id_map[i_true] = i_reco_best;
166  }
167 }
168 
169 void AnaSimDst::FindDimuonRelation(IdMap_t& id_map)
170 {
171  id_map.clear();
172  for (unsigned int i_true = 0; i_true < mi_vec_dim->size(); i_true++) {
173  SQDimuon* dim_true = mi_vec_dim->at(i_true);
174  double mass_true = dim_true->get_mom().M();
175 
176  int i_reco_best = -1;
177  double mass_diff_best;
178  for (int i_reco = 0; i_reco < mi_srec->getNDimuons(); i_reco++) {
179  SRecDimuon dim_reco = mi_srec->getDimuon(i_reco);
180  double mass_diff = fabs(dim_reco.mass - mass_true);
181  if (i_reco_best < 0 || mass_diff < mass_diff_best) {
182  i_reco_best = i_reco;
183  mass_diff_best = mass_diff;
184  }
185  }
186  id_map[i_true] = i_reco_best;
187  }
188 }
int process_event(PHCompositeNode *topNode)
Definition: AnaSimDst.cc:48
int End(PHCompositeNode *topNode)
Called at the end of all processing.
Definition: AnaSimDst.cc:138
int Init(PHCompositeNode *topNode)
Definition: AnaSimDst.cc:20
int InitRun(PHCompositeNode *topNode)
Definition: AnaSimDst.cc:25
virtual const SQDimuon * at(const size_t id) const =0
virtual size_t size() const =0
An SQ interface class to hold one true or reconstructed dimuon.
Definition: SQDimuon.h:8
virtual TLorentzVector get_mom_neg() const =0
Return the momentum of the negative track at vertex.
virtual int get_pdg_id() const =0
Return the GPD ID of parent particle. It is valid only for true dimuon.
virtual TLorentzVector get_mom_pos() const =0
Return the momentum of the positive track at vertex.
virtual TLorentzVector get_mom() const =0
Return the dimuon momentum at vertex.
virtual TVector3 get_pos() const =0
Return the dimuon position at vertex.
virtual bool get_trigger(const SQEvent::TriggerMask i) const =0
Return the trigger bit (fired or not) of the selected trigger channel.
virtual int get_particle_id(const int i) const =0
Return the particle ID of the primary process, where i=0...3 for "0 + 1 -> 2 + 3".
virtual double get_weight() const =0
Return the event weight.
virtual TLorentzVector get_particle_momentum(const int i) const =0
Return the particle momentum of the primary process, where i=0...3 for "0 + 1 -> 2 + 3".
virtual int get_process_id() const =0
Return the primary process ID.
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 TVector3 get_pos_vtx() const =0
Return the track position at vertex.
virtual int get_charge() const =0
Return the charge, i.e. +1 or -1.
virtual TLorentzVector get_mom_vtx() const =0
Return the track momentum at vertex.
TVector3 vtx
3-vector vertex position
Definition: SRecEvent.h:379
TLorentzVector p_neg
Definition: SRecEvent.h:366
Double_t x2
Definition: SRecEvent.h:394
TLorentzVector p_pos
4-momentum of the muon tracks after vertex fit
Definition: SRecEvent.h:365
Double_t x1
Definition: SRecEvent.h:393
Double_t mass
Kinematic variables.
Definition: SRecEvent.h:390
Int_t getRecStatus()
Definition: SRecEvent.h:445
SRecDimuon & getDimuon(Int_t i)
Definition: SRecEvent.h:463
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
Int_t getCharge() const
Gets.
Definition: SRecEvent.h:101
TVector3 getVertex()
Definition: SRecEvent.h:181
TLorentzVector getMomentumVertex()
Get the vertex info.
Definition: SRecEvent.cxx:340
void CalcVar(const SQDimuon *dim, double &mass, double &pT, double &x1, double &x2, double &xF)
Calculate the key kinematic variables of dimuon.
Definition: UtilDimuon.cc:70
TLorentzVector mom_pos
Definition: TreeData.h:37
TLorentzVector mom
Definition: TreeData.h:36
int pdg_id
Definition: TreeData.h:34
double x1
Definition: TreeData.h:41
double phi
Definition: TreeData.h:45
TLorentzVector mom_neg
Definition: TreeData.h:38
double x2
Definition: TreeData.h:42
double mass
Definition: TreeData.h:39
double xF
Definition: TreeData.h:43
double pT
Definition: TreeData.h:40
double costh
Definition: TreeData.h:44
TVector3 pos
Definition: TreeData.h:35
double weight
Definition: TreeData.h:10
int trig_bits
Definition: TreeData.h:11
int proc_id
Definition: TreeData.h:7
int rec_stat
Definition: TreeData.h:12
TLorentzVector par_mom[4]
Definition: TreeData.h:9
int n_dim_true
Definition: TreeData.h:13
int par_id[4]
Definition: TreeData.h:8
int n_dim_reco
Definition: TreeData.h:14
TVector3 pos_vtx
Definition: TreeData.h:24
int charge
Definition: TreeData.h:23
TLorentzVector mom_vtx
Definition: TreeData.h:25