Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
InttCombinedRawDataDecoder.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file InttCombinedRawDataDecoder.cc
2 #include "InttMapping.h"
3 
4 #include <trackbase/InttDefs.h>
5 #include <trackbase/TrkrDefs.h> // for hitkey, hitsetkey
6 #include <trackbase/TrkrHit.h>
7 #include <trackbase/TrkrHitSet.h>
10 #include <trackbase/TrkrHitv2.h>
11 
15 
17 
18 #include <phool/PHCompositeNode.h>
19 #include <phool/PHIODataNode.h> // for PHIODataNode
20 #include <phool/PHNodeIterator.h>
21 #include <phool/getClass.h>
22 #include <phool/phool.h> // for PHWHERE
23 
24 #include <TSystem.h>
25 
26 #include <cstdlib> // for exit
27 #include <iostream> // for operator<<, endl, bas...
28 #include <map> // for _Rb_tree_iterator
29 
31  : SubsysReco(name)
32 {
33  // Do nothing
34  // Consider calling LoadHotChannelMapRemote()
35 }
36 
38 {
39  if (!topNode)
40  {
41  std::cout << "InttCombinedRawDataDecoder::InitRun(PHCompositeNode* topNode)" << std::endl;
42  std::cout << "\tCould not retrieve topNode; doing nothing" << std::endl;
43  exit(1);
44  gSystem->Exit(1);
45 
46  return 1;
47  }
48 
49  PHNodeIterator dst_itr(topNode);
50  PHCompositeNode* dst_node = dynamic_cast<PHCompositeNode*>(dst_itr.findFirst("PHCompositeNode", "DST"));
51  if (!dst_node)
52  {
53  if (Verbosity())
54  {
55  std::cout << "InttCombinedRawDataDecoder::InitRun(PHCompositeNode* topNode)" << std::endl;
56  }
57  if (Verbosity())
58  {
59  std::cout << "\tCould not retrieve dst_node; doing nothing" << std::endl;
60  }
61  exit(1);
62  gSystem->Exit(1);
63 
64  return 1;
65  }
66 
67  PHNodeIterator trkr_itr(dst_node);
68  PHCompositeNode* trkr_node = dynamic_cast<PHCompositeNode*>(trkr_itr.findFirst("PHCompositeNode", "TRKR"));
69  if (!trkr_node)
70  {
71  trkr_node = new PHCompositeNode("TRKR");
72  dst_node->addNode(trkr_node);
73  }
74 
75  TrkrHitSetContainer* trkr_hit_set_container = findNode::getClass<TrkrHitSetContainer>(topNode, "TRKR_HITSET");
76  if (!trkr_hit_set_container)
77  {
78  if (Verbosity())
79  {
80  std::cout << "InttCombinedRawDataDecoder::InitRun(PHCompositeNode* topNode)" << std::endl;
81  }
82  if (Verbosity())
83  {
84  std::cout << "\tMaking TrkrHitSetContainer" << std::endl;
85  }
86 
87  trkr_hit_set_container = new TrkrHitSetContainerv1;
88  PHIODataNode<PHObject>* new_node = new PHIODataNode<PHObject>(trkr_hit_set_container, "TRKR_HITSET", "PHObject");
89  trkr_node->addNode(new_node);
90  }
91 
93 }
94 
96 {
97  TrkrHitSetContainer* trkr_hit_set_container = findNode::getClass<TrkrHitSetContainer>(topNode, "TRKR_HITSET");
98  if (!trkr_hit_set_container)
99  {
100  std::cout << PHWHERE << std::endl;
101  std::cout << "InttCombinedRawDataDecoder::process_event(PHCompositeNode* topNode)" << std::endl;
102  std::cout << "Could not get \"TRKR_HITSET\" from Node Tree" << std::endl;
103  std::cout << "Exiting" << std::endl;
104  gSystem->Exit(1);
105  exit(1);
106 
108  }
109 
110  InttRawHitContainer* inttcont = findNode::getClass<InttRawHitContainer>(topNode, m_InttRawNodeName);
111  if (!inttcont)
112  {
113  std::cout << PHWHERE << std::endl;
114  std::cout << "InttCombinedRawDataDecoder::process_event(PHCompositeNode* topNode)" << std::endl;
115  std::cout << "Could not get \"" << m_InttRawNodeName << "\" from Node Tree" << std::endl;
116  std::cout << "Exiting" << std::endl;
117 
118  gSystem->Exit(1);
119  exit(1);
120  }
121  auto gl1 = findNode::getClass<Gl1RawHit>(topNode, "GL1RAWHIT");
122  if (!gl1)
123  {
124  std::cout << PHWHERE << " no gl1 container, exiting" << std::endl;
126  }
127 
128  uint64_t gl1rawhitbco = gl1->get_bco();
129  // get the last 40 bits by bit shifting left then right to match
130  // to the mvtx bco
131  auto lbshift = gl1rawhitbco << 24U; // clang-tidy: mark as unsigned
132  auto gl1bco = lbshift >> 24U; // clang-tidy: mark as unsigned
133 
134  TrkrDefs::hitsetkey hit_set_key = 0;
135  TrkrDefs::hitkey hit_key = 0;
136  TrkrHitSetContainer::Iterator hit_set_container_itr;
137  TrkrHit* hit = nullptr;
138 
141  for (unsigned int i = 0; i < inttcont->get_nhits(); i++)
142  {
143  InttRawHit* intthit = inttcont->get_hit(i);
144  // uint64_t gtm_bco = intthit->get_bco();
145 
146  InttNameSpace::RawFromHit(raw, intthit);
147  // raw.felix_server = InttNameSpace::FelixFromPacket(intthit->get_packetid());
148  // raw.felix_channel = intthit->get_fee();
149  // raw.chip = (intthit->get_chip_id() + 25) % 26;
150  // raw.channel = intthit->get_channel_id();
151 
152  int adc = intthit->get_adc();
153  // amp = intthit->get_amplitude();
154  // int bco = intthit->get_FPHX_BCO();
155 
156  if (m_HotChannelSet.find(raw) != m_HotChannelSet.end())
157  {
158  continue;
159  }
160 
161  ofl = InttNameSpace::ToOffline(raw);
162  hit_key = InttDefs::genHitKey(ofl.strip_y, ofl.strip_x); // col, row <trackbase/InttDefs.h>
163  hit_set_key = InttDefs::genHitSetKey(ofl.layer, ofl.ladder_z, ofl.ladder_phi, intthit->get_bco() - gl1bco);
164  hit_set_container_itr = trkr_hit_set_container->findOrAddHitSet(hit_set_key);
165  hit = hit_set_container_itr->second->getHit(hit_key);
166 
167  if (hit)
168  {
169  continue;
170  }
171 
172  hit = new TrkrHitv2;
173  hit->setAdc(adc);
174  hit_set_container_itr->second->addHitSpecificKey(hit_key, hit);
175  }
176 
178 }
179 
181 {
182  if (filename.empty())
183  {
184  std::cout << "int InttCombinedRawDataDecoder::LoadHotChannelMapLocal(std::string const& filename)" << std::endl;
185  std::cout << "\tArgument 'filename' is empty string" << std::endl;
186  return 1;
187  }
188  CDBTTree cdbttree(filename);
189  // need to checkt for error exception
190  cdbttree.LoadCalibrations();
191 
192  m_HotChannelSet.clear();
193  Long64_t N = cdbttree.GetSingleIntValue("size");
194  for (Long64_t n = 0; n < N; ++n)
195  {
197  .felix_server = cdbttree.GetIntValue(n, "felix_server"),
198  .felix_channel = cdbttree.GetIntValue(n, "felix_channel"),
199  .chip = cdbttree.GetIntValue(n, "chip"),
200  .channel = cdbttree.GetIntValue(n, "channel")});
201 
202  // if(Verbosity() < 1)
203  // {
204  // continue;
205  // }
206  // std::cout << "Masking channel:\n" << std::endl;
207  // std::cout << "\t" << cdbttree.GetIntValue(n, "felix_server")
208  // << "\t" << cdbttree.GetIntValue(n, "felix_channel")
209  // << "\t" << cdbttree.GetIntValue(n, "chip")
210  // << "\t" << cdbttree.GetIntValue(n, "channel") << std::endl;
211  }
212 
213  return 0;
214 }
215 
217 {
218  if (name.empty())
219  {
220  std::cout << "int InttCombinedRawDataDecoder::LoadHotChannelMapRemote(std::string const& name)" << std::endl;
221  std::cout << "\tArgument 'name' is empty string" << std::endl;
222  return 1;
223  }
224  std::string database = CDBInterface::instance()->getUrl(name);
225  CDBTTree cdbttree(database);
226  cdbttree.LoadCalibrations();
227 
228  m_HotChannelSet.clear();
229  Long64_t N = cdbttree.GetSingleIntValue("size");
230  for (Long64_t n = 0; n < N; ++n)
231  {
233  .felix_server = cdbttree.GetIntValue(n, "felix_server"),
234  .felix_channel = cdbttree.GetIntValue(n, "felix_channel"),
235  .chip = cdbttree.GetIntValue(n, "chip"),
236  .channel = cdbttree.GetIntValue(n, "channel")});
237  }
238 
239  return 0;
240 }
241 
242 /*
243  Packet* p = evt->getPacket(itr->first);
244  if(!p)continue;
245 
246  int N = p->iValue(0, "NR_HITS");
247  full_bco = p->lValue(0, "BCO");
248 
249  if(Verbosity() > 20)std::cout << N << std::endl;
250 
251  for(int n = 0; n < N; ++n)
252  {
253  rawdata = InttNameSpace::RawFromPacket(itr->second, n, p);
254 
255  adc = p->iValue(n, "ADC");
256  //amp = p->iValue(n, "AMPLITUE");
257  bco = p->iValue(n, "FPHX_BCO");
258 
259  offline = InttNameSpace::ToOffline(rawdata);
260 
261  hit_key = InttDefs::genHitKey(offline.strip_y, offline.strip_x); //col, row <trackbase/InttDefs.h>
262  hit_set_key = InttNameSpace::genHitSetKey(offline.layer, offline.ladder_z, offline.ladder_phi, bco);
263 
264  hit_set_container_itr = trkr_hit_set_container->findOrAddHitSet(hit_set_key);
265  hit = hit_set_container_itr->second->getHit(hit_key);
266  if(hit)continue;
267 
268  hit = new TrkrHitv2;
269  hit->setAdc(adc);
270  hit_set_container_itr->second->addHitSpecificKey(hit_key, hit);
271  }
272 
273  delete p;
274  }
275  }
276  if(Verbosity() > 20)
277  {
278  std::cout << std::endl;
279  std::cout << "Identify():" << std::endl;
280  trkr_hit_set_container->identify();
281  std::cout << std::endl;
282  }
283 */