Class Reference for E1039 Core & Analysis Software
AnaEmbeddedData.cc
Go to the documentation of this file.
1 #include <iomanip>
2 #include <TFile.h>
3 #include <TTree.h>
4 #include <geom_svc/GeomSvc.h>
10 #include <ktracker/SRecEvent.h>
12 #include <phool/getClass.h>
13 #include <UtilAna/UtilSQHit.h>
14 #include <UtilAna/UtilDimuon.h>
15 #include "AnaEmbeddedData.h"
16 using namespace std;
17 
18 AnaEmbeddedData::AnaEmbeddedData(const std::string name)
19  : SubsysReco(name)
20 {
21  ;
22 }
23 
25 {
27 }
28 
30 {
31  mi_evt = findNode::getClass<SQEvent >(topNode, "SQEvent");
32  mi_vec_hit = findNode::getClass<SQHitVector >(topNode, "SQHitVector");
33  if (!mi_evt || !mi_vec_hit) {
34  cout << PHWHERE << ": Cannot find SQEvent and/or SQHitVector." << endl;
36  }
37 
38  mi_evt_emb = findNode::getClass<SQEvent >(topNode, "SQEventEmb");
39  mi_evt_mc = findNode::getClass<SQMCEvent >(topNode, "SQMCEvent");
40  mi_evt_mc_emb = findNode::getClass<SQMCEvent >(topNode, "SQMCEventEmb");
41  mi_vec_trk = findNode::getClass<SQTrackVector >(topNode, "SQTruthTrackVector");
42  mi_vec_dim = findNode::getClass<SQDimuonVector>(topNode, "SQTruthDimuonVector");
43  mi_srec = findNode::getClass<SRecEvent >(topNode, "SRecEvent");
44  if (!mi_evt_emb ) cout << "No SQEventEmb node." << endl;
45  if (!mi_evt_mc ) cout << "No SQMCEvent node." << endl;
46  if (!mi_evt_mc_emb) cout << "No SQMCEventEmb node." << endl;
47  if (!mi_vec_trk ) cout << "No SQTruthTrackVector node." << endl;
48  if (!mi_vec_dim ) cout << "No SQTruthDimuonVector node." << endl;
49  if (!mi_srec ) cout << "No SRecEvent node." << endl;
50 
51  mo_file = new TFile("ana_tree.root", "RECREATE");
52  mo_tree = new TTree("tree", "Created by AnaEmbeddedData");
53  mo_tree->Branch("evt" , &mo_evt);
54  mo_tree->Branch("trk_true", &mo_trk_true);
55  mo_tree->Branch("trk_reco", &mo_trk_reco);
56  mo_tree->Branch("dim_true", &mo_dim_true);
57  mo_tree->Branch("dim_reco", &mo_dim_reco);
58 
60 }
61 
63 {
67  mo_evt.job_id = mi_evt->get_spill_id(); // Ppill ID in simulation means job ID
68  mo_evt.event_id = mi_evt->get_event_id();
69  mo_evt.trig_bits = mi_evt->get_trigger();
70  mo_evt.rfp01 = mi_evt->get_qie_rf_intensity(+1);
71  mo_evt.rfp00 = mi_evt->get_qie_rf_intensity( 0);
72  mo_evt.rfm01 = mi_evt->get_qie_rf_intensity(-1);
73 
74  if (mi_evt_mc) {
75  mo_evt.weight = mi_evt_mc->get_weight();
76  }
77  if (mi_srec) {
78  mo_evt.rec_stat = mi_srec->getRecStatus();
79  }
80 
84  map<int, int> list_cnt; // [det ID] -> hit count
85  for (SQHitVector::ConstIter it = mi_vec_hit->begin(); it != mi_vec_hit->end(); it++) {
86  const SQHit* hit = *it;
87  //if (! hit->is_in_time()) continue;
88  int det_id = hit->get_detector_id();
89  list_cnt[det_id]++;
90  }
91 
92  GeomSvc* geom = GeomSvc::instance();
93  mo_evt.n_h1x = list_cnt[geom->getDetectorID("H1T")] + list_cnt[geom->getDetectorID("H1B")];
94  mo_evt.n_h2x = list_cnt[geom->getDetectorID("H2T")] + list_cnt[geom->getDetectorID("H2B")];
95  mo_evt.n_h3x = list_cnt[geom->getDetectorID("H3T")] + list_cnt[geom->getDetectorID("H3B")];
96  mo_evt.n_h4x = list_cnt[geom->getDetectorID("H4T")] + list_cnt[geom->getDetectorID("H4B")];
97  mo_evt.n_d1 = list_cnt[geom->getDetectorID("D0X")];
98  mo_evt.n_d2 = list_cnt[geom->getDetectorID("D2X")];
99  mo_evt.n_d3 = list_cnt[geom->getDetectorID("D3pX")] + list_cnt[geom->getDetectorID("D3mX")];
100 
104  if (mi_vec_trk) {
105  mo_trk_true.clear();
106  for (unsigned int ii = 0; ii < mi_vec_trk->size(); ii++) {
107  SQTrack* trk = mi_vec_trk->at(ii);
108  TrackData td;
109  td.charge = trk->get_charge();
110  td.mom_vtx = trk->get_mom_vtx();
111  mo_trk_true.push_back(td);
112  }
113  }
114  if (mi_vec_dim) {
115  mo_dim_true.clear();
116  for (unsigned int ii = 0; ii < mi_vec_dim->size(); ii++) {
117  SQDimuon* dim = mi_vec_dim->at(ii);
118  DimuonData dd;
119  dd.mom = dim->get_mom();
120  dd.mom_pos = dim->get_mom_pos();
121  dd.mom_neg = dim->get_mom_neg();
122  UtilDimuon::CalcVar(dim, dd.mass, dd.pT, dd.x1, dd.x2, dd.xF);
124  mo_dim_true.push_back(dd);
125  }
126  }
127 
131  if (mi_srec) {
132  mo_trk_reco.clear();
133  for (int ii = 0; ii < mi_srec->getNTracks(); ii++) {
134  SRecTrack* trk = &mi_srec->getTrack(ii);
135  TrackData td;
136  td.charge = trk->getCharge();
137  td.mom_vtx = trk->getMomentumVertex();
138  mo_trk_reco.push_back(td);
139  }
140  mo_dim_reco.clear();
141  for (int ii = 0; ii < mi_srec->getNDimuons(); ii++) {
142  SRecDimuon* dim = &mi_srec->getDimuon(ii);
143  DimuonData dd;
144  dd.mom = dim->p_pos + dim->p_neg;
145  dd.mom_pos = dim->p_pos;
146  dd.mom_neg = dim->p_neg;
147  UtilDimuon::CalcVar(dd.mom_pos, dd.mom_neg, dd.mass, dd.pT, dd.x1, dd.x2, dd.xF);
148  UtilDimuon::Lab2CollinsSoper(dd.mom_pos.Vect(), dd.mom_neg.Vect(), dd.costh_cs, dd.phi_cs);
149  mo_dim_reco.push_back(dd);
150  }
151  }
152 
153  mo_tree->Fill();
155 }
156 
158 {
159  mo_file->cd();
160  mo_file->Write();
161  mo_file->Close();
163 }
164 
166 
178 void AnaEmbeddedData::DivideHitVector(const SQHitVector* vec_in, std::map<int, SQHitVector*>& vec_sep, const bool in_time)
179 {
180  SQHitVector* vec0 = vec_in->Clone();
181  vec0->clear();
182  for (SQHitVector::ConstIter it = mi_vec_hit->begin(); it != mi_vec_hit->end(); it++) {
183  const SQHit* hit = *it;
184  if (in_time && ! hit->is_in_time()) continue;
185  int det_id = hit->get_detector_id();
186  SQHitVector* hv = vec_sep[det_id];
187  if (hv == 0) vec_sep[det_id] = hv = vec0->Clone();
188  hv->push_back(hit);
189  }
190  delete vec0;
191 }
int Init(PHCompositeNode *topNode)
AnaEmbeddedData(const std::string name="AnaEmbeddedData")
int End(PHCompositeNode *topNode)
Called at the end of all processing.
int process_event(PHCompositeNode *topNode)
int InitRun(PHCompositeNode *topNode)
User interface class about the geometry of detector planes.
Definition: GeomSvc.h:164
int getDetectorID(const std::string &detectorName) const
Get the plane position.
Definition: GeomSvc.h:219
static GeomSvc * instance()
singlton instance
Definition: GeomSvc.cxx:212
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 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 int get_qie_rf_intensity(const short i) const =0
Return the i-th QIE RF intensity, where i=-16...+16.
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_spill_id() const =0
Return the spill ID.
virtual int get_event_id() const =0
Return the event ID, which is unique per run.
An SQ interface class to hold a list of SQHit objects.
Definition: SQHitVector.h:32
virtual SQHitVector * Clone() const =0
std::vector< SQHit * >::const_iterator ConstIter
Definition: SQHitVector.h:37
virtual ConstIter end() const =0
virtual ConstIter begin() const =0
virtual void push_back(const SQHit *hit)=0
virtual void clear()=0
An SQ interface class to hold one detector hit.
Definition: SQHit.h:20
virtual bool is_in_time() const
Return 'true' if this hit is in the time window.
Definition: SQHit.h:90
virtual short get_detector_id() const
Return the detector ID of this hit.
Definition: SQHit.h:42
virtual double get_weight() const =0
Return the event weight.
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_charge() const =0
Return the charge, i.e. +1 or -1.
virtual TLorentzVector get_mom_vtx() const =0
Return the track momentum at vertex.
TLorentzVector p_neg
Definition: SRecEvent.h:366
TLorentzVector p_pos
4-momentum of the muon tracks after vertex fit
Definition: SRecEvent.h:365
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
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
void Lab2CollinsSoper(const SQDimuon *dim, double &costh, double &phi)
Convert the momenta of muon pair from Lab frame to Collins-Soper frame.
Definition: UtilDimuon.cc:125
#define PHWHERE
Definition: phool.h:23
TLorentzVector mom_pos
Definition: TreeData.h:37
TLorentzVector mom
Definition: TreeData.h:36
double costh_cs
Definition: TreeData.h:49
double x1
Definition: TreeData.h:41
double phi_cs
Definition: TreeData.h:50
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 weight
Definition: TreeData.h:10
int n_h3x
Definition: TreeData.h:16
int trig_bits
Definition: TreeData.h:11
int n_d2
Definition: TreeData.h:19
int rec_stat
Definition: TreeData.h:12
int n_d1
Definition: TreeData.h:18
int event_id
Definition: TreeData.h:10
int n_h2x
Definition: TreeData.h:15
int rfp00
Definition: TreeData.h:10
int n_h1x
Definition: TreeData.h:14
int n_d3
Definition: TreeData.h:20
int n_h4x
Definition: TreeData.h:17
int rfp01
Definition: TreeData.h:9
int rfm01
Definition: TreeData.h:11
int job_id
Definition: TreeData.h:6
int charge
Definition: TreeData.h:23
TLorentzVector mom_vtx
Definition: TreeData.h:25