Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MvtxCombinedRawDataDecoder.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file MvtxCombinedRawDataDecoder.cc
1 
8 
9 #include <trackbase/MvtxDefs.h>
11 #include <trackbase/TrkrHitSet.h>
13 #include <trackbase/TrkrHitv2.h>
14 
19 
21 
22 #include <phool/PHCompositeNode.h>
23 #include <phool/PHNodeIterator.h>
24 #include <phool/getClass.h>
25 
26 #include <cdbobjects/CDBTTree.h>
27 #include <ffamodules/CDBInterface.h> // for accessing the MVTX hot pixel file from the CDB
28 #include <algorithm>
29 #include <cassert>
30 
31 //_________________________________________________________
33  : SubsysReco(name)
34 {
35 }
36 
37 //_____________________________________________________________________
39 {
41 }
42 
43 //____________________________________________________________________________..
45 {
46  // get dst node
47  PHNodeIterator iter(topNode);
48  PHCompositeNode *dstNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
49  if (!dstNode)
50  {
51  std::cout << "MvtxCombinedRawDataDecoder::InitRun - DST Node missing, doing nothing." << std::endl;
52  exit(1);
53  }
54 
55  // create hitset container if needed
56  hit_set_container = findNode::getClass<TrkrHitSetContainer>(topNode, "TRKR_HITSET");
57  if (!hit_set_container)
58  {
59  // find or create TRKR node
60  auto trkrNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "TRKR"));
61  if (!trkrNode)
62  {
63  trkrNode = new PHCompositeNode("TRKR");
64  dstNode->addNode(trkrNode);
65  }
66 
67  // create container and add to the tree
69  auto newNode = new PHIODataNode<PHObject>(hit_set_container, "TRKR_HITSET", "PHObject");
70  trkrNode->addNode(newNode);
71  }
72 
73  // Check if MVTX event header already exists
75  {
76  auto mvtxNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "MVTX"));
77  if (!mvtxNode)
78  {
79  mvtxNode = new PHCompositeNode("MVTX");
80  dstNode->addNode(mvtxNode);
81  }
82 
83  mvtx_event_header = findNode::getClass<MvtxEventInfo>(mvtxNode, "MVTXEVENTHEADER");
84  if (!mvtx_event_header)
85  {
87  auto newHeader = new PHIODataNode<PHObject>(mvtx_event_header, "MVTXEVENTHEADER", "PHObject");
88  mvtxNode->addNode(newHeader);
89  }
90  }
91 
92  mvtx_raw_event_header = findNode::getClass<MvtxRawEvtHeader>(topNode, m_MvtxRawEvtHeaderNodeName);
94  {
95  std::cout << PHWHERE << "::" << __func__ << ": Could not get \"" << m_MvtxRawEvtHeaderNodeName << "\" from Node Tree" << std::endl;
96  std::cout << "Have you built this yet?" << std::endl;
97  exit(1);
98  }
99 
100  // Mask Hot MVTX Pixels
101  std::string database = CDBInterface::instance()->getUrl("MVTX_HotPixelMap"); // This is specifically for MVTX Hot Pixels
102  CDBTTree *cdbttree = new CDBTTree(database);
103  int NPixel = -1;
104  NPixel = cdbttree->GetSingleIntValue("TotalHotPixels");
105 
106  for (int i = 0; i < NPixel; i++)
107  {
108  int Layer = cdbttree->GetIntValue(i, "layer");
109  int Stave = cdbttree->GetIntValue(i, "stave");
110  int Chip = cdbttree->GetIntValue(i, "chip");
111  int Col = cdbttree->GetIntValue(i, "col");
112  int Row = cdbttree->GetIntValue(i, "row");
113 
114  TrkrDefs::hitsetkey HotPixelHitKey = MvtxDefs::genHitSetKey(Layer, Stave, Chip, 0);
115  TrkrDefs::hitkey HotHitKey = MvtxDefs::genHitKey(Col, Row);
116  m_hotPixelMap.push_back({std::make_pair(HotPixelHitKey, HotHitKey)});
117  }
118 
120 }
121 
122 //___________________________________________________________________________
124 {
125  mvtx_raw_event_header = findNode::getClass<MvtxRawEvtHeader>(topNode, m_MvtxRawEvtHeaderNodeName);
127 
128  hit_set_container = findNode::getClass<TrkrHitSetContainer>(topNode, "TRKR_HITSET");
129 
130  mvtx_hit_container = findNode::getClass<MvtxRawHitContainer>(topNode, m_MvtxRawHitNodeName);
131 
132  if (!mvtx_hit_container)
133  {
134  std::cout << PHWHERE << "::" << __func__ << ": Could not get \"" << m_MvtxRawHitNodeName << "\" from Node Tree" << std::endl;
135  std::cout << "Have you built this yet?" << std::endl;
136  exit(1);
137  }
138  auto gl1 = findNode::getClass<Gl1RawHit>(topNode, "GL1RAWHIT");
139  if (!gl1)
140  {
141  std::cout << PHWHERE << "Could not get gl1 raw hit" << std::endl;
143  }
144  uint64_t gl1rawhitbco = gl1->get_bco();
145  // get the last 40 bits by bit shifting left then right to match
146  // to the mvtx bco
147  auto lbshift = gl1rawhitbco << 24;
148  auto gl1bco = lbshift >> 24;
149 
150  if (Verbosity() >= VERBOSITY_MORE) mvtx_hit_container->identify();
151 
152  uint64_t strobe = -1; // Initialise to -1 for debugging
153  uint8_t layer = 0;
154  uint8_t stave = 0;
155  uint8_t chip = 0;
156  uint16_t row = 0;
157  uint16_t col = 0;
158  std::vector<std::pair<uint64_t, uint32_t>> strobe_bc_pairs;
159  std::set<uint64_t> l1BCOs = mvtx_raw_event_header->getMvtxLvL1BCO();
160  auto mvtxbco = *l1BCOs.begin();
161  if (Verbosity() > 0)
162  {
163  std::cout << "mvtx header bco " << mvtxbco << " and gl1 bco " << gl1bco << std::endl;
164  }
165 
167  {
168  mvtx_event_header = findNode::getClass<MvtxEventInfo>(topNode, "MVTXEVENTHEADER");
170  }
171 
172  // int NMasked = 0;
173 
174  for (unsigned int i = 0; i < mvtx_hit_container->get_nhits(); i++)
175  {
176  mvtx_hit = mvtx_hit_container->get_hit(i);
177  strobe = mvtx_hit->get_bco();
178  layer = mvtx_hit->get_layer_id();
179  stave = mvtx_hit->get_stave_id();
180  chip = mvtx_hit->get_chip_id();
181  row = mvtx_hit->get_row();
182  col = mvtx_hit->get_col();
183 
184  uint64_t bcodiff = gl1bco - strobe;
185  double timeElapsed = bcodiff * 0.106; // 106 ns rhic clock
186  int index = std::floor(timeElapsed / m_strobeWidth);
187 
189 
190  const TrkrDefs::hitsetkey hitsetkey = MvtxDefs::genHitSetKey(layer, stave, chip, index);
191  if (!hitsetkey) continue;
192 
193  // get matching hitset
194  const auto hitset_it = hit_set_container->findOrAddHitSet(hitsetkey);
195 
196  // generate hit key
197  const TrkrDefs::hitkey hitkey = MvtxDefs::genHitKey(col, row);
198 
199  // find existing hit, or create
200  auto hit = hitset_it->second->getHit(hitkey);
201  if (hit)
202  {
203  std::cout << PHWHERE << "::" << __func__ << " - duplicated hit, hitsetkey: " << hitsetkey << " hitkey: " << hitkey << std::endl;
204  continue;
205  }
206 
207  const TrkrDefs::hitsetkey hitsetkeymask = MvtxDefs::genHitSetKey(layer, stave, chip, 0);
208 
209  if (std::find(m_hotPixelMap.begin(), m_hotPixelMap.end(), std::make_pair(hitsetkeymask, hitkey)) == m_hotPixelMap.end())
210  {
211  // create hit and insert in hitset
212  hit = new TrkrHitv2;
213  hitset_it->second->addHitSpecificKey(hitkey, hit);
214  }
215  }
216 
219  {
220  std::set<uint64_t> l1BCOs = mvtx_raw_event_header->getMvtxLvL1BCO();
221  for (auto iter = l1BCOs.begin(); iter != l1BCOs.end(); iter++)
222  {
224  }
226  }
227 
229 }
230 
231 //_____________________________________________________________________
233 {
235 }
236 
237 void MvtxCombinedRawDataDecoder::removeDuplicates(std::vector<std::pair<uint64_t, uint32_t>> &v)
238 {
239  auto end = v.end();
240  for (auto it = v.begin(); it != end; ++it)
241  {
242  end = remove(it + 1, end, *it);
243  }
244  v.erase(end, v.end());
245 }