Class Reference for E1039 Core & Analysis Software
SubsysRecoBG.cc
Go to the documentation of this file.
1 #include <iomanip>
2 #include <fstream>
3 #include <TFile.h>
4 #include <TTree.h>
5 #include <interface_main/SQRun.h>
10 #include <phool/getClass.h>
11 #include <UtilAna/UtilSQHit.h>
12 #include "SubsysRecoBG.h"
13 using namespace std;
14 
15 SubsysRecoBG::SubsysRecoBG(const std::string &name)
16  : SubsysReco(name)
17  , m_mode(DEFAULT)
18  , mi_evt(0)
19  , mi_mc_evt(0)
20  , mi_vec_hit(0)
21  , mo_file(0)
22  , mo_tree(0)
23 {
24  ;
25 }
26 
28 {
30 }
31 
33 {
34  mi_evt = findNode::getClass<SQEvent >(topNode, "SQEvent");
35  mi_vec_hit = findNode::getClass<SQHitVector>(topNode, "SQHitVector");
36  if (!mi_evt || !mi_vec_hit) return Fun4AllReturnCodes::ABORTEVENT;
37 
38  if (m_mode == FULL_BG) {
39  mi_mc_evt = findNode::getClass<SQMCEvent>(topNode, "SQMCEvent");
40  if (!mi_mc_evt) return Fun4AllReturnCodes::ABORTEVENT;
41  }
42 
43  mo_file = new TFile("bg_data.root", "RECREATE");
44  mo_tree = new TTree("bg_tree", "Created by SubsysRecoBG");
45  mo_tree->Branch("bg_data", &mo_bg);
46 
48 }
49 
51 {
52  mo_bg.run = mi_evt->get_run_id();
53  mo_bg.evt = mi_evt->get_event_id();
54  mo_bg.fpga1 = mi_evt->get_trigger(SQEvent::MATRIX1);
55 
56  if (m_mode == FULL_BG) {
57  mo_bg.pot_rfp00 = mi_mc_evt->get_cross_section();
58  }
59  mo_bg.inte_rfp00 = mi_evt->get_qie_rf_intensity(0);
60  mo_bg.inte_max = 0;
61  for (int ii = -8; ii <= 8; ii++) {
62  int inte = mi_evt->get_qie_rf_intensity(ii);
63  if (inte > mo_bg.inte_max) mo_bg.inte_max = inte;
64  }
65 
66  ExtractHits(mi_vec_hit, "H1T", mo_bg.h1t);
67  ExtractHits(mi_vec_hit, "H2T", mo_bg.h2t);
68  ExtractHits(mi_vec_hit, "H3T", mo_bg.h3t);
69  ExtractHits(mi_vec_hit, "H4T", mo_bg.h4t);
70  ExtractHits(mi_vec_hit, "H1B", mo_bg.h1b);
71  ExtractHits(mi_vec_hit, "H2B", mo_bg.h2b);
72  ExtractHits(mi_vec_hit, "H3B", mo_bg.h3b);
73  ExtractHits(mi_vec_hit, "H4B", mo_bg.h4b);
74 
75  mo_tree->Fill();
76 
78 }
79 
81 {
82  mo_file->cd();
83  mo_file->Write();
84  mo_file->Close();
86 }
87 
88 void SubsysRecoBG::ExtractHits(const SQHitVector* hit_vec, const std::string det_name, std::vector<int>& list_ele)
89 {
90  list_ele.clear();
91  shared_ptr<SQHitVector> hv(UtilSQHit::FindFirstHits(hit_vec, det_name, true));
92  for (SQHitVector::Iter it = hv->begin(); it != hv->end(); it++) {
93  list_ele.push_back( (*it)->get_element_id() );
94  }
95 }
virtual int get_run_id() const =0
Return the run ID.
@ MATRIX1
Definition: SQEvent.h:27
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_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
std::vector< SQHit * >::iterator Iter
Definition: SQHitVector.h:38
virtual double get_cross_section() const =0
Return the cross section.
int Init(PHCompositeNode *topNode)
Definition: SubsysRecoBG.cc:27
int End(PHCompositeNode *topNode)
Called at the end of all processing.
Definition: SubsysRecoBG.cc:80
int InitRun(PHCompositeNode *topNode)
Definition: SubsysRecoBG.cc:32
void ExtractHits(const SQHitVector *hit_vec, const std::string det_name, std::vector< int > &list_ele)
Definition: SubsysRecoBG.cc:88
int process_event(PHCompositeNode *topNode)
Definition: SubsysRecoBG.cc:50
SubsysRecoBG(const std::string &name="SubsysRecoBG")
Definition: SubsysRecoBG.cc:15
SQHitVector * FindFirstHits(const SQHitVector *vec_in, const std::string det_name, const bool in_time=false)
Extract a set of first hits that are of the given detector (det_name), where "first" means the earlie...
Definition: UtilSQHit.cc:40
std::vector< int > h1t
Definition: TreeData.h:33
int run
Definition: TreeData.h:27
int inte_max
Definition: TreeData.h:32
std::vector< int > h1b
Definition: TreeData.h:37
int inte_rfp00
In unit of QIE count.
Definition: TreeData.h:31
int evt
Definition: TreeData.h:28
double pot_rfp00
In unit of N of protons.
Definition: TreeData.h:30
std::vector< int > h2t
Definition: TreeData.h:34
std::vector< int > h4b
Definition: TreeData.h:40
bool fpga1
Definition: TreeData.h:29
std::vector< int > h2b
Definition: TreeData.h:38
std::vector< int > h4t
Definition: TreeData.h:36
std::vector< int > h3t
Definition: TreeData.h:35
std::vector< int > h3b
Definition: TreeData.h:39