Class Reference for E1039 Core & Analysis Software
Fun4AllVectEventOutputManager.cc
Go to the documentation of this file.
1 #include <cstdlib>
2 #include <string>
3 #include <iostream>
4 #include <iomanip>
5 #include <TSystem.h>
6 #include <TFile.h>
7 #include <TTree.h>
8 #include <phool/phool.h>
9 #include <phool/getClass.h>
10 #include <phool/PHNode.h>
11 #include <phool/PHNodeIOManager.h>
12 #include <interface_main/SQEvent.h>
16 #include <ktracker/SRecEvent.h>
19 using namespace std;
20 
22  : Fun4AllOutputManager(myname),
23  m_file(0),
24  m_tree(0),
25  m_tree_name("tree"),
26  m_file_name("output.root"),
27  m_evt(0),
28  m_sp_map(0),
29  m_hit_vec(0)
30 {
31  ;
32 }
33 
35  CloseFile();
36 }
37 
39 std::cout << "Fun4AllVectEventOutputManager::OpenFile(): Attempting to open file: " << m_file_name << " with tree: " << m_tree_name << std::endl;
40 m_file = new TFile(m_file_name.c_str(), "RECREATE");
41 
42 if (!m_file || m_file->IsZombie()) {
43  std::cerr << "Error: Could not create file " << m_file_name << std::endl;
44  exit(1);
45 } else {
46  std::cout << "File " << m_file->GetName() << " opened successfully." << std::endl;
47 }
48 
49 m_tree = new TTree(m_tree_name.c_str(), "Tree for storing events");
50 if (!m_tree) {
51  std::cerr << "Error: Could not create tree " << m_tree_name << std::endl;
52  exit(1);
53 } else {
54  std::cout << "Tree " << m_tree->GetName() << " created successfully." << std::endl;
55 }
56 
57 m_tree->Branch("runID", &runID, "runID/I");
58 m_tree->Branch("spillID", &spillID, "spillID/I");
59 m_tree->Branch("eventID", &eventID, "eventID/I");
60 m_tree->Branch("rfID", &rfID, "rfID/I");
61 m_tree->Branch("turnID", &turnID, "turnID/I");
62 m_tree->Branch("rfIntensities", rfIntensities, "rfIntensities[33]/I");
63 m_tree->Branch("fpgaTriggers", fpgaTriggers, "fpgaTriggers[5]/I");
64 m_tree->Branch("nimTriggers", nimTriggers, "nimTriggers[5]/I");
65 
66 m_tree->Branch("detectorIDs", &detectorIDs);
67 m_tree->Branch("elementIDs", &elementIDs);
68 m_tree->Branch("tdcTimes", &tdcTimes);
69 m_tree->Branch("driftDistances", &driftDistances);
70 m_tree->Branch("hitsInTime", &hitsInTime);
71 
72 m_tree->Branch("triggerDetectorIDs", &triggerDetectorIDs);
73 m_tree->Branch("triggerElementIDs", &triggerElementIDs);
74 m_tree->Branch("triggerTdcTimes", &triggerTdcTimes);
75 m_tree->Branch("triggerDriftDistances", &triggerDriftDistances);
76 m_tree->Branch("triggerHitsInTime", &triggerHitsInTime);
77 
78 
79  m_evt = findNode::getClass<SQEvent>(startNode, "SQEvent");
80  m_hit_vec = findNode::getClass<SQHitVector>(startNode, "SQHitVector");
81  m_trig_hit_vec = findNode::getClass<SQHitVector>(startNode, "SQTriggerHitVector");
82 
83  if (!m_evt || !m_hit_vec || !m_trig_hit_vec) return Fun4AllReturnCodes::ABORTEVENT;
85 }
87  if (!m_file || !m_tree) {
88  OpenFile(startNode);
89  }
90  ResetBranches();
91  runID = m_evt->get_run_id();
92  spillID = m_evt->get_spill_id();
93  rfID = m_evt->get_qie_rf_id();
94  eventID = m_evt->get_event_id();
95  turnID = m_evt->get_qie_turn_id();
96 
97  fpgaTriggers[0] = m_evt->get_trigger(SQEvent::MATRIX1);
98  fpgaTriggers[1] = m_evt->get_trigger(SQEvent::MATRIX2);
99  fpgaTriggers[2] = m_evt->get_trigger(SQEvent::MATRIX3);
100  fpgaTriggers[3] = m_evt->get_trigger(SQEvent::MATRIX4);
101  fpgaTriggers[4] = m_evt->get_trigger(SQEvent::MATRIX5);
102 
103  nimTriggers[0] = m_evt->get_trigger(SQEvent::NIM1);
104  nimTriggers[1] = m_evt->get_trigger(SQEvent::NIM2);
105  nimTriggers[2] = m_evt->get_trigger(SQEvent::NIM3);
106  nimTriggers[3] = m_evt->get_trigger(SQEvent::NIM4);
107  nimTriggers[4] = m_evt->get_trigger(SQEvent::NIM5);
108 
109 
110  for (int i = -16; i < 16; ++i) {
111  // cout << "intensity index: i" << i+16 << endl;
112  rfIntensities[i + 16] = m_evt->get_qie_rf_intensity(i);
113 }
114 
115 if (m_hit_vec) {
116  for (int ihit = 0; ihit < m_hit_vec->size(); ++ihit) {
117  SQHit* hit = m_hit_vec->at(ihit);
118  detectorIDs.push_back(hit->get_detector_id());
119  elementIDs.push_back(hit->get_element_id());
120  tdcTimes.push_back(hit->get_tdc_time());
121  driftDistances.push_back(hit->get_drift_distance());
122  // cout << "get drift distance: " << hit->get_drift_distance() << endl;
123  hitsInTime.push_back(hit->is_in_time());
124  }
125 }
126 
127 if (m_trig_hit_vec) {
128  for (int ihit = 0; ihit < m_trig_hit_vec->size(); ++ihit) {
129  SQHit* hit = m_trig_hit_vec->at(ihit);
130  triggerDetectorIDs.push_back(hit->get_detector_id());
131  triggerElementIDs.push_back(hit->get_element_id());
132  triggerTdcTimes.push_back(hit->get_tdc_time());
133  triggerDriftDistances.push_back(hit->get_drift_distance());
134  triggerHitsInTime.push_back(hit->is_in_time());
135  }
136 }
137 
138  m_tree->Fill();
139  return 0;
140 }
141 
143  if (!m_file) return;
144  std::cout << "Fun4AllVectEventOutputManager::CloseFile(): Closing file: " << m_file_name << std::endl;
145  m_file->Write();
146  m_file->Close();
147  delete m_file;
148  m_file = nullptr;
149 }
150 
152  detectorIDs.clear();
153  elementIDs.clear();
154  tdcTimes.clear();
155  driftDistances.clear();
156  hitsInTime.clear();
157 
158  triggerDetectorIDs.clear();
159  triggerElementIDs.clear();
160  triggerTdcTimes.clear();
161  triggerDriftDistances.clear();
162  triggerHitsInTime.clear();
163 }
164 
virtual int Write(PHCompositeNode *startNode)
write starting from given node
Fun4AllVectEventOutputManager(const std::string &myname="UNIVERSALOUT")
int OpenFile(PHCompositeNode *startNode)
virtual int get_run_id() const =0
Return the run ID.
@ MATRIX2
Definition: SQEvent.h:28
@ NIM5
Definition: SQEvent.h:26
@ MATRIX1
Definition: SQEvent.h:27
@ NIM4
Definition: SQEvent.h:25
@ NIM2
Definition: SQEvent.h:23
@ NIM1
Definition: SQEvent.h:22
@ MATRIX3
Definition: SQEvent.h:29
@ NIM3
Definition: SQEvent.h:24
@ MATRIX4
Definition: SQEvent.h:30
@ MATRIX5
Definition: SQEvent.h:31
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_qie_turn_id() const =0
Return the QIE turn ID.
virtual int get_qie_rf_id() const =0
Return the QIE RF ID.
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.
virtual const SQHit * at(const size_t idkey) const =0
virtual size_t size() const =0
An SQ interface class to hold one detector hit.
Definition: SQHit.h:20
virtual float get_drift_distance() const
Return the drift distance of this hit. Probably the value is not properly set at present....
Definition: SQHit.h:57
virtual bool is_in_time() const
Return 'true' if this hit is in the time window.
Definition: SQHit.h:90
virtual short get_element_id() const
Return the element ID of this hit.
Definition: SQHit.h:45
virtual float get_tdc_time() const
Return the TDC time (nsec) of this hit.
Definition: SQHit.h:54
virtual short get_detector_id() const
Return the detector ID of this hit.
Definition: SQHit.h:42