ATLAS Offline Software
Loading...
Searching...
No Matches
FPGATrackSimMapMakerAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
7
8#include "TH2.h"
9#include "TFile.h"
10
11#include "GaudiKernel/IEventProcessor.h"
12
13
15// Initialize
16
17FPGATrackSimMapMakerAlg::FPGATrackSimMapMakerAlg (const std::string& name, ISvcLocator* pSvcLocator) :
18 AthAlgorithm(name, pSvcLocator)
19{
20}
21
22
24{
25
26 m_monitorFile.reset (new TFile((m_outFileName.value() + "region" + std::to_string(m_region) + ".root").c_str(), "RECREATE"));
27 std::stringstream ss(m_description);
28 std::string line;
29 ATH_MSG_INFO("Tag config:");
30 if (m_description.value() != "")
31 while (std::getline(ss, line, '\n'))
32 ATH_MSG_INFO('\t' << line);
33
34 // reset Hit and Module vectors
35 m_pbHits.clear(); m_peHits.clear(); m_sbHits.clear(); m_seHits.clear(); m_allHits.clear();
36 m_track2modules.clear(); m_track2slice.clear(); m_slice2modules.clear();
37 m_key_etamods.clear(); m_key_etamods2.clear(); m_usedTracks.clear();
38 m_radii.clear(); m_keylayer.clear(); m_keylayer2.clear();
39 m_modules.clear();
40
41 // Set up the module relabel tool.
43
44 // We need to select this before calling parseKeyString() in case the user
45 // has specified a plane (logical layer) as the keystring.
46 // If running with the new regions. Don't use the hardcoded layer assignments. Use either a generic one
47 // or one passed in via a flag.
48 ATH_MSG_INFO("Using Python configuration for regions. All pixel layers -> first stage, all strip layers -> second stage");
49 m_planes = &(m_overridePlanes.value());
50 m_planes2 = &(m_overridePlanes2.value());
51
53 ATH_CHECK(m_evtSel.retrieve());
54
55 ATH_CHECK(m_hitSGInputTool.retrieve(EnableTool{!m_hitSGInputTool.empty()}));
56 ATH_CHECK(m_hitInputTool.retrieve(EnableTool{!m_hitInputTool.empty()}));
57
58 ATH_MSG_DEBUG("initialize() Instantiating root objects");
59 ATH_MSG_DEBUG("initialize() Finished");
60
61
62 return StatusCode::SUCCESS;
63}
64
66// MAIN EXECUTE ROUTINE //
68
70{
71 // Read inputs
72 bool done = false;
73 ATH_CHECK(readInputs(done));
74
75 if (done) {
76 SmartIF<IEventProcessor> appMgr{service("ApplicationMgr")};
77 if (!appMgr) {
78 ATH_MSG_ERROR("Failed to retrieve ApplicationMgr as IEventProcessor");
79 return StatusCode::FAILURE;
80 }
81 return appMgr->stopRun();
82 }
83
84 // Reset data pointers
85 m_eventHeader.reset();
86
87 return StatusCode::SUCCESS;
88}
89
90
92// INPUT READING AND PROCESSING //
94
95
97{
98 // Read primary input
99 if ( !m_hitSGInputTool.empty()) {
100 ATH_CHECK(m_hitSGInputTool->readData(&m_eventHeader, Gaudi::Hive::currentContext()));
101 ATH_MSG_DEBUG("Loaded " << m_eventHeader.nHits() << " hits in event header from SG");
102 }
103 else {
104 ATH_CHECK(m_hitInputTool->readData(&m_eventHeader, done));
105 if (done)
106 {
107 ATH_MSG_INFO("Cannot read more events from file, returning");
108 return StatusCode::SUCCESS; // end of loop over events
109 }
110 }
111
112 // Ask the event selection service if this event really falls within the region.
113 if (!m_evtSel->selectEvent(&m_eventHeader))
114 {
115 ATH_MSG_DEBUG("Event skipped by: " << m_evtSel->name());
116 return StatusCode::SUCCESS;
117 }
118
119 FPGATrackSimEventInfo eventinfo = m_eventHeader.event();
120 ATH_MSG_DEBUG ("Getting Event " << eventinfo);
121
122 for (auto hit: m_eventHeader.hits()) // fill track to modules map and hit vectors
123 {
124 // barcode cut: these hits are not associated with our single muon tracks
125 if (hit.getBarcodePt() == 0) continue;
126
127 SiliconTech det = hit.getDetType();
128 DetectorZone bec = hit.getDetectorZone();
129 int lyr = hit.getPhysLayer();
130 int eta = hit.getEtaModule();
131 int phi = hit.getPhiModule();
132
133 if (hit.isPixel() && hit.isBarrel()) {
134 m_pbHits.push_back(hit);
135 if (lyr > m_pbmax) m_pbmax = lyr;
136 }
137 else if (hit.isPixel() && !hit.isBarrel()) {
138 // Perform the layer + module number remapping using our tool.
139 m_moduleRelabel->remap(hit);
140 eta = hit.getEtaModule();
141 lyr = hit.getPhysLayer();
142
143 m_peHits.push_back(hit);
144 if (bec == DetectorZone::posEndcap && lyr > m_pemax[0]) m_pemax[0] = lyr;
145 if (bec == DetectorZone::negEndcap && lyr > m_pemax[1]) m_pemax[1] = lyr;
146 }
147 else if (!hit.isPixel() && hit.isBarrel()) {
148 m_sbHits.push_back(hit);
149 if (lyr > m_sbmax) m_sbmax = lyr;
150 }
151 else if (!hit.isPixel() && !hit.isBarrel()) {
152 m_seHits.push_back(hit);
153 if (bec == DetectorZone::posEndcap && lyr > m_semax[0]) m_semax[0] = lyr;
154 if (bec == DetectorZone::negEndcap && lyr > m_semax[1]) m_semax[1] = lyr;
155 }
156 else ATH_MSG_ERROR("Invalid Region");
157
158 // keep track of all etamods on the key layer for slicing
159 if (isOnKeyLayer(1,det,bec,lyr)) {
160 if (m_key_etamods.count(eta) == 0) {
161 m_key_etamods.insert(std::pair<int, int>(eta, 1));
162 } else {
163 m_key_etamods[eta] += 1;
164 }
165 }
166
167 if (m_key2) // if doing 2D slicing
168 if (isOnKeyLayer(2,det,bec,lyr)) m_key_etamods2.insert(eta);
169
170 // Keep a global record of all modules as we see them, indexed by the ID tuple.
171 Module mod = Module(det, bec, lyr, eta, phi);
172 if (m_modules.count(mod.moduleId()) == 0) {
173 m_modules.insert(std::pair<FPGATrackSimModuleId, Module>(mod.moduleId(), mod));
174 }
175 Module* modref = &(m_modules[mod.moduleId()]);
176
177 if (std::find(m_track2modules[hit.getEventIndex()].begin(), m_track2modules[hit.getEventIndex()].end(), modref) == m_track2modules[hit.getEventIndex()].end())
178 m_track2modules[hit.getEventIndex()].push_back(modref);
179 m_allHits.push_back(hit);
180 }
181
182 return StatusCode::SUCCESS;
183}
184
185
186
187
189// OUTPUT WRITING AND MONITORING //
191StatusCode FPGATrackSimMapMakerAlg::writePmapAndRmap(std::vector<FPGATrackSimHit> const & pbHits, std::vector<FPGATrackSimHit> const & peHits, std::vector<FPGATrackSimHit> const & sbHits, std::vector<FPGATrackSimHit> const & seHits, int reg)
192{
193 // Plane Map
194
195 // to avoid vector _M_range_check errors, print at least one line for each DetectorZone
196 if (m_pbmax == -1) m_pbmax = 0;
197 if (m_sbmax == -1) m_sbmax = 0;
198 if (m_pemax[0] == -1) m_pemax[0] = 0;
199 if (m_pemax[1] == -1) m_pemax[1] = 0;
200 if (m_semax[0] == -1) m_semax[0] = 0;
201 if (m_semax[1] == -1) m_semax[1] = 0;
202
203 std::string pmap_path = m_outFileName.value() + "region" + std::to_string(m_region) + ".pmap";
204 ATH_MSG_INFO("Creating pmap: " << pmap_path);
205 m_pmap.open(pmap_path, std::ofstream::out);
206
207 m_pmap << m_geoTag.value() << "\n" << m_planes->size() << " logical_s1\n" << m_planes2->size() << " logical_s2\n";
208 m_pmap << m_pbmax+1 << " pixel barrel \n" << m_pemax[0]+1 << " pixel endcap+ \n" << m_pemax[1]+1 << " pixel endcap- \n";
209 m_pmap << m_sbmax+1 << " SCT barrel \n" << m_semax[0]+1 << " SCT endcap+\n" << m_semax[1]+1 << " SCT endcap-\n";
210 m_pmap << "! silicon endCap physDisk physLayer ['stereo' stripSide <strip only>] 'plane1' logiLayer1 'plane2' logiLayer2\n";
211 m_pmap << "\nregion " << reg << "\n";
212
213 int p1,p2;
214 for (int lyr = 0; lyr <= m_pbmax; lyr++) { // Pixel Barrel
215 p1 = findPlane(m_planes, "pb" + std::to_string(lyr));
216 p2 = findPlane(m_planes2, "pb" + std::to_string(lyr));
217 m_pmap << "pixel 0 -1 " << lyr << " plane1 " << p1 << " plane2 " << p2 << "\n";
218 }
219 for (int lyr = 0; lyr <= m_pemax[0]; lyr++) { // Pixel Postive Endap
220 p1 = findPlane(m_planes, "pe" + std::to_string(lyr) + "+");
221 p2 = findPlane(m_planes2, "pe" + std::to_string(lyr) + "+");
222 m_pmap << "pixel 1 " << lyr << " " << lyr << " plane1 " << p1 << " plane2 " << p2 << "\n";
223 }
224 for (int lyr = 0; lyr <= m_pemax[1]; lyr++) { // Pixel Negative Endcap
225 p1 = findPlane(m_planes, "pe" + std::to_string(lyr) + "-");
226 p2 = findPlane(m_planes2, "pe" + std::to_string(lyr) + "-");
227 m_pmap << "pixel 2 " << lyr << " " << lyr << " plane1 " << p1 << " plane2 " << p2 << "\n";
228 }
229 for (int lyr = 0; lyr <= m_sbmax; lyr++) { // Strip Barrel
230 p1 = findPlane(m_planes, "sb" + std::to_string(lyr));
231 p2 = findPlane(m_planes2, "sb" + std::to_string(lyr));
232 m_pmap << "SCT 0 -1 " << lyr << " stereo " << lyr % 2 << " plane1 " << p1 << " plane2 " << p2 << "\n";
233 }
234 for (int lyr = 0; lyr <= m_semax[0]; lyr++) { // Strip Positive Endcap
235 p1 = findPlane(m_planes, "se" + std::to_string(lyr) + "+");
236 p2 = findPlane(m_planes2, "se" + std::to_string(lyr) + "+");
237 m_pmap << "SCT 1 " << lyr/2 << " " << lyr << " stereo " << lyr % 2 << " plane1 " << p1 << " plane2 " << p2 << "\n";
238 }
239 for (int lyr = 0; lyr <= m_semax[1]; lyr++) { // Strip Negative Endcap
240 p1 = findPlane(m_planes, "se" + std::to_string(lyr) + "-");
241 p2 = findPlane(m_planes2, "se" + std::to_string(lyr) + "-");
242 m_pmap << "SCT 2 " << lyr/2 << " " << lyr << " stereo " << lyr % 2 << " plane1 " << p1 << " plane2 " << p2 << "\n";
243 }
244
245 // Region Map
246 std::string rmap_path = m_outFileName.value() + "region" + std::to_string(m_region) + ".rmap";
247 ATH_MSG_INFO("Creating rmap: " << rmap_path);
248 m_rmap.open(rmap_path, std::ofstream::out);
249 m_rmap << "towers 1 phi 16\n\n0\n";
250
254
258
259 m_pmap.close();
260 m_rmap.close();
261
262 // Step 1: Read the entire contents of the file
263 std::ifstream inputFile(pmap_path);
264 if (!inputFile) {
265 ATH_MSG_ERROR("Error: Unable to open file for reading: " << pmap_path);
266 return StatusCode::FAILURE;
267 }
268
269 std::ostringstream buffer;
270 buffer << inputFile.rdbuf(); // Reading the entire file into the stringstream
271 std::string fileContent = buffer.str();
272 inputFile.close();
273
274 // Step 2: Concatenate the content n times
275 std::string newContent;
276 for (int i = 0; i < m_nSlices; ++i) {
277 newContent += fileContent;
278 newContent += '\n';
279 }
280
281 // Step 3: Write the new content back to the same file
282 std::ofstream outputFile(pmap_path);
283 if (!outputFile) {
284 ATH_MSG_ERROR("Error: Unable to open file for writing: " << pmap_path);
285 return StatusCode::FAILURE;
286 }
287
288 outputFile << newContent;
289 outputFile.close();
290
291 return StatusCode::SUCCESS;
292}
293
294StatusCode FPGATrackSimMapMakerAlg::writeSubrmap(std::vector<FPGATrackSimHit> const & allHits)
295{
296 /* ---------- Create z-slices ---------- */
297
298 // BEFORE DOING ANYTHING ELSE, apply global trimming to m_key_etamods.
299 std::set<int> key_etamods;
300 // First, sum up the total number of hits in the key layer. There should be a one-liner for this...
301 float total_hits = 0;
302 for (auto const &etamod : m_key_etamods) {
303 total_hits += etamod.second;
304 }
305 ATH_MSG_INFO("Found " << total_hits << " hits in the key layer, applying global trim factor of " << m_globalTrim << "%");
306
307 // Then, do the trim.
308 for (auto const &etamod : m_key_etamods) {
309 if (m_globalTrim == 0 || ((etamod.second / total_hits) >= m_globalTrim*0.01)) {
310 key_etamods.insert(etamod.first);
311 } else {
312 ATH_MSG_INFO("Eta module " << etamod.first << " only contains " << etamod.second << " out of " << total_hits << " hits, excluding from slices.");
313 }
314 }
315
316 if (m_nSlices.value() == -1) m_nSlices.value() = key_etamods.size(); //default is full granularity slicing
317 float etasPerSlice = float(key_etamods.size())/m_nSlices.value();
318 std::vector<std::vector<int>> key_modules_for_slices; // indexed by slice, holds module eta values
319
320 // easier to use vector than set, convert m_key_etamods into key_etas
321 std::vector<int> key_etas;
322 std::vector<int> key_etas2; // used if 2D slicing
323 key_etas.insert(key_etas.end(), key_etamods.begin(), key_etamods.end());
324
325 for (unsigned i = 0; i < key_etas.size(); i++)
326 {
327 if (i >= (key_modules_for_slices.size() * etasPerSlice)) key_modules_for_slices.push_back(std::vector<int>());
328 key_modules_for_slices.back().push_back(key_etas[i]);
329 }
330
331 std::map<int, int> keymod2slice;
332 for (unsigned s = 0; s < key_modules_for_slices.size(); s++)
333 for (unsigned e = 0; e < key_modules_for_slices[s].size(); e++)
334 keymod2slice[key_modules_for_slices[s][e]] = s;
335
336 std::string key = m_keystring.value();
337 key.erase(std::remove(key.begin(), key.end(), ','), key.end());
338 key.erase(std::remove(key.begin(), key.end(), ' '), key.end());
339 std::string key2 = m_keystring2.value();
340 key2.erase(std::remove(key2.begin(), key2.end(), ','), key2.end());
341 key2.erase(std::remove(key2.begin(), key2.end(), ' '), key2.end());
342
343 ATH_MSG_INFO("Doing z-slicing");
344 if (key_etamods.size() == 0) ATH_MSG_ERROR("Found 0 slices using the keystring: '" << m_keystring << "'");
345 ATH_MSG_INFO("Nslices = " << std::to_string(m_nSlices.value()) << ":");
346
347 std::stringstream eta_slices;
348 for (unsigned s = 0; s < key_modules_for_slices.size(); s++){
349 for (unsigned e = 0; e < key_modules_for_slices[s].size(); e++){
350 eta_slices << key_modules_for_slices[s][e] << " ";
351 }
352 eta_slices << ", ";
353 }
354 ATH_MSG_INFO(eta_slices.str());
355
356 // 2D key layer slicing
357 if (m_key2)
358 { // make new slices from combinations of keylayer1 slices and keylayer2 slices
359
360 /*------------- setup keylayer2 keymodule to slice map -------------- */
361 std::vector<std::vector<int>> key_modules_for_slices2;
362 float etasPerSlice2 = float(m_key_etamods2.size())/m_nSlices.value();
363 key_etas2.insert(key_etas2.end(), m_key_etamods2.begin(), m_key_etamods2.end());
364 for (unsigned i = 0; i < key_etas2.size(); i++)
365 {
366 if (i >= (key_modules_for_slices2.size() * etasPerSlice2)) key_modules_for_slices2.push_back(std::vector<int>());
367 key_modules_for_slices2.back().push_back(key_etas2[i]);
368 }
369
370 std::map<int, int> keymod2slice2;
371 for (unsigned s = 0; s < key_modules_for_slices2.size(); s++)
372 for (unsigned e = 0; e < key_modules_for_slices2[s].size(); e++)
373 keymod2slice2[key_modules_for_slices2[s][e]] = s;
374 /*----------------------------------------------------------------- */
375
376 int new_nSlice = m_nSlices.value();
377 for (int s1 = 0; s1 < m_nSlices.value(); s1++){ // loop over keylayer1's slices
378 for (int s2 = 0; s2 < m_nSlices.value(); s2++){ // loop over keylayer2's slices
379
380 for (auto& pair: m_track2modules) // track loop
381 {
382 if (m_usedTracks.find(pair.first) != m_usedTracks.end()) continue; // skip if already assigned a slice to this track
383 bool key1 = false;
384 bool key2 = false;
385 for (auto& m: pair.second) // module loop
386 {
387 if (isOnKeyLayer(1,m->det,m->bec,m->lyr))
388 if (keymod2slice[m->eta] == s1) key1 = true;
389
390 if (isOnKeyLayer(2,m->det,m->bec,m->lyr))
391 if (keymod2slice2[m->eta] == s2) key2 = true;
392 }
393
394 if (key1 && key2)
395 {
396 int newSlice = m_nSlices.value()*s1 + s2;
397 m_track2slice[pair.first] = newSlice;
398 m_usedTracks.insert(pair.first);
399 if (newSlice + 1 > new_nSlice) new_nSlice = newSlice + 1; // find max slice, to set new total number of slices
400 }
401 }
402 }
403 }
404 m_nSlices.value() = new_nSlice;
405 ATH_MSG_INFO("These slices were further divided based on the key layer 2: '" << key2 << "'. Now nSlices = " << m_nSlices.value());
406 ATH_MSG_INFO("Using " << m_usedTracks.size() << " tracks out of " << m_maxEvents << ". The rest were missing a hit in one or both of the key layers");
407 }
408
409
410 // 1D key layer slicing
411 else
412 {
413 for (const FPGATrackSimHit& hit: allHits) // Fill the track to slice map
414 {
415 if (m_usedTracks.find(hit.getEventIndex()) != m_usedTracks.end()) continue; // skip if already done a hit from this track
416 if (isOnKeyLayer(1,hit.getDetType(),hit.getDetectorZone(), hit.getPhysLayer()))
417 { // if hit is in key layer, add it's barcode to the map
418 if (keymod2slice.count(hit.getEtaModule()) > 0) {
419 int s = keymod2slice[hit.getEtaModule()];
420 m_track2slice[hit.getEventIndex()] = s;
421 m_usedTracks.insert(hit.getEventIndex());
422 }
423 }
424 }
425 ATH_MSG_INFO("Using " << m_usedTracks.size() << " tracks out of " << m_maxEvents << ". The rest missed the key layer");
426 }
427
428 std::stringstream trim;
429 std::stringstream gtrim;
430 trim << std::fixed << std::setprecision(3) << m_trim;
431 gtrim << std::fixed << std::setprecision(3) << m_globalTrim;
432 std::string str_trim = trim.str();
433 std::string str_gtrim = gtrim.str();
434 int dot = str_trim.find_last_of(".");
435 str_trim.replace(dot,1,"p");
436 str_gtrim.replace(str_gtrim.find_last_of("."), 1, "p");
437
438 std::string subrmap_path = m_outFileName.value() + "region" + std::to_string(m_region) + ".subrmap";
439
440 ATH_MSG_INFO("Creating subrmap: " << subrmap_path);
441 m_subrmap.open(subrmap_path, std::ofstream::out);
442 m_subrmap << "towers " << m_nSlices.value() << " phi 16\n\n";
443
444 // Resize numTracks vector to be equal to the number of slices
445 // Now that this just stores module pointers we could loop over m_modules instead.
446 for (auto& pair: m_track2modules) {
447 for (Module* m: pair.second) {
448 if (m->numTracks.empty()) {
449 m->numTracks.resize(m_nSlices.value(), 0);
450 }
451 }
452 }
453
454
455 // Count tracks per slice
456 ATH_MSG_INFO("Counting number of tracks per slice.");
457 std::vector<std::vector<int>> slicedTracks (m_nSlices.value()); // vector of tracks, indexed by slice
458 for (auto trk: m_usedTracks) {
459 int s = m_track2slice[trk];
460 slicedTracks[s].push_back(trk);
461 std::vector<Module*> mods = m_track2modules[trk];
462 for (Module* mod: mods) {
463 mod->numTracks[s]++;
464 }
465 }
466
467 if (m_drawSlices)
468 drawSlices(allHits);
469
470 // Now do trimming and Fill slice2module map
471 int trimmed = 0;
472 for (int s = 0; s < m_nSlices.value(); s++) {
473 ATH_MSG_INFO("Applying local trimming in slice " << s);
474 for (auto trk : slicedTracks[s]) {
475 auto it = std::remove_if (m_track2modules[trk].begin(),
476 m_track2modules[trk].end(),
477 [&] (const Module* m) {
478 return 100 * ( float(m->numTracks[s]) / float(slicedTracks[s].size()) ) < m_trim;
479 });
480 trimmed += m_track2modules[trk].end() - it;
481 m_track2modules[trk].erase (it, m_track2modules[trk].end());
482
483 ATH_MSG_DEBUG("About to query trk2slice");
484 int s = m_track2slice[trk];
485 ATH_MSG_DEBUG("Queried trk2slice.");
486 // add all modules from track to slices
487 if (m_track2modules[trk].size() > 0) {
488 ATH_MSG_DEBUG("About to insert trk2modules");
489 m_slice2modules[s].insert( m_slice2modules[s].end(), m_track2modules[trk].begin(), m_track2modules[trk].end() );
490 ATH_MSG_DEBUG("Inserted trk2modules.");
491 }
492 }
493 }
494 ATH_MSG_INFO("Trimmed off " << trimmed << " modules that were hit by less than " << m_trim << "% of tracks");
495
496 /* ---------- Print z-slice map ---------- */
497
498 for (int s = 0; s < m_nSlices.value(); s++)
499 {
500 m_subrmap << s << "\n";
504
508 m_subrmap << "\n\n";
509 }
510
511 m_subrmap.close();
512 return StatusCode::SUCCESS;
513}
514
516{
517 std::stringstream trim;
518 std::stringstream gtrim;
519 trim << std::fixed << std::setprecision(3) << m_trim;
520 gtrim << std::fixed << std::setprecision(3) << m_globalTrim;
521 std::string str_trim = trim.str();
522 std::string str_gtrim = gtrim.str();
523 int dot = str_trim.find_last_of(".");
524 str_trim.replace(dot,1,"p");
525 str_gtrim.replace(str_gtrim.find_last_of("."), 1, "p");
526
527 std::string slicingType = "";
528 if (m_key2) slicingType = "2D";
529
530 std::string etapat_path = m_outFileName.value() + "region" + std::to_string(m_region) + ".patt";
531
532 ATH_MSG_INFO("Creating eta patterns file: " << etapat_path);
533 m_etapat.open(etapat_path, std::ofstream::out);
534
535 // assign logical layer to each module
536 for (auto& pair: m_track2modules) {
537 for (auto& m: pair.second)
538 {
539 if (m->det == SiliconTech::pixel && m->bec == DetectorZone::barrel) m->plane = findPlane(m_planes, "pb" + std::to_string(m->lyr));
540 if (m->det == SiliconTech::pixel && m->bec == DetectorZone::posEndcap) m->plane = findPlane(m_planes, "pe" + std::to_string(m->lyr) + "+");
541 if (m->det == SiliconTech::pixel && m->bec == DetectorZone::negEndcap) m->plane = findPlane(m_planes, "pe" + std::to_string(m->lyr) + "-");
542
543 if (m->det == SiliconTech::strip && m->bec == DetectorZone::barrel) m->plane = findPlane(m_planes, "sb" + std::to_string(m->lyr));
544 if (m->det == SiliconTech::strip && m->bec == DetectorZone::posEndcap) m->plane = findPlane(m_planes, "se" + std::to_string(m->lyr) + "+");
545 if (m->det == SiliconTech::strip && m->bec == DetectorZone::negEndcap) m->plane = findPlane(m_planes, "se" + std::to_string(m->lyr) + "-");
546 }
547 }
548
549 for(auto trk: m_usedTracks)
550 {
551 std::stringstream track_etapatts;
552 unsigned planesDone = 0;
553 for (unsigned p = 0; p < (m_planes)->size(); p++)
554 {
555 for (const Module* m : m_track2modules[trk]) {
556 if (m->plane == static_cast<int>(p))
557 {
558 track_etapatts << std::to_string(static_cast<int>(m->det)) << "\t" << std::to_string(static_cast<int>(m->bec)) << "\t" << std::to_string(m->eta) << "\t\t";
559 planesDone++;
560 break;
561 }
562 }
563 }
564 if (planesDone == (m_planes)->size())
565 m_etapat << track_etapatts.str() << "\n";
566
567 }
568 m_etapat.close();
569 return StatusCode::SUCCESS;
570}
571
572StatusCode FPGATrackSimMapMakerAlg::writeRadiiFile(std::vector<FPGATrackSimHit> const & allHits)
573{
574 // calculate mean radii.
575 m_radii.resize(m_nSlices.value(), std::vector<std::vector<float>>(m_planes2->size(),std::vector<float>(0)));
576 for (const auto& hit: allHits)
577 {
578 SiliconTech det = hit.getDetType();
579 DetectorZone bec = hit.getDetectorZone();
580 int lyr = hit.getPhysLayer();
581 int slice = m_track2slice[hit.getEventIndex()];
582 int plane = -1;
583 if (det == SiliconTech::pixel && bec == DetectorZone::barrel) plane = findPlane(m_planes2, "pb" + std::to_string(lyr));
584 if (det == SiliconTech::pixel && bec == DetectorZone::posEndcap) plane = findPlane(m_planes2, "pe" + std::to_string(lyr) + "+");
585 if (det == SiliconTech::pixel && bec == DetectorZone::negEndcap) plane = findPlane(m_planes2, "pe" + std::to_string(lyr) + "-");
586 if (det == SiliconTech::strip && bec == DetectorZone::barrel) plane = findPlane(m_planes2, "sb" + std::to_string(lyr));
587 if (det == SiliconTech::strip && bec == DetectorZone::posEndcap) plane = findPlane(m_planes2, "se" + std::to_string(lyr) + "+");
588 if (det == SiliconTech::strip && bec == DetectorZone::negEndcap) plane = findPlane(m_planes2, "se" + std::to_string(lyr) + "-");
589
590 if (plane != -1) {
591 m_radii[slice][plane].push_back(hit.getR());
592 }
593 }
594
595 // print file
596 std::string radii_path = m_outFileName.value() + "region" + std::to_string(m_region) + "_radii.txt";
597
598 ATH_MSG_INFO("Creating radii file: " << radii_path);
599 m_radfile.open(radii_path, std::ofstream::out);
600 for (int s = 0; s < m_nSlices.value(); s++){
601 m_radfile << std::to_string(s) << " ";
602 for (unsigned p = 0; p < (m_planes2)->size(); p++){
603 if (m_radii[s][p].size() != 0){
604 // "If left to type inference, op operates on values of the same type as
605 // init which can result in unwanted casting of the iterator elements."
606 // https://en.cppreference.com/w/cpp/algorithm/accumulate
607 float avg = std::accumulate(m_radii[s][p].begin(), m_radii[s][p].end(), 0.0f) / float(m_radii[s][p].size());
608 m_radfile << std::setprecision(3) << std::fixed << avg << " ";
609 } else {
610 int avg = -1;
611 m_radfile << avg << " ";
612 }
613 }
614 m_radfile << "\n";
615 }
616
617 // Calculate global mean radii by reversing the order of the above two loops.
618 m_radfile << -1 << " ";
619 for (unsigned p = 0; p < (m_planes2)->size(); p++) {
620 float avg = 0;
621 int count = 0;
622 for (int s = 0; s < m_nSlices.value(); s++) {
623 if (m_radii[s][p].size() != 0) {
624 avg = std::accumulate(m_radii[s][p].begin(), m_radii[s][p].end(), avg);
625 count += m_radii[s][p].size();
626 }
627 }
628 if (count > 0) {
629 avg /= float(count);
630 m_radfile << std::setprecision(3) << std::fixed << avg << " ";
631 } else {
632 m_radfile << -1 << " ";
633 }
634 }
635 m_radfile << std::endl;
636
637 m_radfile.close();
638
639 return StatusCode::SUCCESS;
640}
641
642StatusCode FPGATrackSimMapMakerAlg::writeMedianZFile(std::vector<FPGATrackSimHit> const & allHits)
643{
644 // calculate median z. We do this globally and slice-by-slice.
645 m_z.resize(m_nSlices.value(), std::vector<std::vector<float>>((m_planes2)->size(),std::vector<float>(0)));
646 for (const auto& hit: allHits)
647 {
648 SiliconTech det = hit.getDetType();
649 DetectorZone bec = hit.getDetectorZone();
650 int lyr = hit.getPhysLayer();
651 int slice = m_track2slice[hit.getEventIndex()];
652 int plane = -1;
653 if (det == SiliconTech::pixel && bec == DetectorZone::barrel) plane = findPlane(m_planes2, "pb" + std::to_string(lyr));
654 if (det == SiliconTech::pixel && bec == DetectorZone::posEndcap) plane = findPlane(m_planes2, "pe" + std::to_string(lyr) + "+");
655 if (det == SiliconTech::pixel && bec == DetectorZone::negEndcap) plane = findPlane(m_planes2, "pe" + std::to_string(lyr) + "-");
656 if (det == SiliconTech::strip && bec == DetectorZone::barrel) plane = findPlane(m_planes2, "sb" + std::to_string(lyr));
657 if (det == SiliconTech::strip && bec == DetectorZone::posEndcap) plane = findPlane(m_planes2, "se" + std::to_string(lyr) + "+");
658 if (det == SiliconTech::strip && bec == DetectorZone::negEndcap) plane = findPlane(m_planes2, "se" + std::to_string(lyr) + "-");
659
660 if (plane != -1) {
661 m_z[slice][plane].push_back(hit.getZ());
662 }
663 }
664
665 // print file
666 std::string zed_path = m_outFileName.value() + "region" + std::to_string(m_region) + "_z.txt";
667
668 ATH_MSG_INFO("Creating median z file: " << zed_path);
669 m_zedfile.open(zed_path, std::ofstream::out);
670 for (int s = 0; s < m_nSlices.value(); s++){
671 m_zedfile << std::to_string(s) << " ";
672 for (unsigned p = 0; p < (m_planes2)->size(); p++){
673 if (m_z[s][p].size() != 0){
674 float minZ = *std::min_element(m_z[s][p].begin(), m_z[s][p].end());
675 float maxZ = *std::max_element(m_z[s][p].begin(), m_z[s][p].end());
676 float median = (minZ + maxZ)/2;
677 m_zedfile << std::setprecision(3) << std::fixed << median << " ";
678 } else {
679 int median = -1;
680 m_zedfile << median << " ";
681 }
682 }
683 m_zedfile << std::endl;
684 }
685
686 // Now do this globally. Note: should this be meanZ instead of medianZ in the forward region?
687 m_zedfile << -1 << " ";
688 for (unsigned p = 0; p < (m_planes2)->size(); p++) {
689 float minZ = 0;
690 float maxZ = 0;
691 bool doneInitial = false;
692 for (int s = 0; s < m_nSlices.value(); s++) {
693 if (m_z[s][p].size() != 0) {
694 float newMinZ = *std::min_element(m_z[s][p].begin(), m_z[s][p].end());
695 float newMaxZ = *std::max_element(m_z[s][p].begin(), m_z[s][p].end());
696 // slightly clunky, but z can be positive and negative so there's not a sane initial/placeholder value.
697 if (doneInitial) {
698 minZ = std::min(minZ, newMinZ);
699 maxZ = std::max(maxZ, newMaxZ);
700 } else {
701 minZ = newMinZ;
702 maxZ = newMaxZ;
703 doneInitial = true;
704 }
705 }
706 }
707 if (doneInitial) {
708 float median = (minZ + maxZ)/2;
709 m_zedfile << std::setprecision(3) << std::fixed << median << " ";
710 } else {
711 int median = -1;
712 m_zedfile << median << " ";
713 }
714 }
715 m_zedfile << std::endl;
716
717 m_zedfile.close();
718
719 return StatusCode::SUCCESS;
720}
721
722// Finalize
724{
725 // Write the output
731 m_monitorFile.reset();
732 delete m_moduleRelabel;
733 m_moduleRelabel = nullptr;
734 return StatusCode::SUCCESS;
735}
736
737
739// Helpers
740
741void FPGATrackSimMapMakerAlg::drawSlices(std::vector<FPGATrackSimHit> const & allHits)
742{
743 m_monitorFile->cd();
744
745 std::vector<TH2F*> h_slicemap;
746 char *hname = new char[20];
747
748 for (unsigned i = 0; i < (unsigned)m_nSlices.value(); i++)
749 {
750 sprintf(hname,"rz_slice%u",i);
751 // This should just default to the entire range, I think.
752 // The user can reduce the binning.
753 TH2F *h = new TH2F(hname,hname,7000,-3500,3500,1200,0,1200);
754 h_slicemap.push_back(h);
755 }
756
757 for (const auto& hit: allHits)
758 {
759 if (m_usedTracks.find(hit.getEventIndex()) == m_usedTracks.end()) continue; // skip if we don't use this track
760 int s = m_track2slice[hit.getEventIndex()];
761 h_slicemap[s]->Fill(hit.getZ(),hit.getR());
762 }
763
764 for (int i = 0; i < m_nSlices.value(); i++)
765 h_slicemap[i]->Write();
766
767 delete [] hname;
768}
769
771{
772 int det = static_cast<int>(t_det);
773 int bec = static_cast<int>(t_bec);
774 if (keynum == 1)
775 if (m_keylayer["det"].find(det) != m_keylayer["det"].end() && m_keylayer["bec"].find(bec) != m_keylayer["bec"].end() && m_keylayer["lyr"].find(lyr) != m_keylayer["lyr"].end())
776 return true;
777
778 if (keynum == 2)
779 if (m_keylayer2["det"].find(det) != m_keylayer2["det"].end() && m_keylayer2["bec"].find(bec) != m_keylayer2["bec"].end() && m_keylayer2["lyr"].find(lyr) != m_keylayer2["lyr"].end())
780 return true;
781
782 return false;
783}
784
785int FPGATrackSimMapMakerAlg::findPlane(const std::vector<std::vector<std::string>>* planes, const std::string& test) // find what plane a layer is assigned to.
786{
787 int pcounter = 0;
788 for (auto& plane : *planes) {
789 for (auto& layer : plane) {
790 if (test == layer) return pcounter;
791 }
792 pcounter ++;
793 }
794 return -1;
795}
796
797std::string FPGATrackSimMapMakerAlg::makeRmapLines(std::vector<FPGATrackSimHit> const & hits, SiliconTech det, DetectorZone bec, int max)
798{
799 std::stringstream rmap_line;
800 std::set<int> etas, phis;
801
802 for(int lyr = 0; lyr <= max; lyr++)
803 {
804 etas.clear();
805 phis.clear();
806 for (const auto& hit: hits)
807 {
808 if(static_cast<int>(hit.getPhysLayer()) == lyr && hit.getDetectorZone() == bec) // cast from uint to int just to remove Wsign-compare warnings
809 {
810 etas.insert(hit.getEtaModule());
811 phis.insert(hit.getPhiModule());
812 }
813 }
814 if (etas.size() != 0) rmap_line << static_cast<int>(det) << " " << static_cast<int>(bec) << " " << lyr << " " << *phis.begin() << " " << *phis.rbegin() << " " << phis.size() << " " << *etas.begin() << " " << *etas.rbegin() << " " << etas.size() << "\n";
815 else rmap_line << static_cast<int>(det) << " " << static_cast<int>(bec) << " " << lyr << " 0 0 0 0 0 0\n";
816
817 }
818
819 return rmap_line.str();
820}
821
823{
824 // Parse keystring and define the Key Layer
825 std::string delimiter = ",";
826 std::string s = m_keystring.value();
827 std::string det, bec, lyr;
828 std::map <std::string, std::vector<std::string>> abrevs = { {"pb",{"pixel","barrel"}}, {"pe",{"pixel","endcap"}}, {"sb",{"strip","barrel"}}, {"se",{"strip","endcap"}} };
829 if( s.find(delimiter) != std::string::npos) // keylayer format is 'strip,barrel,4'
830 {
831 try {
832 std::string det = s.substr(0, s.find(delimiter));
833 s.erase(0, s.find(delimiter) + delimiter.length());
834 std::string bec = s.substr(0, s.find(delimiter));
835 s.erase(0, s.find(delimiter) + delimiter.length());
836 std::string lyr = s.substr(0, s.find(delimiter));
837 s.erase(0, s.find(delimiter) + delimiter.length());
838 m_keylayer["det"].insert(static_cast<int>(m_det2tech[det]));
839 m_keylayer["bec"].insert(static_cast<int>(m_bec2zone[bec]));
840 m_keylayer["lyr"].insert(std::stoi(lyr));
841 ATH_MSG_INFO("Using key layer: " << m_keystring.value());
842 } catch (...){
843 ATH_MSG_ERROR("Invalid KeyString: '" << m_keystring.value() << "'." << "Accepted formats are 'strip,posEndcap,2', 'pixel,barrel,3', or 'plane 0'");
844 }
845 }
846 else // keylayer format is 'plane 0'
847 {
848 std::string delimiter = " ";
849 try {
850 s.erase(0, s.find(delimiter) + delimiter.length());
851 std::string plane = s.substr(0, s.find(delimiter));
852 std::vector<std::string> s = (m_planes)->at(std::stoi(plane));
853 for (unsigned i = 0; i < s.size(); i++){
854 std::string reg = s[i].substr(0, 2);
855 std::vector<std::string> zone = abrevs[reg];
856 if (s[i].back() == '+') zone[1] = "posEndcap";
857 if (s[i].back() == '-') zone[1] = "negEndcap";
858 s[i].erase(std::remove(s[i].begin(), s[i].end(), '+'), s[i].end());
859 s[i].erase(std::remove(s[i].begin(), s[i].end(), '-'), s[i].end());
860 std::string lyr = s[i].substr(2);
861 m_keylayer["det"].insert(static_cast<int>(m_det2tech[zone[0]]));
862 m_keylayer["bec"].insert(static_cast<int>(m_bec2zone[zone[1]]));
863 m_keylayer["lyr"].insert(std::stoi(lyr));
864 }
865 ATH_MSG_INFO("Using key layer (in plane X format) " << m_keystring.value());
866 } catch (...){
867 ATH_MSG_ERROR("Invalid KeyString: '" << m_keystring.value() << "'." << "Accepted formats are 'strip,posEndcap,2', 'pixel,barrel,3', or 'plane 0'");
868 }
869 }
870
871 // if 2D slicing
872 if (m_keystring2.value() != "")
873 {
874 m_key2 = true;
875 std::string s = m_keystring2.value();
876 if( s.find(delimiter) != std::string::npos) // keylayer format is 'strip,barrel,8'
877 {
878 try {
879 std::string det = s.substr(0, s.find(delimiter));
880 s.erase(0, s.find(delimiter) + delimiter.length());
881 std::string bec = s.substr(0, s.find(delimiter));
882 s.erase(0, s.find(delimiter) + delimiter.length());
883 std::string lyr = s.substr(0, s.find(delimiter));
884 s.erase(0, s.find(delimiter) + delimiter.length());
885 m_keylayer2["det"].insert(static_cast<int>(m_det2tech[det]));
886 m_keylayer2["bec"].insert(static_cast<int>(m_bec2zone[bec]));
887 m_keylayer2["lyr"].insert(std::stoi(lyr));
888 ATH_MSG_INFO("Using key layer (2D): " << m_keystring2.value());
889 } catch (...){
890 ATH_MSG_ERROR("Invalid KeyString2: '" << m_keystring2.value() << "'." << "Accepted formats are 'strip,posEndcap,2', 'pixel,barrel,3', or 'plane 0'");
891 }
892 }
893 else // keylayer format is 'plane 0'
894 {
895 std::string delimiter = " ";
896 try {
897 s.erase(0, s.find(delimiter) + delimiter.length());
898 std::string plane = s.substr(0, s.find(delimiter));
899 std::vector<std::string> s = (m_planes)->at(std::stoi(plane));
900 for (unsigned i = 0; i < s.size(); i++){
901 std::string reg = s[i].substr(0, 2);
902 std::vector<std::string> zone = abrevs[reg];
903 if (s[i].back() == '+') zone[1] = "posEndcap";
904 if (s[i].back() == '-') zone[1] = "negEndcap";
905 s[i].erase(std::remove(s[i].begin(), s[i].end(), '+'), s[i].end());
906 s[i].erase(std::remove(s[i].begin(), s[i].end(), '-'), s[i].end());
907 std::string lyr = s[i].substr(2);
908 m_keylayer2["det"].insert(static_cast<int>(m_det2tech[zone[0]]));
909 m_keylayer2["bec"].insert(static_cast<int>(m_bec2zone[zone[1]]));
910 m_keylayer2["lyr"].insert(std::stoi(lyr));
911 }
912 ATH_MSG_INFO("Using key layer (2D) (in plane X format): " << m_keystring2.value());
913 } catch (...){
914 ATH_MSG_ERROR("Invalid KeyString2: '" << m_keystring2.value() << "'." << "Accepted formats are 'strip,posEndcap,2', 'pixel,barrel,3', or 'plane 0'");
915 }
916 }
917 }
918}
919
920std::string FPGATrackSimMapMakerAlg::makeSubrmapLines(std::vector<Module*> const & allmods, SiliconTech det, DetectorZone bec, int max)
921{
922 std::stringstream subrmap_line;
923 std::set<int> etas, phis;
924
925 std::vector<Module*> mods;
926 for (auto* mod: allmods) // just want modules in pb/pe/sb/se etc, not all at once
927 if (mod->det == det && mod->bec == bec) mods.push_back(mod);
928
929 for(int lyr = 0; lyr <= max; lyr++)
930 {
931 etas.clear();
932 phis.clear();
933 for (auto mod: mods)
934 {
935 if(mod->lyr == lyr)
936 {
937 etas.insert(mod->eta);
938 phis.insert(mod->phi);
939 }
940 }
941 if (etas.size() != 0) subrmap_line << static_cast<int>(det) << " " << static_cast<int>(bec) << " " << lyr << " " << *phis.begin() << " " << *phis.rbegin() << " " << phis.size() << " " << *etas.begin() << " " << *etas.rbegin() << " " << etas.size() << "\n";
942 else subrmap_line << static_cast<int>(det) << " " << static_cast<int>(bec) << " " << lyr << " 0 0 0 0 0 0\n";
943
944 }
945
946 return subrmap_line.str();
947}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
SiliconTech
DetectorZone
static Double_t ss
string trim(string s)
#define max(a, b)
Definition cfImp.cxx:41
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
Header file for AthHistogramAlgorithm.
std::unique_ptr< TFile > m_monitorFile
StatusCode writePmapAndRmap(std::vector< FPGATrackSimHit > const &pbHits, std::vector< FPGATrackSimHit > const &peHits, std::vector< FPGATrackSimHit > const &sbHits, std::vector< FPGATrackSimHit > const &seHits, int region)
Gaudi::Property< float > m_trim
std::string makeSubrmapLines(std::vector< Module * > const &allmods, SiliconTech det, DetectorZone bec, int max)
StatusCode writeRadiiFile(std::vector< FPGATrackSimHit > const &allHits)
ToolHandle< IFPGATrackSimInputTool > m_hitSGInputTool
std::vector< std::vector< std::vector< float > > > m_radii
FPGATrackSimMapMakerAlg(const std::string &name, ISvcLocator *pSvcLocator)
std::string makeRmapLines(std::vector< FPGATrackSimHit > const &hits, SiliconTech det, DetectorZone bec, int max)
bool isOnKeyLayer(int keynum, SiliconTech det, DetectorZone bec, int lyr)
std::map< std::string, std::set< int > > m_keylayer2
std::map< std::string, std::set< int > > m_keylayer
Gaudi::Property< std::string > m_description
std::vector< FPGATrackSimHit > m_sbHits
Gaudi::Property< float > m_globalTrim
std::map< FPGATrackSimModuleId, Module > m_modules
FPGATrackSimEventInputHeader m_eventHeader
std::map< std::string, SiliconTech > m_det2tech
std::vector< FPGATrackSimHit > m_seHits
std::vector< FPGATrackSimHit > m_pbHits
Gaudi::Property< std::string > m_keystring2
Gaudi::Property< std::string > m_outFileName
const std::vector< std::vector< std::string > > * m_planes
std::map< int, int > m_track2slice
ToolHandle< IFPGATrackSimEventInputHeaderTool > m_hitInputTool
int findPlane(const std::vector< std::vector< std::string > > *planes, const std::string &test)
Gaudi::Property< std::string > m_keystring
std::map< int, std::vector< Module * > > m_track2modules
Gaudi::Property< std::vector< std::vector< std::string > > > m_overridePlanes
StatusCode writeMedianZFile(std::vector< FPGATrackSimHit > const &allHits)
Gaudi::Property< bool > m_remapModules
const std::vector< std::vector< std::string > > * m_planes2
Gaudi::Property< bool > m_drawSlices
FPGATrackSimModuleRelabel * m_moduleRelabel
Gaudi::Property< int > m_region
std::vector< std::vector< std::vector< float > > > m_z
std::map< int, int > m_key_etamods
Gaudi::Property< std::vector< std::vector< std::string > > > m_overridePlanes2
std::vector< FPGATrackSimHit > m_peHits
Gaudi::Property< std::string > m_geoTag
std::map< int, std::vector< Module * > > m_slice2modules
std::map< std::string, DetectorZone > m_bec2zone
std::vector< FPGATrackSimHit > m_allHits
StatusCode readInputs(bool &done)
StatusCode writeSubrmap(std::vector< FPGATrackSimHit > const &allHits)
Gaudi::Property< int > m_maxEvents
ServiceHandle< IFPGATrackSimEventSelectionSvc > m_evtSel
void drawSlices(std::vector< FPGATrackSimHit > const &allHits)
Gaudi::Property< int > m_nSlices
STL class.
std::string find(const std::string &s)
return a remapped string
Definition hcg.cxx:138
int count(std::string s, const std::string &regx)
count how many occurances of a regx are in a string
Definition hcg.cxx:146
Definition dot.py:1
DataModel_detail::iterator< DVL > remove(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end, const T &value)
Specialization of remove for DataVector/List.
DataModel_detail::iterator< DVL > remove_if(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end, Predicate pred)
Specialization of remove_if for DataVector/List.