ATLAS Offline Software
Loading...
Searching...
No Matches
RpcDigitizationTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6//
7// RpcDigitizationTool
8// ------------
9// Authors:
10// Andrea Di Simone <Andrea.Di.Simone@cern.ch>
11// Gabriele Chiodini <gabriele.chiodini@le.infn.it>
12// Stefania Spagnolo <stefania.spagnolo@le.infn.it>
14
16
17// Inputs
18#include "GaudiKernel/SystemOfUnits.h"
19#include "GaudiKernel/PhysicalConstants.h"
21#include "GeoModelHelpers/TransformToStringConverter.h"
24
25// Geometry
30
31// run n. from geometry DB
36
37// Truth
41#include "GeoModelKernel/throwExcept.h"
42// Random Numbers
44#include "CLHEP/Random/RandExponential.h"
45#include "CLHEP/Random/RandFlat.h"
46#include "CLHEP/Random/RandGaussZiggurat.h"
47
48// Core includes
49#include <TString.h> // for Form
50
51#include <atomic>
52#include <fstream>
53#include <iostream>
54#include <sstream>
55#include <utility>
56
59
60// 12 charge points, 15 BetaGamma points, 180 efficiency points for fcp search
61namespace {
62 constexpr int N_Charge = 12;
63 constexpr int N_Velocity = 15;
64 constexpr std::array<double, N_Charge> Charge{0.1, 0.2, 0.3, 0.33, 0.4, 0.5, 0.6, 0.66, 0.7, 0.8, 0.9, 1.0};
65 constexpr std::array<double, N_Velocity> Velocity{0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 2.0, 3.0, 10.0, 100.0, 1000.0};
66 constexpr double Eff_garfield[N_Charge][N_Velocity] = {
67 {0.8648, 0.3476, 0.1407, 0.0618, 0.0368, 0.0234, 0.0150, 0.0120, 0.0096, 0.0079, 0.0038, 0.0041, 0.0035, 0.0049, 0.0054},
68 {0.9999, 0.9238, 0.6716, 0.4579, 0.3115, 0.2238, 0.1727, 0.1365, 0.1098, 0.0968, 0.0493, 0.0451, 0.0528, 0.0694, 0.0708},
69 {1.0000, 0.9978, 0.9517, 0.8226, 0.6750, 0.5611, 0.4674, 0.3913, 0.3458, 0.3086, 0.1818, 0.1677, 0.1805, 0.2307, 0.2421},
70 {1.0000, 0.9994, 0.9758, 0.8918, 0.7670, 0.6537, 0.5533, 0.4856, 0.4192, 0.3852, 0.2333, 0.2186, 0.2479, 0.2957, 0.2996},
71 {1.0000, 1.0000, 0.9972, 0.9699, 0.9022, 0.8200, 0.7417, 0.6660, 0.6094, 0.5622, 0.3846, 0.3617, 0.3847, 0.4578, 0.4583},
72 {1.0000, 1.0000, 0.9998, 0.9956, 0.9754, 0.9479, 0.9031, 0.8604, 0.8126, 0.7716, 0.5827, 0.5545, 0.5865, 0.6834, 0.6706},
73 {1.0000, 1.0000, 1.0000, 0.9997, 0.9968, 0.9876, 0.9689, 0.9464, 0.9221, 0.8967, 0.7634, 0.7385, 0.7615, 0.8250, 0.8309},
74 {1.0000, 1.0000, 1.0000, 1.0000, 0.9995, 0.9952, 0.9866, 0.9765, 0.9552, 0.9427, 0.8373, 0.8127, 0.8412, 0.8899, 0.8891},
75 {1.0000, 1.0000, 1.0000, 1.0000, 0.9995, 0.9981, 0.9918, 0.9803, 0.9754, 0.9602, 0.8730, 0.8564, 0.8746, 0.9178, 0.9261},
76 {1.0000, 1.0000, 1.0000, 1.0000, 1.0000, 0.9993, 0.9990, 0.9951, 0.9935, 0.9886, 0.9419, 0.9277, 0.9422, 0.9686, 0.9700},
77 {1.0000, 1.0000, 1.0000, 1.0000, 1.0000, 1.0000, 0.9998, 0.9996, 0.9980, 0.9966, 0.9786, 0.9718, 0.9748, 0.9875, 0.9882},
78 {1.0000, 1.0000, 1.0000, 1.0000, 1.0000, 1.0000, 0.9998, 1.0000, 0.9991, 0.9988, 0.9913, 0.9872, 0.9917, 0.9970, 0.9964}};
79 bool
80 validIndex(int idx, int arraySize){
81 return (idx>=0) and (idx<arraySize);
82 }
83} // namespace
84
85using namespace MuonGM;
86namespace {
87 constexpr double SIG_VEL = 4.8;
88}
89
90
91RpcDigitizationTool::RpcDigitizationTool(const std::string& type, const std::string& name, const IInterface* pIID) :
92 PileUpToolBase(type, name, pIID) {}
93
94// member function implementation
95//--------------------------------------------
97 ATH_MSG_DEBUG("RpcDigitizationTool:: in initialize()");
98 ATH_MSG_DEBUG("Configuration RpcDigitizationTool ");
99
100 ATH_MSG_DEBUG("InputObjectName " << m_inputHitCollectionName);
101 ATH_MSG_DEBUG("OutputObjectName " << m_outputDigitCollectionKey.key());
102 ATH_MSG_DEBUG("OutputSDOName " << m_outputSDO_CollectionKey.key());
103 ATH_MSG_DEBUG("WindowLowerOffset " << m_timeWindowLowerOffset);
104 ATH_MSG_DEBUG("WindowUpperOffset " << m_timeWindowUpperOffset);
105 ATH_MSG_DEBUG("DeadTime " << m_deadTime);
106 ATH_MSG_DEBUG("RndmSvc " << m_rndmSvc);
107 ATH_MSG_DEBUG("PatchForRpcTime " << m_patch_for_rpc_time);
108 ATH_MSG_DEBUG("RpcTimeShift " << m_rpc_time_shift);
109 ATH_MSG_DEBUG("RPC_TimeSchema " << m_RPC_TimeSchema);
110 ATH_MSG_DEBUG("RPCSDOareRPCDigits " << m_sdoAreOnlyDigits);
111
112 ATH_MSG_DEBUG("IgnoreRunDependentConfig " << m_ignoreRunDepConfig);
113 ATH_MSG_DEBUG("turnON_efficiency " << m_turnON_efficiency);
114 ATH_MSG_DEBUG("Efficiency_fromCOOL " << m_Efficiency_fromCOOL);
115 ATH_MSG_DEBUG("Efficiency_BIS78_fromCOOL" << m_Efficiency_BIS78_fromCOOL);
116 ATH_MSG_DEBUG("turnON_clustersize " << m_turnON_clustersize);
117 ATH_MSG_DEBUG("ClusterSize_fromCOOL " << m_ClusterSize_fromCOOL);
118 ATH_MSG_DEBUG("ClusterSize_BIS78_fromCOOL" << m_ClusterSize_BIS78_fromCOOL);
119 ATH_MSG_DEBUG("FirstClusterSizeInTail " << m_FirstClusterSizeInTail);
120 ATH_MSG_DEBUG("ClusterSize1_2uncorr " << m_ClusterSize1_2uncorr);
121 ATH_MSG_DEBUG("BOG_BOF_DoubletR2_OFF " << m_BOG_BOF_DoubletR2_OFF);
122 ATH_MSG_DEBUG("CutMaxClusterSize " << m_CutMaxClusterSize);
123 ATH_MSG_DEBUG("CutProjectedTracks " << m_CutProjectedTracks);
124 ATH_MSG_DEBUG("ValidationSetup " << m_validationSetup);
125 ATH_MSG_DEBUG("IncludePileUpTruth " << m_includePileUpTruth);
126 ATH_MSG_DEBUG("VetoPileUpTruthLinks " << m_vetoPileUpTruthLinks);
127
128 ATH_CHECK(m_detMgrKey.initialize());
129 if (m_onlyUseContainerName) { ATH_CHECK(m_mergeSvc.retrieve()); }
130 ATH_CHECK(detStore()->retrieve(m_idHelper));
131 // check the identifiers
132
133 ATH_MSG_INFO("Max Number of RPC Gas Gaps for these Identifiers = " << m_idHelper->gasGapMax());
134
135 // check the input object name
136 if (m_hitsContainerKey.key().empty()) {
137 ATH_MSG_FATAL("Property InputObjectName not set !");
138 return StatusCode::FAILURE;
139 }
141 ATH_MSG_DEBUG("Input objects in container : '" << m_inputHitCollectionName << "'");
142
143 // Initialize ReadHandleKey
144 ATH_CHECK(m_hitsContainerKey.initialize());
145
146 // initialize the output WriteHandleKeys
150 ATH_MSG_DEBUG("Output digits: '" << m_outputDigitCollectionKey.key() << "'");
151
152 // set the configuration based on run1/run2
154
155 ATH_MSG_DEBUG("Ready to read parameters for cluster simulation from file");
156
157 ATH_CHECK(m_rndmSvc.retrieve());
158
159 // fill the taginfo information
161
163
165 // m_turnON_clustersize=false;
166 m_BOF_id = m_idHelper->stationNameIndex("BOF");
167 m_BOG_id = m_idHelper->stationNameIndex("BOG");
168 m_BOS_id = m_idHelper->stationNameIndex("BOS");
169 m_BIL_id = m_idHelper->stationNameIndex("BIL");
170 m_BIS_id = m_idHelper->stationNameIndex("BIS");
172
173 return StatusCode::SUCCESS;
174}
175
177 // TODO This should all be in a conditions Alg
178 // Retrieve geometry config information from the database (RUN1, RUN2, etc...)
179 SmartIF<IGeoModelSvc> geoModel{Gaudi::svcLocator()->service("GeoModelSvc")};
180 if ( !geoModel ) {
181 ATH_MSG_ERROR("Could not locate GeoModelSvc");
182 return StatusCode::FAILURE;
183 }
184
185 // check the DetDescr version
186 std::string atlasVersion = geoModel->atlasVersion();
187
188 SmartIF<IRDBAccessSvc> rdbAccess{Gaudi::svcLocator()->service("RDBAccessSvc")};
189 if ( !rdbAccess ) {
190 ATH_MSG_ERROR("Could not locate RDBAccessSvc");
191 return StatusCode::FAILURE;
192 }
193
194 enum DataPeriod {Unknown, Run1, Run2, Run3, Run4 };
195 DataPeriod run = Unknown;
196
197 std::string configVal = "";
198
199 IRDBRecordset_ptr atlasCommonRec = rdbAccess->getRecordsetPtr("AtlasCommon", atlasVersion, "ATLAS");
200 if (atlasCommonRec->size() == 0) {
201 run = Run1;
202 } else {
203 configVal = (*atlasCommonRec)[0]->getString("CONFIG");
204 ATH_MSG_INFO("From DD Database, Configuration is " << configVal);
205 if (configVal == "RUN1") {
206 run = Run1;
207 } else if (configVal == "RUN2") {
208 run = Run2;
209 } else if (configVal == "RUN3") {
210 run = Run3;
211 } else if (configVal == "RUN4") {
212 run = Run4;
213 }
214 if (run == DataPeriod::Unknown) {
215 ATH_MSG_FATAL("Unexpected value for geometry config read from the database: " << configVal);
216 return StatusCode::FAILURE;
217 }
218 }
219 if (run == Run3 && m_idHelper->gasGapMax() < 3)
220 ATH_MSG_WARNING("Run3, configVal = " << configVal << " and GasGapMax =" << m_idHelper->gasGapMax());
221
222 if (run == Run1)
223 ATH_MSG_INFO("From Geometry DB: MuonSpectrometer configuration is: RUN1 or MuonGeometry = R.06");
224 else if (run == Run2)
225 ATH_MSG_INFO("From Geometry DB: MuonSpectrometer configuration is: RUN2 or MuonGeometry = R.07");
226 else if (run == Run3)
227 ATH_MSG_INFO("From Geometry DB: MuonSpectrometer configuration is: RUN3 or MuonGeometry = R.09");
228 else if (run == Run4)
229 ATH_MSG_INFO("From Geometry DB: MuonSpectrometer configuration is: RUN4 or MuonGeometry = R.10");
230
231 if (m_ignoreRunDepConfig == false) {
233 m_Efficiency_fromCOOL = false;
235 m_RPCInfoFromDb = false;
236 m_kill_deadstrips = false;
237 if (run == Run1) {
238 // m_BOG_BOF_DoubletR2_OFF = true
239 // m_Efficiency_fromCOOL = true
240 // m_ClusterSize_fromCOOL = true
242 if (configVal == "RUN1") { // MC12 setup
245 m_RPCInfoFromDb = true;
246 m_kill_deadstrips = true;
248 }
249 } else {
250 // m_BOG_BOF_DoubletR2_OFF = false # do not turn off at digitization the hits in the dbR=2 chambers in the feet
251 // m_Efficiency_fromCOOL = false # use common average values in python conf.
252 // m_ClusterSize_fromCOOL = false # use common average values in python conf.
254 if (run == Run2) { // MC15c setup
257 m_RPCInfoFromDb = true;
258 m_kill_deadstrips = false;
260 } else {
261 ATH_MSG_INFO("Run3/4: configuration parameter not from COOL");
262 m_Efficiency_fromCOOL = false;
264 m_RPCInfoFromDb = false;
265 m_kill_deadstrips = false;
266 }
267 }
268 ATH_MSG_INFO("RPC Run1/2/3-dependent configuration is enforced");
269 } else {
270 ATH_MSG_WARNING("Run1/2/3-dependent configuration is bypassed; be careful with option settings");
271 }
272
273 ATH_MSG_DEBUG("......RPC Efficiency_fromCOOL " << m_Efficiency_fromCOOL);
274 ATH_MSG_DEBUG("......RPC ClusterSize_fromCOOL " << m_ClusterSize_fromCOOL);
275 ATH_MSG_DEBUG("......RPC BOG_BOF_DoubletR2_OFF " << m_BOG_BOF_DoubletR2_OFF);
276 ATH_MSG_DEBUG("......RPC RPCInfoFromDb " << m_RPCInfoFromDb);
277 ATH_MSG_DEBUG("......RPC KillDeadStrips " << m_kill_deadstrips);
278 ATH_MSG_DEBUG("......RPC CutProjectedTracks " << m_CutProjectedTracks);
279
280
281 return StatusCode::SUCCESS;
282}
283
284template <class CondType>
285StatusCode RpcDigitizationTool::retrieveCondData(const EventContext& ctx,
287 const CondType* & condPtr) const {
288
289 if (key.empty()) {
290 ATH_MSG_DEBUG("No key has been configured for object "<<typeid(CondType).name()<<". Clear pointer");
291 condPtr = nullptr;
292 return StatusCode::SUCCESS;
293 }
294 SG::ReadCondHandle<CondType> readHandle{key, ctx};
295 if (!readHandle.isValid()){
296 ATH_MSG_FATAL("Failed to load conditions object "<<key.fullKey()<<".");
297 return StatusCode::FAILURE;
298 }
299 condPtr = readHandle.cptr();
300 return StatusCode::SUCCESS;
301
302}
303//--------------------------------------------
304StatusCode RpcDigitizationTool::prepareEvent(const EventContext& /*ctx*/, unsigned int) {
305 ATH_MSG_DEBUG("RpcDigitizationTool::in prepareEvent()");
306
307 // John's Hacks START
308 m_RPCHitCollList.clear();
309 m_thpcRPC = std::make_unique<TimedHitCollection<RPCSimHit>>();
310 // John's Hacks END
311
312 return StatusCode::SUCCESS;
313}
314
315//--------------------------------------------
316StatusCode RpcDigitizationTool::processBunchXing(int bunchXing, SubEventIterator bSubEvents, SubEventIterator eSubEvents) {
317 ATH_MSG_DEBUG("RpcDigitizationTool::in processBunchXing()");
318
320 TimedHitCollList hitCollList;
321
322 if (!(m_mergeSvc->retrieveSubSetEvtData(m_inputHitCollectionName, hitCollList, bunchXing, bSubEvents, eSubEvents).isSuccess()) &&
323 hitCollList.empty()) {
324 ATH_MSG_ERROR("Could not fill TimedHitCollList");
325 return StatusCode::FAILURE;
326 } else {
327 ATH_MSG_VERBOSE(hitCollList.size() << " RPCSimHitCollection with key " << m_inputHitCollectionName << " found");
328 }
329
330 TimedHitCollList::iterator iColl(hitCollList.begin());
331 TimedHitCollList::iterator endColl(hitCollList.end());
332
333 // Iterating over the list of collections
334 for (; iColl != endColl; ++iColl) {
335 RPCSimHitCollection* hitCollPtr = new RPCSimHitCollection(*iColl->second);
336 PileUpTimeEventIndex timeIndex(iColl->first);
337
338 ATH_MSG_DEBUG("RPCSimHitCollection found with " << hitCollPtr->size() << " hits");
339 ATH_MSG_VERBOSE("time index info. time: " << timeIndex.time() << " index: " << timeIndex.index() << " type: " << timeIndex.type());
340
341 m_thpcRPC->insert(timeIndex, hitCollPtr);
342 m_RPCHitCollList.emplace_back(hitCollPtr);
343 }
344
345 return StatusCode::SUCCESS;
346}
347
348//--------------------------------------------
349// Get next event and extract collection of hit collections:
350StatusCode RpcDigitizationTool::getNextEvent(const EventContext& ctx) {
351 ATH_MSG_DEBUG("RpcDigitizationTool::getNextEvent()");
352
353 // initialize pointer
354 m_thpcRPC.reset();
355
356 // get the container(s)
358
359 // In case of single hits container just load the collection using read handles
362 if (!hitCollection.isValid()) {
363 ATH_MSG_ERROR("Could not get RPCSimHitCollection container " << hitCollection.name() << " from store "
364 << hitCollection.store());
365 return StatusCode::FAILURE;
366 }
367
368 // create a new hits collection
369 m_thpcRPC = std::make_unique<TimedHitCollection<RPCSimHit>>(1);
370 m_thpcRPC->insert(0, hitCollection.cptr());
371 ATH_MSG_DEBUG("RPCSimHitCollection found with " << hitCollection->size() << " hits");
372
373 return StatusCode::SUCCESS;
374 }
375 // this is a list<pair<time_t, DataLink<RPCSimHitCollection> > >
376 TimedHitCollList hitCollList;
377
378 if (!(m_mergeSvc->retrieveSubEvtsData(m_inputHitCollectionName, hitCollList).isSuccess())) {
379 ATH_MSG_ERROR("Could not fill TimedHitCollList");
380 return StatusCode::FAILURE;
381 }
382 if (hitCollList.empty()) {
383 ATH_MSG_ERROR("TimedHitCollList has size 0");
384 return StatusCode::FAILURE;
385 } else {
386 ATH_MSG_DEBUG(hitCollList.size() << " RPCSimHitCollections with key " << m_inputHitCollectionName << " found");
387 }
388
389 // create a new hits collection
390 m_thpcRPC = std::make_unique<TimedHitCollection<RPCSimHit>>();
391 // now merge all collections into one
392 TimedHitCollList::iterator iColl(hitCollList.begin());
393 TimedHitCollList::iterator endColl(hitCollList.end());
394 while (iColl != endColl) {
395 const RPCSimHitCollection* p_collection(iColl->second);
396 m_thpcRPC->insert(iColl->first, p_collection);
397 // if ( m_debug ) ATH_MSG_DEBUG ( "RPCSimHitCollection found with "
398 // << p_collection->size() << " hits" ); // loop on the hit collections
399 ++iColl;
400 }
401 return StatusCode::SUCCESS;
402}
403
404//--------------------------------------------
405StatusCode RpcDigitizationTool::mergeEvent(const EventContext& ctx) {
406 StatusCode status = StatusCode::SUCCESS;
407
408 ATH_MSG_DEBUG("RpcDigitizationTool::in mergeEvent()");
409 // create and record the Digit container in StoreGate
411 ATH_CHECK(digitContainer.record(std::make_unique<RpcDigitContainer>(m_idHelper->module_hash_max())));
412 ATH_MSG_DEBUG("RpcDigitContainer recorded in StoreGate.");
413
414 // Create and record the SDO container in StoreGate
416 ATH_CHECK(sdoContainer.record(std::make_unique<MuonSimDataCollection>()));
417 ATH_MSG_DEBUG("RpcSDOCollection recorded in StoreGate.");
418
420 m_sdo_tmp_map.clear();
422
423 Collections_t collections;
424 status = doDigitization(ctx, collections, sdoContainer.ptr());
425 if (status.isFailure()) { ATH_MSG_ERROR("doDigitization Failed"); }
426 for (size_t coll_hash = 0; coll_hash < collections.size(); ++coll_hash) {
427 if (collections[coll_hash]) {
428 ATH_CHECK( digitContainer->addCollection (collections[coll_hash].release(), coll_hash) );
429 }
430 }
431
432 // Clean-up
433 m_RPCHitCollList.clear();
434
435 return status;
436}
437
438//--------------------------------------------
439StatusCode RpcDigitizationTool::processAllSubEvents(const EventContext& ctx) {
440 StatusCode status = StatusCode::SUCCESS;
441
442 // merging of the hit collection in getNextEvent method
443
444 ATH_MSG_DEBUG("RpcDigitizationTool::in digitize()");
445
446 // create and record the Digit container in StoreGate
448 ATH_CHECK(digitContainer.record(std::make_unique<RpcDigitContainer>(m_idHelper->module_hash_max())));
449 ATH_MSG_DEBUG("RpcDigitContainer recorded in StoreGate.");
450
451 // Create and record the SDO container in StoreGate
453 ATH_CHECK(sdoContainer.record(std::make_unique<MuonSimDataCollection>()));
454 ATH_MSG_DEBUG("RpcSDOCollection recorded in StoreGate.");
455
457 m_sdo_tmp_map.clear();
459
460 if (!m_thpcRPC) {
461 status = getNextEvent(ctx);
462 if (StatusCode::FAILURE == status) {
463 ATH_MSG_INFO("There are no RPC hits in this event");
464 return status; // there are no hits in this event
465 }
466 }
467
468 Collections_t collections;
469 ATH_CHECK(doDigitization(ctx, collections, sdoContainer.ptr()));
470 for (size_t coll_hash = 0; coll_hash < collections.size(); ++coll_hash) {
471 if (collections[coll_hash]) {
472 ATH_CHECK( digitContainer->addCollection (collections[coll_hash].release(), coll_hash) );
473 }
474 }
475
476 return status;
477}
478
479//--------------------------------------------
480StatusCode RpcDigitizationTool::doDigitization(const EventContext& ctx,
481 Collections_t& collections,
482 MuonSimDataCollection* sdoContainer) {
483 ATHRNG::RNGWrapper* rngWrapper = m_rndmSvc->getEngine(this);
484 rngWrapper->setSeed(name(), ctx);
485 CLHEP::HepRandomEngine* rndmEngine = rngWrapper->getEngine(ctx);
486
487 const MuonGM::MuonDetectorManager* detMgr{nullptr};
489
490
491 std::unique_ptr<RPCSimHitCollection> inputSimHitColl{std::make_unique<RPCSimHitCollection>("RPC_Hits")};
492
493
494 // get the iterator pairs for this DetEl
495 // iterate over hits
497
498 // Perform null check on m_thpcRPC
499 if (!m_thpcRPC) {
500 ATH_MSG_ERROR("m_thpcRPC is null");
501 return StatusCode::FAILURE;
502 }
503
504 struct SimDataContent {
505 Identifier channelId{};
506 std::vector<MuonSimData::Deposit> deposits;
507 Amg::Vector3D gpos{Amg::Vector3D::Zero()};
508 double simTime{0.};
509 };
510
511 while (m_thpcRPC->nextDetectorElement(i, e)) {
512 // to store the a single
513
514 std::map<Identifier, SimDataContent> channelSimDataMap;
515
516 // Loop over the hits:
517 while (i != e) {
518 ATH_MSG_DEBUG("RpcDigitizationTool::loop over the hits");
519
520 TimedHitPtr<RPCSimHit> phit(*i++);
521
522 // the hit
523 const RPCSimHit& hit(*phit);
524 // the hit id
525 const int idHit = hit.RPCid();
526 // the global time (G4 time + bunch time)
527 const double globalHitTime{hitTime(phit)};
528 // the G4 time or TOF from IP
529 const double G4Time{hit.globalTime()};
530 // the bunch time
531 const double bunchTime{globalHitTime - hit.globalTime()};
532
533 ATH_MSG_DEBUG("Global time " << globalHitTime << " G4 time " << G4Time << " Bunch time " << bunchTime);
534
535 if (!m_simHitValidKey.empty()) {
536 ATH_MSG_VERBOSE("Validation: globalHitTime, G4Time, BCtime = " << globalHitTime << " " << G4Time << " " << bunchTime);
537 inputSimHitColl->Emplace(idHit, globalHitTime, hit.localPosition(),
538 HepMcParticleLink::getRedirectedLink(phit->particleLink(), phit.eventId(), ctx), // This link should now correctly resolve to the TruthEvent McEventCollection in the main StoreGateSvc.
539 hit.postLocalPosition(),
540 hit.energyDeposit(), hit.stepLength(), hit.particleEncoding(), hit.kineticEnergy());
541 }
542
543 // convert sim id helper to offline id
544 const std::string stationName = m_muonHelper->GetStationName(idHit);
545 const int stationEta = m_muonHelper->GetZSector(idHit);
546 const int stationPhi = m_muonHelper->GetPhiSector(idHit);
547 const int doubletR = m_muonHelper->GetDoubletR(idHit);
548 const int doubletZ = m_muonHelper->GetDoubletZ(idHit);
549 const int doubletPhi = m_muonHelper->GetDoubletPhi(idHit);
550 int gasGap = m_muonHelper->GetGasGapLayer(idHit);
551
552 if (m_muonHelper->GetMeasuresPhi(idHit)) continue; // Skip phi strip . To be created after efficiency evaluation
553
554
555 bool isValid{false};
556 const Identifier elementID = m_idHelper->elementID(stationName,stationEta,stationPhi,doubletR, isValid);
557 if (!isValid) {
558 ATH_MSG_WARNING("Failed to construct the element ID from "<<stationName
559 <<", stationEta: "<<stationEta<<", stationPhi: "<<stationPhi<<", doubletR: "<<doubletR);
560 continue;
561 }
562 // construct Atlas identifier from components
563 ATH_MSG_DEBUG("creating id for hit in element:"
564 << " stationName " << stationName << " stationEta " << stationEta << " stationPhi " << stationPhi << " doubletR "
565 << doubletR << " doubletZ " << doubletZ << " doubletPhi " << doubletPhi << " gasGap " << gasGap);
566 const Identifier detElId{m_idHelper->channelID(elementID, doubletZ, doubletPhi, 1,0, 1, isValid)};
567 if (!isValid) {
568 continue;
569 }
570 const RpcReadoutElement* reEle = detMgr->getRpcReadoutElement(detElId);
572 if (false && reEle->rotatedRpcModule()) {
573 gasGap = gasGap == 1 ? 2 : 1;
574 }
575
576
577 bool isValidEta{false}, isValidPhi{false};
578 const Identifier idpaneleta = m_idHelper->channelID(elementID, doubletZ, doubletPhi, gasGap, 0, 1, isValidEta);
579 const Identifier idpanelphi = m_idHelper->channelID(elementID, doubletZ, doubletPhi, gasGap, 1, 1, isValidPhi);
580 if (!isValidEta || !isValidPhi) {
581 ATH_MSG_WARNING("Found an invalid identifier "
582 << " stationName " << stationName << " stationEta " << stationEta << " stationPhi " << stationPhi
583 << " doubletR " << doubletR << " doubletZ " << doubletZ << " doubletPhi " << doubletPhi << " gasGap "
584 << gasGap);
585 continue;
586 }
587 // loop on eta and phi to apply correlated efficiency between the two views
588
590 const double tmp_CorrJitter = m_idHelper->stationName(idpaneleta) < 2 ? m_CorrJitter_BIS78 : m_CorrJitter;
592 const double corrtimejitter = tmp_CorrJitter > 0.01 ?
593 CLHEP::RandGaussZiggurat::shoot(rndmEngine, 0., tmp_CorrJitter) : 0.; // correlated jitter
594 // handle here the special case where eta panel is dead => phi strip status (dead or eff.) cannot be resolved;
595 // measured panel eff. will be used in that case and no phi strip killing will happen
596
597
598 // Extrapolate the hit to the gas gap centre located at x=0
599 const Amg::Vector3D hitDir{(hit.postLocalPosition() - hit.localPosition()).unit()};
600 const Amg::Vector3D gapCentre = hit.localPosition() +
601 Amg::intersect<3>(hit.localPosition(), hitDir, Amg::Vector3D::UnitX(), 0).value_or(0) * hitDir;
602
603 std::array<int, 3> pcseta = physicalClusterSize(ctx, reEle, idpaneleta, gapCentre, rndmEngine); // set to one for new algorithms
604 ATH_MSG_VERBOSE("Simulated cluster on eta panel: size/first/last= " << pcseta[0] << "/" << pcseta[1] << "/" << pcseta[2]);
605 std::array<int, 3> pcsphi = physicalClusterSize(ctx, reEle, idpanelphi, gapCentre, rndmEngine); // set to one for new algorithms
606 ATH_MSG_VERBOSE("Simulated cluster on phi panel: size/first/last= " << pcsphi[0] << "/" << pcsphi[1] << "/" << pcsphi[2]);
607
608
609
610 // create Identifiers
611 const Identifier atlasRpcIdeta = m_idHelper->channelID(elementID, doubletZ, doubletPhi, gasGap, 0, pcseta[1], isValidEta);
612 const Identifier atlasRpcIdphi = m_idHelper->channelID(elementID, doubletZ, doubletPhi, gasGap, 1, pcsphi[1], isValidPhi);
613
614 const HepMcParticleLink particleLink = HepMcParticleLink::getRedirectedLink(phit->particleLink(), phit.eventId(), ctx); // This link should now correctly resolve to the TruthEvent McEventCollection in the main StoreGateSvc.
615 const auto [etaStripOn, phiStripOn] = detectionEfficiency(ctx, idpaneleta, idpanelphi, rndmEngine, particleLink);
616 ATH_MSG_DEBUG("SetPhiOn " << phiStripOn << " SetEtaOn " << etaStripOn);
617
618 for (bool imeasphi : {false, true}) {
619 if (!imeasphi && (!etaStripOn || !isValidEta)) continue;
620 if (imeasphi && (!phiStripOn || !isValidPhi)) continue;
621
622
623 // get Identifier and list of clusters for this projection
624 const Identifier& atlasId = !imeasphi ? atlasRpcIdeta : atlasRpcIdphi;
625 std::array<int, 3> pcs{!imeasphi ? pcseta : pcsphi};
626
627 ATH_MSG_DEBUG("SetOn: stationName " << stationName << " stationEta " << stationEta << " stationPhi " << stationPhi
628 << " doubletR " << doubletR << " doubletZ " << doubletZ << " doubletPhi " << doubletPhi
629 << " gasGap " << gasGap << " measphi " << imeasphi);
630
631 // pcs contains the cluster size, the first strip number and the last strip number of the cluster
632 pcs = TurnOnStrips(reEle, std::move(pcs), atlasId);
633 if (pcs[2] < 0){
634 continue;
635 }
636
637 ATH_MSG_DEBUG("Simulated cluster1: size/first/last= " << pcs[0] << "/" << pcs[1] << "/" << pcs[2]);
638
639
640 const Amg::Vector3D pos = fromSimHitToLayer(reEle, atlasId) * hit.localPosition();
641 const Amg::Vector3D gpos = reEle->transform(atlasId) * pos;
642
643 ATH_MSG_VERBOSE(" evt: "<<ctx.eventID().event_number()
644 <<" hit "<<m_idHelper->print_to_string(atlasId)
645 <<" local simHit "<<Amg::toString(hit.localPosition())
646 <<" corrected: "<<Amg::toString(pos)
647 <<" transform: "<<GeoTrf::toString(fromSimHitToLayer(reEle, atlasId))
648 <<" local strip: "<<Amg::toString(reEle->localToGlobalTransf(atlasId).inverse()*reEle->stripPos(atlasId))
649 <<" local strip (II): "<<Amg::toString(reEle->transform(atlasId).inverse()*reEle->stripPos(atlasId))
650 <<" global: "<<Amg::toString(gpos)
651 <<" strip Pos: "<<Amg::toString(reEle->stripPos(atlasId)));
652
653 // Calculate propagation time along readout strip in seconds
654 double proptime = PropagationTime(reEle, atlasId, gpos);
655
656 double tns = G4Time + proptime + corrtimejitter; // the time is in nanoseconds
657 ATH_MSG_VERBOSE("TOF+propagation time " << tns << " /s where proptime " << proptime << "/s");
658
659 double time = tns + bunchTime;
660 ATH_MSG_VERBOSE("final time in ns: BC+TOF+prop " << time << " /ns");
661
662 // pack propagation time along strip, bunch time and local hit position
663 long long int packedMCword = PackMCTruth(proptime, bunchTime, pos.y(), pos.z());
664 //cppcheck-suppress invalidPointerCast
665 double* b = reinterpret_cast<double*>(&packedMCword);
666
668 // create here deposit for MuonSimData
669 // MuonMCData first word is the packing of : proptime, bunchTime, posy, posz
670 // MuonMCData second word is the total hit time: bunchcTime+tof+proptime+correlatedJitter / ns
671 MuonSimData::Deposit deposit(particleLink, MuonMCData((*b), time)); // store tof+strip_propagation+corr.jitter
672 // MuonMCData((*b),G4Time+bunchTime+proptime )); // store tof+strip_propagation
673
674 // Do not store pile-up truth information
676 if (std::abs(hit.particleEncoding()) == 13 || hit.particleEncoding() == 0) {
677 if (channelSimDataMap.find(atlasId) == channelSimDataMap.end()) {
678 SimDataContent& content = channelSimDataMap[atlasId];
679 content.channelId = atlasId;
680 content.deposits.push_back(deposit);
681 content.gpos = reEle->transform(atlasId)*
682 fromSimHitToLayer(reEle,atlasId) * gapCentre;
683 content.simTime = hitTime(phit);
684 ATH_MSG_VERBOSE("adding SDO entry: r " << content.gpos.perp() << " z " << content.gpos.z());
685 }
686 }
687 }
688
689
690 //---------------------------------------------------------------------
691 // construct new digit and store it in the respective digit collection
692 // --------------------------------------------------------------------
693
694 // we create one digit-vector/deposit for each strip in the cluster
695 bool isValid{false};
696 for (int clus = pcs[1]; clus <= pcs[2]; ++clus) {
697 Identifier newId = m_idHelper->channelID(stationName, stationEta, stationPhi, doubletR, doubletZ,
698 doubletPhi, gasGap, imeasphi, clus, isValid);
699 if (!isValid) {
700 ATH_MSG_WARNING(__FILE__<<":"<<__LINE__<< "Channel "<< stationName<<" "<<stationEta<<" "<<stationPhi<<" "<< doubletR<<" "<<doubletZ
701 <<" "<< doubletPhi<<" "<< gasGap <<" "<< imeasphi<<" "<< clus<<" is invalid");
702 continue;
703 }
704
705 if (!m_idHelper->valid(newId)) {
706 if (stationName.find("BI") != std::string::npos) {
707 ATH_MSG_WARNING("Temporary skipping creation of RPC digit for stationName="
708 << stationName << ", eta=" << stationEta << ", phi=" << stationPhi << ", doubletR=" << doubletR
709 << ", doubletZ=" << doubletZ << ", doubletPhi=" << doubletPhi << ", gasGap=" << gasGap
710 << ", measuresPhi=" << imeasphi << ", strip=" << clus << ", cf. ATLASRECTS-6124");
711 return StatusCode::SUCCESS;
712 } else {
713 ATH_MSG_ERROR("Created an invalid id, aborting!");
714 m_idHelper->print(newId);
715 return StatusCode::FAILURE;
716 }
717 }
718
722 // One identifier but several deposits // name m_sdo_tmp_map is wrong call it m_sdo_map
723 if (m_sdo_tmp_map.find(newId) == m_sdo_tmp_map.end()) {
724 std::vector<MuonSimData::Deposit> newdeps;
725 newdeps.push_back(deposit);
726 m_sdo_tmp_map.insert(std::map<Identifier, std::vector<MuonSimData::Deposit>>::value_type(newId, newdeps));
727 } else {
728 m_sdo_tmp_map[newId].push_back(deposit);
729 }
730 } // end for cluster
731 } // loop on eta and phi
732 } // end loop hits
733
734 if (m_muonOnlySDOs) {
735 for (auto it = channelSimDataMap.begin(); it != channelSimDataMap.end(); ++it) {
736 MuonSimData simData(it->second.deposits, 0);
737 simData.setPosition(it->second.gpos);
738 simData.setTime(it->second.simTime);
739 auto insertResult = sdoContainer->insert(std::make_pair(it->first, simData));
740 if (!insertResult.second)
741 ATH_MSG_WARNING("Attention: this sdo is not recorded, since the identifier already exists in the sdoContainer map");
742 }
743 }
744
745 } // end loop detector elements
746
748
749 std::map<Identifier, std::vector<MuonSimData::Deposit>>::iterator map_iter = m_sdo_tmp_map.begin();
750 ATH_MSG_DEBUG("Start the digit map loop");
751
752 for (; map_iter != m_sdo_tmp_map.end(); ++map_iter) {
753 // Identifier
754 const Identifier theId = (*map_iter).first;
755 ATH_MSG_DEBUG("in the map loop: id " << m_idHelper->show_to_string(theId));
756 // Deposit
757 const std::vector<MuonSimData::Deposit> theDeps = (*map_iter).second;
758
759 // store the SDO from the muon
760 MuonSimData::Deposit theMuon; // useful beacuse it sorts the digits in ascending time.
761 std::multimap<double, MuonSimData::Deposit> times; // extract here time info from deposits.
762
763 // loop on the vector deposit
764 for (unsigned int k = 0; k < theDeps.size(); k++) {
765 double time = theDeps[k].second.secondEntry();
766 times.insert(std::multimap<double, MuonSimData::Deposit>::value_type(time, theDeps[k]));
767 }
768
769 // now iterate again over the multimap entries and store digits after dead time applied
770
771 IdContext rpcContext = m_idHelper->module_context(); // work on chamber context
772
773 std::multimap<double, MuonSimData::Deposit>::iterator map_dep_iter = times.begin();
774
775 // loop to suppress digits too close in time (emulate Front-End and CMA dead time)
776 double last_time = -10000; // init to high value
777 for (; map_dep_iter != times.end(); ++map_dep_iter) {
778 double currTime = (*map_dep_iter).first;
779 ATH_MSG_VERBOSE("deposit with time " << currTime);
780
782 // store (before any cut: all G4 hits) in the SDO container
783 // Identifier sdo and digit are the same
784 if (sdoContainer->find(theId) != sdoContainer->end()) // Identifier exist -> increase deposit
785 {
786 std::map<Identifier, MuonSimData>::const_iterator it = sdoContainer->find(theId);
787 std::vector<MuonSimData::Deposit> deps = ((*it).second).getdeposits();
788 deps.push_back((*map_dep_iter).second);
789 } else // Identifier does not exist -> create (Id,deposit)
790 {
791 std::vector<MuonSimData::Deposit> deposits;
792 deposits.push_back((*map_dep_iter).second);
793 std::pair<std::map<Identifier, MuonSimData>::iterator, bool> insertResult =
794 sdoContainer->insert(std::make_pair(theId, MuonSimData(deposits, 0)));
795 if (!insertResult.second)
797 "Attention TEMP: this sdo is not recorded, since the identifier already exists in the sdoContainer map");
798 }
799 }
800 // apply dead time
801 if (std::abs(currTime - last_time) > (m_deadTime)) {
802 ATH_MSG_DEBUG("deposit with time " << currTime << " is distant enough from previous (if any) hit on teh same strip");
803 last_time = (*map_dep_iter).first;
804
805 // first add time jitter to the time:
806 double uncorrjitter = 0;
807 double tmp_UncorrJitter = m_UncorrJitter;
808 if (m_idHelper->stationName(theId) < 2) tmp_UncorrJitter = m_UncorrJitter_BIS78;
809 if (tmp_UncorrJitter > 0.01) uncorrjitter = CLHEP::RandGaussZiggurat::shoot(rndmEngine, 0., tmp_UncorrJitter);
810 // Historically patch for the cavern background
811 // Now we subtract TOF from IP to assume full time calibrated detector (t=0 for particle from IP at light speed)
812 // We add a time shift to emulate FE global offset
813
814 const RpcReadoutElement* ele = detMgr->getRpcReadoutElement(theId);
815 Amg::Vector3D posi = ele->stripPos(theId);
816 double tp = m_patch_for_rpc_time ? posi.mag() / Gaudi::Units::c_light : 0.;
817 // Calculate propagation time for a hit at the center of the strip, to be subtructed as well as the nominal TOF
818 double propTimeFromStripCenter = PropagationTime(ele, theId, posi);
819 double newDigit_time = currTime + uncorrjitter + m_rpc_time_shift - tp - propTimeFromStripCenter;
820
821 double digi_ToT = -1.; // Time over threshold, for Narrow-gap RPCs only
822 if (m_idHelper->stationName(theId) < 2) digi_ToT = timeOverThreshold(rndmEngine); //mn
823
824 ATH_MSG_VERBOSE("last_time=currTime " << last_time << " jitter " << uncorrjitter << " TOFcorrection " << tp << " shift "
825 << m_rpc_time_shift << " newDigit_time " << newDigit_time);
826
827 // Apply readout window (sensitive detector time window)
828 bool outsideDigitizationWindow = outsideWindow(newDigit_time);
829 if (outsideDigitizationWindow) {
830 ATH_MSG_VERBOSE("hit outside digitization window - do not produce digits");
831 ATH_MSG_DEBUG("Hit outside time window!!"
832 << " hit time (ns) = " << newDigit_time << " timeWindow = " << m_timeWindowLowerOffset << " / "
834
835 continue;
836 }
837 // ok, let's store this digit
838 // this is an accepted hit to become digit
839 last_time = (*map_dep_iter).first;
840
841 std::unique_ptr<RpcDigit> newDigit = std::make_unique<RpcDigit>(theId, newDigit_time, digi_ToT, false);
842
843 Identifier elemId = m_idHelper->elementID(theId);
844 RpcDigitCollection* digitCollection = nullptr;
845
846 IdentifierHash coll_hash;
847 if (m_idHelper->get_hash(elemId, coll_hash, &rpcContext)) {
848 ATH_MSG_ERROR("Unable to get RPC hash id from RPC Digit collection "
849 << "context begin_index = " << rpcContext.begin_index()
850 << " context end_index = " << rpcContext.end_index() << " the identifier is \n"<<elemId);
851 }
852
853 // make new digit
854 ATH_MSG_DEBUG("Digit Id = " << m_idHelper->show_to_string(theId) << " digit time " << newDigit_time);
855
856 // remember new collection.
857 if (coll_hash >= collections.size()) {
858 collections.resize (coll_hash+1);
859 }
860 digitCollection = collections[coll_hash].get();
861 if (!digitCollection) {
862 collections[coll_hash] = std::make_unique<RpcDigitCollection>(elemId, coll_hash);
863 digitCollection = collections[coll_hash].get();
864 }
865 digitCollection->push_back(std::move(newDigit));
866
868 // put SDO collection in StoreGate
869 if (sdoContainer->find(theId) != sdoContainer->end()) {
870 std::map<Identifier, MuonSimData>::const_iterator it = sdoContainer->find(theId);
871 std::vector<MuonSimData::Deposit> deps = ((*it).second).getdeposits();
872 deps.push_back((*map_dep_iter).second);
873 } else {
874 std::vector<MuonSimData::Deposit> deposits;
875 deposits.push_back((*map_dep_iter).second);
876 std::pair<std::map<Identifier, MuonSimData>::iterator, bool> insertResult =
877 sdoContainer->insert(std::make_pair(theId, MuonSimData(deposits, 0)));
878 if (!insertResult.second)
880 "Attention: this sdo is not recorded, since teh identifier already exists in the sdoContainer map");
881 }
882 }
883
884 } else
885 ATH_MSG_DEBUG("discarding digit due to dead time: " << (*map_dep_iter).first << " " << last_time);
886 }
887
888 } // loop to suppress digits too close in time ended
889
890 // reset the pointer if it not null
891 m_thpcRPC.reset();
892 if (!m_simHitValidKey.empty()) {
894 ATH_CHECK(validHandle.record(std::move(inputSimHitColl)));
895 }
896
897 return StatusCode::SUCCESS;
898}
900 const Identifier& layerId) const {
901
902 Amg::Vector3D lGasGapPos = reEle->localGasGapPos(layerId);
903 if (reEle->NphiStripPanels() != reEle->nGasGapPerLay()) {
904 lGasGapPos.y() =0.;
905 }
906
910 const bool flip = reEle->numberOfLayers() == 2 &&
911 (m_idHelper->gasGap(layerId) == 2) != reEle->rotatedRpcModule();
912 const Amg::Transform3D fromHitToGap{reEle->transform(layerId).inverse() *
913 reEle->absTransform() * Amg::getTranslate3D(lGasGapPos) *
914 (flip ? Amg::getRotateY3D(180.*Gaudi::Units::deg) : Amg::Transform3D::Identity())};
915 ATH_MSG_VERBOSE("Transformation to go from hit to gap restframe "<<m_idHelper->print_to_string(layerId)
916 <<" "<<Amg::toString(fromHitToGap));
917 return fromHitToGap;
918}
919
920//--------------------------------------------
921std::array<int, 3> RpcDigitizationTool::physicalClusterSize(const EventContext& ctx,
922 const RpcReadoutElement* ele,
923 const Identifier& id,
924 const Amg::Vector3D& gapCentre,
925 CLHEP::HepRandomEngine* rndmEngine) const {
926
927
928 std::array<int, 3> result{};
929
930 const Amg::Vector3D position = fromSimHitToLayer(ele, id) * gapCentre;
931
932 const int doubletPhi = m_idHelper->doubletPhi(id);
933 const int gasGap = m_idHelper->gasGap(id);
934 const bool measuresPhi = m_idHelper->measuresPhi(id);
935 const double pitch= ele->StripPitch(measuresPhi);
936
937
938 const int nstrip = ele->stripNumber(position.block<2,1>(0,0), id);
939 const int numStrips = ele->Nstrips(measuresPhi);
940
941 result[1] = nstrip;
942 result[2] = nstrip;
943
944 if (nstrip < 1 || nstrip > numStrips) {
945 return make_array<int, 3>(-1);
946 }
947 const Amg::Vector3D locStripPos = ele->transform(id).inverse()*ele->stripPos(doubletPhi, gasGap, measuresPhi, nstrip);
948 float xstripnorm = (locStripPos -position).x() / pitch ;
949 result[0] = determineClusterSize(ctx, id, xstripnorm, rndmEngine);
950
951 //
952
953
954 if (m_turnON_clustersize == false) result[0] = 1;
955
956 return result;
957}
958
959//--------------------------------------------
961 std::array<int, 3>&& pcs,
962 const Identifier& id) const {
963
964
965 const int nstrips = ele->Nstrips(m_idHelper->measuresPhi(id));
966
967 if (pcs[0] == -2) {
968 pcs[1] = pcs[2] - 1;
969 } else if (pcs[0] == 2) {
970 pcs[2] = pcs[1] + 1;
971 } else if (pcs[0] > 2) {
972 pcs[1] = pcs[1] - pcs[0] / 2;
973 if (fmod(pcs[0], 2) == 0) pcs[1] = pcs[1] + 1;
974 pcs[2] = pcs[1] + pcs[0] - 1;
975 } else if (pcs[0] < -2) {
976 pcs[1] = pcs[1] + pcs[0] / 2;
977 pcs[2] = pcs[1] - pcs[0] - 1;
978 }
979
980 // cut the clusters at the beginning and at the end of the chamber
981
982 pcs[1] = std::clamp(pcs[1], 1, nstrips);
983 pcs[2] = std::clamp(pcs[2], 1, nstrips);
984
985 pcs[0] = pcs[2] - pcs[1] + 1;
986
987 return pcs;
988}
989
990//--------------------------------------------
992 const Identifier& id,
993 const Amg::Vector3D& globPos) const {
994
995 double distance{0.};
996 if (m_idHelper->measuresPhi(id)) {
997 distance = ele->distanceToPhiReadout(globPos);
998 } else {
999 distance = ele->distanceToEtaReadout(globPos);
1000 }
1001
1002 // distance in mm, SIG_VEL in ns/m
1003 return std::abs(distance * SIG_VEL * 1.e-3);
1004}
1005
1006//--------------------------------------------
1007long long int RpcDigitizationTool::PackMCTruth(float proptime, float bctime, float posy, float posz) const {
1008 // start with proptime: it is usually ~ns. It comes in ns. We express it in ns/10. use only 8 bits
1009 if (proptime < 0) {
1010 ATH_MSG_WARNING("A poblem: packing a propagation time <0 " << proptime << " redefine it as 0");
1011 proptime = 0.;
1012 }
1013 long long int new_proptime = int(proptime * 10) & 0xff;
1014
1015 // now tof. it is ~100ns. comes in ns. express it in ns/10. 16 bits needed (0-32768)
1016 // now BC time: it is ~100ns. comes in ns. express it in ns/10. 16 bits needed (0-32768)
1017 // can be negative (=> add 300 ns)
1018
1019 long long int new_bctime = int((bctime + 300.) * 10.) & 0xffff;
1020
1021 // posy: ~1000mm comes in mm, write it in mm*10. need 16 bits (0-32768)
1022 // can be negative (=>add 1500 mm)
1023
1024 long long int new_posy = int((posy + 1500.) * 10.) & 0xffff;
1025
1026 // posz: ~1000mm comes in mm, write it in mm*10. need 16 bits (0-32768)
1027 // can be negative (=>add 1500 mm)
1028
1029 long long int new_posz = int((posz + 1500.) * 10.) & 0xffff;
1030
1031 return (new_proptime + (new_bctime << 8) + (new_posy << 24) + (new_posz << 40));
1032}
1033
1034//--------------------------------------------
1035void RpcDigitizationTool::UnPackMCTruth(double theWord, float& proptime, float& bctime, float& posy, float& posz) {
1036 // int64_t is just a shorter way of writing long long int
1037 using Repacker = union
1038
1039 {
1040 double dWord;
1041
1042 int64_t iWord;
1043 };
1044 Repacker MCTruth;
1045 MCTruth.dWord = theWord;
1046 proptime = ((MCTruth.iWord) & 0x00000000000000ffLL) / 10.;
1047 bctime = (((MCTruth.iWord) & 0x0000000000ffff00LL) >> 8) / 10.;
1048 posy = (((MCTruth.iWord) & 0x000000ffff000000LL) >> 24) / 10.;
1049 posz = (((MCTruth.iWord) & 0x00ffff0000000000LL) >> 40) / 10.;
1050
1051 //
1052 bctime = bctime - 300.;
1053 posy = posy - 1500.;
1054 posz = posz - 1500.;
1055}
1056
1057//--------------------------------------------
1059 // get TagInfoMgr
1060 SmartIF<ITagInfoMgr> tagInfoMgr{Gaudi::svcLocator()->service("TagInfoMgr")}; // Tag Info Manager
1061 if (!tagInfoMgr) { return StatusCode::FAILURE; }
1062
1063 std::string RpctimeSchema = "";
1064 std::stringstream RpctimeShift;
1065 RpctimeShift << (int)m_rpc_time_shift;
1066
1068 RpctimeSchema = "Datalike_TOFoff_TimeShift" + RpctimeShift.str() + "nsec";
1069 } else {
1070 RpctimeSchema = "G4like_TOFon_TimeShift" + RpctimeShift.str() + "nsec";
1071 }
1072
1073 StatusCode sc = tagInfoMgr->addTag(m_RPC_TimeSchema, RpctimeSchema);
1074
1075 if (sc.isFailure()) {
1076 ATH_MSG_WARNING(m_RPC_TimeSchema << " " << RpctimeSchema << " not added to TagInfo ");
1077 return sc;
1078 }
1079
1080 ATH_MSG_DEBUG(m_RPC_TimeSchema << " " << RpctimeSchema << " added to TagInfo ");
1081
1082 return StatusCode::SUCCESS;
1083}
1084
1085
1086//--------------------------------------------
1087std::pair<bool,bool> RpcDigitizationTool::detectionEfficiency(const EventContext& ctx,
1088 const Identifier& IdEta,
1089 const Identifier& IdPhi,
1090 CLHEP::HepRandomEngine* rndmEngine,
1091 const HepMcParticleLink& trkParticle) const {
1092
1093
1094
1095 ATH_MSG_DEBUG("RpcDigitizationTool::in DetectionEfficiency");
1096
1097 ATH_MSG_DEBUG("EtaPanelId to look for Eff is " << m_idHelper->show_to_string(IdEta));
1098 ATH_MSG_DEBUG("PhiPanelId to look for Eff is " << m_idHelper->show_to_string(IdPhi));
1099
1100
1101 // dead spacers are not simulated in GEANT4 => their effect must be emulated in the digitizer as an effective max. efficiency = 99%
1102 // (spacers are 1x1cm^2 over a grid of 10x10cm^2 =? geometrical ineff. introduced is 1% for normal incidence)
1103 float maxGeomEff{0.99}, PhiAndEtaEff{0.99}, OnlyEtaEff{0.f}, OnlyPhiEff{0.f};
1104
1105 // 2=BML,3=BMS,4=BOL,5=BOS,8=BMF,9=BOF,10=BOG
1106 int stationName = m_idHelper->stationName(IdEta);
1107 int stationEta = m_idHelper->stationEta(IdEta);
1108 int doubletR = m_idHelper->doubletR(IdEta);
1109
1110 // remove feet extension. driven by joboption
1111 if (m_BOG_BOF_DoubletR2_OFF && (stationName == m_BOF_id || stationName == m_BOG_id) && doubletR == 2) {
1112 return std::make_pair(false, false);
1113 }
1114
1115
1116 if (!m_turnON_efficiency) {
1117 return std::make_pair(true, true);
1118 }
1119 bool etaStripOn{true}, phiStripOn{true};
1120
1121 // int stripetadead = 0 ; // not used
1122 // int stripphidead = 0 ; // not used
1123
1124 unsigned int index = stationName - 2;
1125 // BML and BMS, BOL and BOS come first (stationName= 2 and 3, 4 and 5 -> index 0-3)
1126 if (stationName > 5 && stationName < 50) index = index - 2;
1127 // BMF, BOF and BOG are 8,9,10 => must be 4,5 and 6
1128 else if (stationName > 50)
1129 index = index - 44;
1130 // BME and BOE 53 and 54 are at indices 7 and 8
1131
1132 if (!m_Efficiency_fromCOOL && stationName >= 2) {
1133 if (index > m_PhiAndEtaEff_A.size() || index > m_OnlyEtaEff_A.size() || index > m_OnlyPhiEff_A.size()) {
1134 THROW_EXCEPTION("Index out of array in Detection Efficiency SideA " << index << " stationName = " << stationName);
1135 }
1136
1137 PhiAndEtaEff = m_PhiAndEtaEff_A[index];
1138 OnlyEtaEff = m_OnlyEtaEff_A[index];
1139 OnlyPhiEff = m_OnlyPhiEff_A[index];
1140
1141 if (stationEta < 0) {
1142 if (index > m_PhiAndEtaEff_C.size() || index > m_OnlyEtaEff_C.size() || index > m_OnlyPhiEff_C.size()) {
1143 THROW_EXCEPTION("Index out of array in Detection Efficiency SideC " << index << " stationName = " << stationName);
1144 }
1145 PhiAndEtaEff = m_PhiAndEtaEff_C[index];
1146 OnlyEtaEff = m_OnlyEtaEff_C[index];
1147 OnlyPhiEff = m_OnlyPhiEff_C[index];
1148 }
1149 } else if (stationName < 2 && (!m_Efficiency_fromCOOL || !m_Efficiency_BIS78_fromCOOL)) { // BIS
1150 PhiAndEtaEff = m_PhiAndEtaEff_BIS78;
1151 OnlyEtaEff = m_OnlyEtaEff_BIS78;
1152 OnlyPhiEff = m_OnlyPhiEff_BIS78;
1153 } else { // Efficiency from Cool
1154
1155 const RpcCondDbData* readCdo{nullptr};
1156 if(!retrieveCondData(ctx, m_readKey, readCdo).isSuccess()){
1157 THROW_EXCEPTION("Failed to retrieve conditions object");
1158 }
1159
1160 ATH_MSG_DEBUG("Efficiencies and cluster size + dead strips will be extracted from COOL");
1161
1162 double FracDeadStripEta{0.}, FracDeadStripPhi{0.};
1163 double EtaPanelEfficiency{1.}, PhiPanelEfficiency{1.}, GapEfficiency{1.};
1164 int RPC_ProjectedTracksEta = 0;
1165
1166 std::optional<double> fracDeadStripEtaFromCOOL = readCdo->getFracDeadStrip(IdEta);
1167 std::optional<double> fracDeadStripPhiFromCOOL = readCdo->getFracDeadStrip(IdPhi);
1168
1169 bool noEntryInDb = !fracDeadStripEtaFromCOOL || !fracDeadStripPhiFromCOOL;
1170
1171 FracDeadStripEta = fracDeadStripEtaFromCOOL.value_or(0.);
1172 FracDeadStripPhi = fracDeadStripPhiFromCOOL.value_or(0.);
1173 RPC_ProjectedTracksEta = readCdo->getProjectedTrack(IdEta).value_or(0);
1174
1175 EtaPanelEfficiency = readCdo->getEfficiency(IdEta).value_or(1.);
1176 PhiPanelEfficiency = readCdo->getEfficiency(IdPhi).value_or(1.);
1177 GapEfficiency = readCdo->getGapEfficiency(IdEta).value_or(1.);
1178
1179 if (std::abs(FracDeadStripEta - 1.) < 0.001) {
1180 ATH_MSG_DEBUG("Watch out: SPECIAL CASE: Read from Cool: FracDeadStripEta/Phi "
1181 << FracDeadStripEta << "/" << FracDeadStripPhi << " RPC_ProjectedTracksEta " << RPC_ProjectedTracksEta
1182 << " Eta/PhiPanelEfficiency " << EtaPanelEfficiency << "/" << PhiPanelEfficiency << " gapEff " << GapEfficiency
1183 << " for gas gap " << m_idHelper->show_to_string(IdEta) << " id " << IdEta.get_identifier32().get_compact());
1184 // dead eta panel => cannot determine the strip status for phi strips
1185 // FracDeadStripPhi must be reset to 0. and undefinedPhiStripStatus = true
1186 FracDeadStripPhi = 0.;
1187 ATH_MSG_VERBOSE("Watch out: SPECIAL CASE: Resetting FracDeadStripPhi " << FracDeadStripPhi << " ignoring phi dead strips ");
1188 }
1189
1190 // special test
1191 // here redefining the efficiencies:
1192 // EtaPanelEfficiency = 0.92;
1193 // PhiPanelEfficiency = 0.85;
1194 // GapEfficiency = 0.97;
1195 bool changing = false;
1196 ATH_MSG_DEBUG("Read from Cool: FracDeadStripEta/Phi " << FracDeadStripEta << "/" << FracDeadStripPhi << " RPC_ProjectedTracksEta "
1197 << RPC_ProjectedTracksEta << " Eta/PhiPanelEfficiency " << EtaPanelEfficiency
1198 << "/" << PhiPanelEfficiency << " gapEff " << GapEfficiency);
1199 // if ((1.-FracDeadStripEta)<EtaPanelEfficiency)
1200 if ((maxGeomEff - FracDeadStripEta) - EtaPanelEfficiency < -0.011) {
1201 ATH_MSG_DEBUG("Ineff. from dead strips on Eta Panel larger that measured efficiency: deadFrac="
1202 << FracDeadStripEta << " Panel Eff=" << EtaPanelEfficiency << " for Panel " << m_idHelper->show_to_string(IdEta));
1203 ATH_MSG_DEBUG("... see the corresponding report from RpcDetectorStatusDbTool");
1204 // EtaPanelEfficiency = 1.-FracDeadStripEta;
1205 EtaPanelEfficiency = maxGeomEff - FracDeadStripEta;
1206 changing = true;
1207 }
1208 // if ((1.-FracDeadStripPhi)<PhiPanelEfficiency)
1209 if ((maxGeomEff - FracDeadStripPhi) - PhiPanelEfficiency < -0.011) {
1210 ATH_MSG_DEBUG("Ineff. from dead strips on Phi Panel larger that measured efficiency: deadFrac="
1211 << FracDeadStripPhi << " Panel Eff=" << PhiPanelEfficiency << " for Panel " << m_idHelper->show_to_string(IdPhi));
1212 ATH_MSG_DEBUG("... see the corresponding report among the warnings of RpcDetectorStatusDbTool");
1213 // PhiPanelEfficiency = 1.-FracDeadStripPhi;
1214 PhiPanelEfficiency = maxGeomEff - FracDeadStripPhi;
1215 changing = true;
1216 }
1217 // if ((1.-FracDeadStripEta*FracDeadStripPhi)<GapEfficiency)
1218 if ((maxGeomEff - FracDeadStripEta * FracDeadStripPhi) - GapEfficiency < -0.011) {
1219 ATH_MSG_DEBUG("Ineff. from dead strips on Eta/Phi Panels larger that measured EtaORPhi efficiency: deadFrac="
1220 << FracDeadStripEta * FracDeadStripPhi << " EtaORPhi Eff=" << GapEfficiency << " for GasGap "
1221 << m_idHelper->show_to_string(IdEta));
1222 ATH_MSG_DEBUG("... see the corresponding report among the warnings of RpcDetectorStatusDbTool");
1223 // GapEfficiency = 1.-FracDeadStripEta*FracDeadStripPhi;
1224 GapEfficiency = maxGeomEff - FracDeadStripEta * FracDeadStripPhi;
1225 changing = true;
1226 }
1227 if (changing)
1228 ATH_MSG_DEBUG("Rinormalized Values from Cool: FracDeadStripEta/Phi "
1229 << FracDeadStripEta << "/" << FracDeadStripPhi << " RPC_ProjectedTracksEta " << RPC_ProjectedTracksEta
1230 << " Eta/PhiPanelEfficiency " << EtaPanelEfficiency << "/" << PhiPanelEfficiency << " gapEff " << GapEfficiency);
1231
1232 // gabriele //..stefania - if there are dead strips renormalize the eff. to the active area
1233 if (m_kill_deadstrips) {
1234 if ((FracDeadStripEta > 0.0 && FracDeadStripEta < 1.0) || (FracDeadStripPhi > 0.0 && FracDeadStripPhi < 1.0) || (noEntryInDb)) {
1235 EtaPanelEfficiency = EtaPanelEfficiency / (maxGeomEff - FracDeadStripEta);
1236 PhiPanelEfficiency = PhiPanelEfficiency / (maxGeomEff - FracDeadStripPhi);
1237 GapEfficiency = GapEfficiency / (maxGeomEff - FracDeadStripEta * FracDeadStripPhi);
1238
1239 if (EtaPanelEfficiency > maxGeomEff) EtaPanelEfficiency = maxGeomEff;
1240 if (PhiPanelEfficiency > maxGeomEff) PhiPanelEfficiency = maxGeomEff;
1241 if (GapEfficiency > maxGeomEff) GapEfficiency = maxGeomEff;
1242
1243 if (EtaPanelEfficiency > GapEfficiency) GapEfficiency = EtaPanelEfficiency;
1244 if (PhiPanelEfficiency > GapEfficiency) GapEfficiency = PhiPanelEfficiency;
1245 ATH_MSG_DEBUG("Eff Redefined (to correct for deadfrac): FracDeadStripEta/Phi "
1246 << " Eta/PhiPanelEfficiency " << EtaPanelEfficiency << "/" << PhiPanelEfficiency << " gapEff "
1247 << GapEfficiency);
1248 }
1249 }
1250
1251 // values from COOLDB (eventually overwritten later)
1252 PhiAndEtaEff = float(EtaPanelEfficiency + PhiPanelEfficiency - GapEfficiency);
1253 if (PhiAndEtaEff < 0.) PhiAndEtaEff = 0.;
1254 OnlyEtaEff = float(EtaPanelEfficiency - PhiAndEtaEff);
1255 if (OnlyEtaEff < 0.) OnlyEtaEff = 0.;
1256 OnlyPhiEff = float(PhiPanelEfficiency - PhiAndEtaEff);
1257 if (OnlyPhiEff < 0.) OnlyPhiEff = 0.;
1258
1259 // special patch to be true only when m_Efficiency_fromCOOL=true and /RPC/DQMF/ELEMENT_STATUS tag is
1260 // RPCDQMFElementStatus_2012_Jaunuary_26
1261 bool applySpecialPatch = false;
1263 if (m_idHelper->stationName(IdEta) == 3)
1264 {
1265 if (std::abs(m_idHelper->stationEta(IdEta)) == 6 && m_idHelper->doubletR(IdEta) == 1 &&
1266 m_idHelper->doubletZ(IdEta) == 2 && m_idHelper->doubletPhi(IdEta) == 1) {
1267 applySpecialPatch = true;
1269 "Applying special patch for BMS at |eta|=6 lowPt plane -dbbZ=2 and dbPhi=1 ... will use default eff. for Id "
1270 << m_idHelper->show_to_string(IdEta));
1272 "Applying special patch: THIS HAS TO BE DONE IF /RPC/DQMF/ELEMENT_STATUS tag is "
1273 "RPCDQMFElementStatus_2012_Jaunuary_2");
1274 }
1275 }
1276 }
1277
1278 // if projected tracks number too low or inconsistent values get efficiencies from joboption and overwrite previous values
1279 if (applySpecialPatch || RPC_ProjectedTracksEta < m_CutProjectedTracks || RPC_ProjectedTracksEta > 10000000 ||
1280 EtaPanelEfficiency > 1 || EtaPanelEfficiency < 0 || PhiPanelEfficiency > 1 || PhiPanelEfficiency < 0 || GapEfficiency > 1 ||
1281 GapEfficiency < 0) {
1282 if (index > m_PhiAndEtaEff_A.size() || index > m_OnlyEtaEff_A.size() || index > m_OnlyPhiEff_A.size()) {
1283 THROW_EXCEPTION("Index out of array in Detection Efficiency SideA COOLDB" << index << " stationName = " << stationName);
1284 }
1285 if (RPC_ProjectedTracksEta < m_CutProjectedTracks)
1286 ATH_MSG_DEBUG("# of proj tracks = " << RPC_ProjectedTracksEta << " < cut = " << m_CutProjectedTracks
1287 << " resetting eff. from cool with default(python) values ");
1288
1289 PhiAndEtaEff = m_PhiAndEtaEff_A[index];
1290 OnlyEtaEff = m_OnlyEtaEff_A[index];
1291 OnlyPhiEff = m_OnlyPhiEff_A[index];
1292
1293 if (stationEta < 0) {
1294 if (index > m_PhiAndEtaEff_C.size() || index > m_OnlyEtaEff_C.size() || index > m_OnlyPhiEff_C.size()) {
1295 THROW_EXCEPTION("Index out of array in Detection Efficiency SideC COOLDB" << index << " stationName = " << stationName);
1296 }
1297 PhiAndEtaEff = m_PhiAndEtaEff_C[index];
1298 OnlyEtaEff = m_OnlyEtaEff_C[index];
1299 OnlyPhiEff = m_OnlyPhiEff_C[index];
1300 }
1301
1302 // if (m_applyEffThreshold) {
1303 // gabriele Set efficiency from dead strip fraction instead of nominal value
1304 float effgap = PhiAndEtaEff + OnlyEtaEff + OnlyPhiEff;
1305 float s_EtaPanelEfficiency = 1. - FracDeadStripEta;
1306 float s_PhiPanelEfficiency = 1. - FracDeadStripPhi;
1307 float s_PhiAndEtaEff = s_EtaPanelEfficiency * s_PhiPanelEfficiency / effgap;
1308 if (s_PhiAndEtaEff < PhiAndEtaEff) PhiAndEtaEff = s_PhiAndEtaEff;
1309 float s_OnlyEtaEff = s_EtaPanelEfficiency - PhiAndEtaEff;
1310 float s_OnlyPhiEff = s_PhiPanelEfficiency - PhiAndEtaEff;
1311
1312 if (s_OnlyEtaEff < OnlyEtaEff) OnlyEtaEff = s_OnlyEtaEff;
1313 if (s_OnlyPhiEff < OnlyPhiEff) OnlyPhiEff = s_OnlyPhiEff;
1314 // }
1315 }
1316
1317 float VolEff = PhiAndEtaEff + OnlyEtaEff + OnlyPhiEff;
1318 if (VolEff > maxGeomEff) {
1319 PhiAndEtaEff = (PhiAndEtaEff / VolEff) * maxGeomEff;
1320 OnlyEtaEff = (OnlyEtaEff / VolEff) * maxGeomEff;
1321 OnlyPhiEff = (OnlyPhiEff / VolEff) * maxGeomEff;
1322 }
1323
1324 } // End eff from COOL
1325
1326 // Efficiency correction factor for fractional-charged particles(added by Quanyin Li: quli@cern.ch)
1327 // link to truth particles and calculate the charge and betagamma
1328 HepMC::ConstGenParticlePtr genparticle = trkParticle.cptr();
1329 if (genparticle) {
1330 // only apply efficiency correction to fractional-charged particles based on pdgId betagamma
1331 if (MC::isGenericMultichargedParticle(genparticle)) {
1332 const double eff_sf = FCPEfficiency(genparticle);
1333 // Apply scale factor to the 3 Eff.
1334 PhiAndEtaEff = PhiAndEtaEff * eff_sf;
1335 OnlyEtaEff = OnlyEtaEff * eff_sf;
1336 OnlyPhiEff = OnlyPhiEff * eff_sf;
1337 }
1338 }
1339
1340 float I0 = PhiAndEtaEff;
1341 float I1 = PhiAndEtaEff + OnlyEtaEff;
1342 float ITot = PhiAndEtaEff + OnlyEtaEff + OnlyPhiEff;
1343
1344 float GapEff = ITot ;
1345 float PhiEff = PhiAndEtaEff + OnlyPhiEff;
1346 float EtaEff = PhiAndEtaEff + OnlyEtaEff;
1347
1348 ATH_MSG_DEBUG("DetectionEfficiency: Final Efficiency Values applied for "
1349 << m_idHelper->show_to_string(IdEta) << " are " << PhiAndEtaEff << "=PhiAndEtaEff " << OnlyEtaEff
1350 << "=OnlyEtaEff " << OnlyPhiEff << "=OnlyPhiEff " << GapEff << "=GapEff " << EtaEff << "=EtaEff " << PhiEff
1351 << "=PhiEff ");
1352
1353 float rndmEff = CLHEP::RandFlat::shoot(rndmEngine, 1);
1354
1355 if (rndmEff < I0) {
1356 phiStripOn = true;
1357 etaStripOn = true;
1358 } else if ((I0 <= rndmEff) && (rndmEff < I1)) {
1359 phiStripOn = false;
1360 etaStripOn = true;
1361 } else if ((I1 <= rndmEff) && (rndmEff <= ITot)) {
1362 phiStripOn = true;
1363 etaStripOn = false;
1364 } else {
1365 phiStripOn = false;
1366 etaStripOn = false;
1367 }
1368
1369 return std::make_pair(etaStripOn, phiStripOn);
1370}
1371
1372//--------------------------------------------
1374 const Identifier& idRpcStrip,
1375 double xstripnorm,
1376 CLHEP::HepRandomEngine* rndmEngine) const {
1377 ATH_MSG_DEBUG("RpcDigitizationTool::in determineClusterSize");
1378
1379 ATH_MSG_DEBUG("Digit Id = " << m_idHelper->show_to_string(idRpcStrip));
1380
1381 int ClusterSize = 1;
1382
1383 double FracClusterSize1{1.}, FracClusterSize2{0.}, MeanClusterSize{1.},
1384 FracClusterSizeTail{0.}, MeanClusterSizeTail{1.},
1385 FracClusterSize2norm{0.};
1386
1387 // 2=BML,3=BMS,4=BOL,5=BOS,8=BMF,9=BOF,10=BOG
1388 int stationName = m_idHelper->stationName(idRpcStrip);
1389 int stationEta = m_idHelper->stationEta(idRpcStrip);
1390 int measuresPhi = m_idHelper->measuresPhi(idRpcStrip);
1391
1392 unsigned int index = stationName - 2;
1393 // BML and BMS, BOL and BOS come first (stationName= 2 and 3, 4 and 5 -> index 0-3)
1394 if (stationName > 5 && stationName < 50) index = index - 2;
1395 // BMF, BOF and BOG are 8,9,10 => must be 4,5 and 6
1396 else if (stationName > 50)
1397 index = index - 44;
1398 // BME and BOE 53 and 54 are at indices 7 and 8
1399
1400 if (!m_ClusterSize_fromCOOL && stationName >= 2) {
1401 index += m_FracClusterSize1_A.size() / 2 * measuresPhi;
1402 if (index >= m_FracClusterSize1_A.size() ||
1403 index >= m_FracClusterSize2_A.size() ||
1404 index >= m_FracClusterSizeTail_A.size() ||
1405 index >= m_MeanClusterSizeTail_A.size()) {
1406 ATH_MSG_ERROR("Index out of array in determineClusterSize SideA " << index << " statName " << stationName);
1407 return 1;
1408 }
1409 FracClusterSize1 = m_FracClusterSize1_A[index];
1410 FracClusterSize2 = m_FracClusterSize2_A[index];
1411 FracClusterSizeTail = m_FracClusterSizeTail_A[index];
1412 MeanClusterSizeTail = m_MeanClusterSizeTail_A[index];
1413
1414 if (stationEta < 0) {
1415 index += m_FracClusterSize1_C.size() / 2 * measuresPhi - m_FracClusterSize1_A.size() / 2 * measuresPhi;
1416 if (index >= m_FracClusterSize1_C.size() ||
1417 index >= m_FracClusterSize2_C.size() ||
1418 index >= m_FracClusterSizeTail_C.size() ||
1419 index >= m_MeanClusterSizeTail_C.size()) {
1420 ATH_MSG_ERROR("Index out of array in determineClusterSize SideC " << index << " statName " << stationName);
1421 return 1;
1422 }
1423 FracClusterSize1 = m_FracClusterSize1_C[index];
1424 FracClusterSize2 = m_FracClusterSize2_C[index];
1425 FracClusterSizeTail = m_FracClusterSizeTail_C[index];
1426 MeanClusterSizeTail = m_MeanClusterSizeTail_C[index];
1427 }
1428 } else if (stationName < 2 && (!m_ClusterSize_fromCOOL || !m_ClusterSize_BIS78_fromCOOL)) { // BIS78
1429 FracClusterSize1 = m_FracClusterSize1_BIS78;
1430 FracClusterSize2 = m_FracClusterSize2_BIS78;
1431 FracClusterSizeTail = m_FracClusterSizeTail_BIS78;
1432 MeanClusterSizeTail = m_MeanClusterSizeTail_BIS78;
1433 } else { // Cluster size from COOL
1434 const RpcCondDbData* readCdo{nullptr};
1435 retrieveCondData(ctx, m_readKey, readCdo).ignore();
1436
1437 Identifier Id = m_idHelper->panelID(idRpcStrip);
1438
1439 int RPC_ProjectedTracks = readCdo->getProjectedTrack(Id).value_or(0);
1440 FracClusterSize1 = readCdo->getFracClusterSize1(Id).value_or(1.);
1441 FracClusterSize2 = readCdo->getFracClusterSize2(Id).value_or(0.);
1442 MeanClusterSize = readCdo->getMeanClusterSize(Id).value_or(1.);
1443
1444
1445 ATH_MSG_DEBUG("FracClusterSize1 and 2 " << FracClusterSize1 << " " << FracClusterSize2);
1446
1447 FracClusterSizeTail = 1. - FracClusterSize1 - FracClusterSize2;
1448
1449 MeanClusterSizeTail = MeanClusterSize - FracClusterSize1 - 2 * FracClusterSize2;
1450
1451 ATH_MSG_DEBUG("MeanClusterSizeTail and FracClusterSizeTail " << MeanClusterSizeTail << " " << FracClusterSizeTail);
1452
1453 // if clustersize have anomalous values set to the average cluster size from joboption
1454 if (RPC_ProjectedTracks < m_CutProjectedTracks || RPC_ProjectedTracks > 10000000 || MeanClusterSize > m_CutMaxClusterSize ||
1455 MeanClusterSize <= 1 || FracClusterSizeTail < 0 || FracClusterSize1 < 0 || FracClusterSize2 < 0 || FracClusterSizeTail > 1 ||
1456 FracClusterSize1 > 1 || FracClusterSize2 > 1) {
1457 if (stationName >= 2) {
1458 index += m_FracClusterSize1_A.size() / 2 * measuresPhi;
1459 if (index >= m_FracClusterSize1_A.size() ||
1460 index >= m_FracClusterSize2_A.size() ||
1461 index >= m_FracClusterSizeTail_A.size() ||
1462 index >= m_MeanClusterSizeTail_A.size()) {
1463 ATH_MSG_ERROR("Index out of array in determineClusterSize SideA " << index << " statName " << stationName);
1464 return 1;
1465 }
1466 FracClusterSize1 = m_FracClusterSize1_A[index];
1467 FracClusterSize2 = m_FracClusterSize2_A[index];
1468 FracClusterSizeTail = m_FracClusterSizeTail_A[index];
1469 MeanClusterSizeTail = m_MeanClusterSizeTail_A[index];
1470
1471 if (stationEta < 0) {
1472 index += m_FracClusterSize1_C.size() / 2 * measuresPhi - m_FracClusterSize1_A.size() / 2 * measuresPhi;
1473 if (index > m_FracClusterSize1_C.size() || index > m_FracClusterSize2_C.size() ||
1475 ATH_MSG_ERROR("Index out of array in determineClusterSize SideC " << index << " statName " << stationName);
1476 return 1;
1477 }
1478
1479 FracClusterSize1 = m_FracClusterSize1_C[index];
1480 FracClusterSize2 = m_FracClusterSize2_C[index];
1481 FracClusterSizeTail = m_FracClusterSizeTail_C[index];
1482 MeanClusterSizeTail = m_MeanClusterSizeTail_C[index];
1483 }
1484 } else {
1485 FracClusterSize1 = m_FracClusterSize1_BIS78;
1486 FracClusterSize2 = m_FracClusterSize2_BIS78;
1487 FracClusterSizeTail = m_FracClusterSizeTail_BIS78;
1488 MeanClusterSizeTail = m_MeanClusterSizeTail_BIS78;
1489 }
1490 }
1491 }
1492 FracClusterSize1 = std::min(FracClusterSize1, 1.);
1493 FracClusterSize2 = std::min(FracClusterSize2, 1.);
1494 FracClusterSizeTail = std::min(FracClusterSizeTail, 1.);
1495 float FracTot = FracClusterSize1 + FracClusterSize2 + FracClusterSizeTail;
1496 if (FracTot != 1. && FracTot > 0) {
1497 FracClusterSize1 = FracClusterSize1 / FracTot;
1498 FracClusterSize2 = FracClusterSize2 / FracTot;
1499 FracClusterSizeTail = FracClusterSizeTail / FracTot;
1500 }
1501 if (MeanClusterSizeTail < 0 || MeanClusterSizeTail > 10) MeanClusterSizeTail = 1;
1502
1503 ATH_MSG_VERBOSE("ClusterSize Final " << FracClusterSize1 << " FracClusterSize1 " << FracClusterSize2 << " FracClusterSize2 "
1504 << FracClusterSizeTail << " " << FracClusterSizeTail << " MeanClusterSizeTail "
1505 << MeanClusterSizeTail);
1506
1507 float FracClusterSize1plus2 = FracClusterSize1 + FracClusterSize2;
1508 float ITot = FracClusterSize1 + FracClusterSize2 + FracClusterSizeTail;
1509
1510 if (FracClusterSize1plus2 != 0) {
1511 // FracClusterSize1norm = FracClusterSize1 / FracClusterSize1plus2 ; // not used
1512 FracClusterSize2norm = FracClusterSize2 / FracClusterSize1plus2;
1513 }
1514
1515 float rndmCS = CLHEP::RandFlat::shoot(rndmEngine, ITot);
1516
1517 if (stationName >= 2) { // Legacy RPCs
1518 // Expanded CS2 of 1.3 to match average CS1 and CS2 (to be investigate)
1519 if (rndmCS < FracClusterSize1plus2) {
1520 // deterministic assignment of CS 1 or 2
1521 if (xstripnorm <= FracClusterSize2norm / 2. * 1.3) {
1522 ClusterSize = -2;
1523 } else if ((1.0 - FracClusterSize2norm / 2. * 1.3) <= xstripnorm) {
1524 ClusterSize = 2;
1525 } else {
1526 ClusterSize = 1;
1527 }
1529 float rndmCS1_2 = CLHEP::RandFlat::shoot(rndmEngine, 1);
1530 ClusterSize = 1 + (rndmCS1_2 < FracClusterSize2norm);
1531 }
1532
1533 } else if ((FracClusterSize1plus2 <= rndmCS) && (rndmCS <= ITot)) {
1534 ClusterSize = m_FirstClusterSizeInTail;
1535 ClusterSize += int(CLHEP::RandExponential::shoot(rndmEngine, MeanClusterSizeTail));
1536 float rndmLR = CLHEP::RandFlat::shoot(rndmEngine, 1.0);
1537 if (rndmLR > 0.5) ClusterSize = -ClusterSize;
1538 } else {
1539 ClusterSize = 1;
1540 }
1541
1542 } else { // NRPCs
1543 if (rndmCS < FracClusterSize1) {
1544 ClusterSize = 1;
1545 } else if (rndmCS < FracClusterSize1 + FracClusterSize2) {
1546 ClusterSize = 2;
1547 } else {
1548 ClusterSize = int(CLHEP::RandExponential::shoot(rndmEngine, MeanClusterSizeTail));
1549 }
1550 ClusterSize = std::max(ClusterSize, 1);
1551 if (ClusterSize > 1) {
1552 float rndmLR = CLHEP::RandFlat::shoot(rndmEngine, 1.0);
1553 if (rndmLR > 0.5) ClusterSize = -ClusterSize;
1554 }
1555 }
1556
1557 // negative CS correspond to left asymmetric cluster with respect to nstrip
1558 return ClusterSize;
1559}
1561 double qcharge = 1.;
1562 const int particlePdgId = genParticle->pdg_id();
1563 // charge calculation
1564 qcharge = (static_cast<double>((std::abs(particlePdgId) / 1000) % 100)) / (static_cast<double>((std::abs(particlePdgId) / 10) % 100));
1565 qcharge = ((static_cast<double>((static_cast<int>(qcharge * 100))))) / 100;
1566 if (particlePdgId < 0.0) qcharge = -qcharge;
1567 // BetaGamma calculation
1568 const double QPx = genParticle->momentum().px();
1569 const double QPy = genParticle->momentum().py();
1570 const double QPz = genParticle->momentum().pz();
1571 const double QE = genParticle->momentum().e();
1572 const double QM2 = std::pow(QE, 2) - std::pow(QPx, 2) - std::pow(QPy, 2) - std::pow(QPz, 2);
1573 const double QP = std::hypot(QPx, QPy, QPz);
1574 const double QM = QM2 >=0 ? std::sqrt(QM2) : -1.;
1575
1576 const double qbetagamma = QM > 0. ? QP / QM : -1.;
1577
1578 // find the i in the array
1579 int i_e = -1;
1580 for (int i = 0; i < 12; i++) {
1581 if (Charge[i] == std::abs(qcharge)) {
1582 i_e = i;
1583 break;
1584 }
1585 }
1586 int i_v = -99, j_v = 99;
1587 if (qbetagamma != -1) {
1588 for (int i = 0; i < 15; i++) {
1589 if (Velocity[i] <= qbetagamma) { i_v = i; }
1590 }
1591 for (int i = 14; i >= 0; i--) {
1592 if (Velocity[i] >= qbetagamma) { j_v = i; }
1593 }
1594 }
1595 // calculate the efficiency according to charge and velocity. Using linear function to calculate efficiency of a specific velocity
1596 // between velocity1 and velocity2
1597 double eff_fcp = 1.0, eff_muon = 1.0;
1598 if (i_e >= 0 && i_e <= 11) {
1599 if (validIndex(j_v, N_Velocity) && validIndex(i_v, N_Velocity) && (j_v - i_v) == 1) {
1600 const double delta_v = Velocity[i_v] - Velocity[j_v];
1601 eff_fcp = (Eff_garfield[i_e][i_v] - Eff_garfield[i_e][j_v]) / delta_v * qbetagamma +
1602 (Eff_garfield[i_e][j_v] * Velocity[i_v] - Eff_garfield[i_e][i_v] * Velocity[j_v]) / delta_v;
1603 eff_muon = (Eff_garfield[11][i_v] - Eff_garfield[11][j_v]) / delta_v * qbetagamma +
1604 (Eff_garfield[11][j_v] * Velocity[i_v] - Eff_garfield[11][i_v] * Velocity[j_v]) / delta_v;
1605 } else if (i_v == 14 && j_v == 99) {
1606 eff_fcp = Eff_garfield[i_e][14];
1607 eff_muon = Eff_garfield[11][14];
1608 } else if (i_v == -99 && j_v == 0) {
1609 eff_fcp = Eff_garfield[i_e][0];
1610 eff_muon = Eff_garfield[11][0];
1611 } else {
1612 ATH_MSG_WARNING("Wrong particle with unknown velocity! Scale factor is set to be 1.");
1613 }
1614 } else {
1615 ATH_MSG_WARNING("Wrong particle with unknown charge! Scale factor is set to be 1.");
1616 }
1617 // A scale factor is calculated by efficiency of fcp / efficiency of muon(charge==1.0
1618 const double eff_SF = eff_fcp / eff_muon;
1619 return eff_SF;
1620}
1621
1622double RpcDigitizationTool::timeOverThreshold(CLHEP::HepRandomEngine* rndmEngine) {
1623 //mn Time-over-threshold modeled as a narrow and a wide gaussian
1624 //mn based on the fit documented in https://its.cern.ch/jira/browse/ATLASRECTS-7820
1625 constexpr double tot_mean_narrow = 16.;
1626 constexpr double tot_sigma_narrow = 2.;
1627 constexpr double tot_mean_wide = 15.;
1628 constexpr double tot_sigma_wide = 4.5;
1629
1630 double thetot = 0.;
1631
1632 if (CLHEP::RandFlat::shoot(rndmEngine)<0.75) {
1633 thetot = CLHEP::RandGaussZiggurat::shoot(rndmEngine, tot_mean_narrow, tot_sigma_narrow);
1634 } else {
1635 thetot = CLHEP::RandGaussZiggurat::shoot(rndmEngine, tot_mean_wide, tot_sigma_wide);
1636 }
1637
1638 return (thetot > 0.) ? thetot : 0.;
1639}
float hitTime(const AFP_SIDSimHit &hit)
constexpr std::array< T, N > make_array(const T &def_val)
Helper function to initialize in-place arrays with non-zero values.
Definition ArrayHelper.h:10
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
bool isValid(const T &p)
Av: we implement here an ATLAS-sepcific convention: all particles which are 99xxxxx are fine.
Definition AtlasPID.h:878
ATLAS-specific HepMC functions.
std::vector< xAOD::EventInfo::SubEvent >::const_iterator SubEventIterator
Definition IPileUpTool.h:22
Definition of the abstract IRDBAccessSvc interface.
std::shared_ptr< IRDBRecordset > IRDBRecordset_ptr
Definition of the abstract IRDBRecord interface.
Definition of the abstract IRDBRecordset interface.
static Double_t sc
AtlasHitsVector< RPCSimHit > RPCSimHitCollection
@ Unknown
Definition TruthClasses.h:9
#define x
A wrapper class for event-slot-local random engines.
Definition RNGWrapper.h:56
void setSeed(const std::string &algName, const EventContext &ctx)
Set the random seed using a string (e.g.
Definition RNGWrapper.h:169
CLHEP::HepRandomEngine * getEngine(const EventContext &ctx) const
Retrieve the random engine corresponding to the provided EventContext.
Definition RNGWrapper.h:134
size_type size() const
const T * get(size_type n) const
Access an element, as an rvalue.
value_type push_back(value_type pElem)
Add an element to the end of the collection.
virtual unsigned int size() const =0
This class saves the "context" of an expanded identifier (ExpandedIdentifier) for compact or hash ver...
Definition IdContext.h:26
size_type begin_index() const
Definition IdContext.h:45
size_type end_index() const
Definition IdContext.h:46
value_type get_compact() const
Get the compact id.
This is a "hash" representation of an Identifier.
Identifier32 get_identifier32() const
Get the 32-bit version Identifier, will be invalid if >32 bits needed.
virtual const Amg::Transform3D & transform() const override
Return local to global transform.
The MuonDetectorManager stores the transient representation of the Muon Spectrometer geometry and pro...
const RpcReadoutElement * getRpcReadoutElement(const Identifier &id) const
access via extended identifier (requires unpacking)
An RpcReadoutElement corresponds to a single RPC module; therefore typicaly a barrel muon station con...
int nGasGapPerLay() const
returns the number of gasgaps
int NphiStripPanels() const
returns the number of phi strip panels (1 or 2)
virtual int stripNumber(const Amg::Vector2D &pos, const Identifier &id) const override final
strip number corresponding to local position.
Amg::Vector3D localGasGapPos(const Identifier &id) const
Returns the position of the gasGap w.r.t. rest frame of the chamber.
Amg::Transform3D localToGlobalTransf(const Identifier &id) const
int Nstrips(bool measphi) const
returns the number of strips for the phi or eta plane
virtual int numberOfLayers(bool measphi=true) const override final
number of layers in phi/eta projection, same for eta/phi planes
double StripPitch(bool measphi) const
returns the strip pitch for the phi or eta plane
std::pair< HepMcParticleLink, MuonMCData > Deposit
Definition MuonSimData.h:66
Gaudi::Property< int > m_vetoPileUpTruthLinks
PileUpToolBase(const std::string &type, const std::string &name, const IInterface *parent)
int particleEncoding() const
Definition RPCSimHit.h:61
double stepLength() const
Definition RPCSimHit.h:60
const Amg::Vector3D & postLocalPosition() const
Definition RPCSimHit.h:57
int RPCid() const
Definition RPCSimHit.h:63
const Amg::Vector3D & localPosition() const
Definition RPCSimHit.h:55
double kineticEnergy() const
Definition RPCSimHit.h:62
double globalTime() const
Definition RPCSimHit.h:54
double energyDeposit() const
Definition RPCSimHit.h:58
std::optional< double > getFracClusterSize2(const Identifier &) const
std::optional< int > getProjectedTrack(const Identifier &) const
std::optional< double > getFracDeadStrip(const Identifier &) const
std::optional< double > getFracClusterSize1(const Identifier &) const
std::optional< double > getGapEfficiency(const Identifier &) const
std::optional< double > getEfficiency(const Identifier &) const
std::optional< double > getMeanClusterSize(const Identifier &) const
std::pair< bool, bool > detectionEfficiency(const EventContext &ctx, const Identifier &ideta, const Identifier &idphi, CLHEP::HepRandomEngine *rndmEngine, const HepMcParticleLink &trkParticle) const
Evaluate detection efficiency.
Gaudi::Property< float > m_MeanClusterSizeTail_BIS78
Gaudi::Property< double > m_CorrJitter_BIS78
std::vector< std::unique_ptr< RpcDigitCollection > > Collections_t
bool outsideWindow(double time) const
const RpcHitIdHelper * m_muonHelper
int determineClusterSize(const EventContext &ctx, const Identifier &id, double xstripnorm, CLHEP::HepRandomEngine *rndmEngine) const
virtual StatusCode mergeEvent(const EventContext &ctx) override final
When being run from PileUpToolsAlgs, this method is called at the end of the subevts loop.
virtual StatusCode initialize() override final
Initialize.
long long int PackMCTruth(float proptime, float tof, float posx, float posz) const
Gaudi::Property< bool > m_turnON_efficiency
double FCPEfficiency(const HepMC::ConstGenParticlePtr &genParticle) const
Gaudi::Property< bool > m_EfficiencyPatchForBMShighEta
Gaudi::Property< bool > m_patch_for_rpc_time
std::array< int, 3 > physicalClusterSize(const EventContext &ctx, const MuonGM::RpcReadoutElement *reEle, const Identifier &id, const Amg::Vector3D &posAtCentre, CLHEP::HepRandomEngine *rndmEngine) const
Cluster simulation: first step.
RpcDigitizationTool(const std::string &type, const std::string &name, const IInterface *pIID)
Gaudi::Property< int > m_FirstClusterSizeInTail
Gaudi::Property< std::vector< double > > m_FracClusterSize2_A
ServiceHandle< IAthRNGSvc > m_rndmSvc
Gaudi::Property< double > m_timeWindowLowerOffset
Gaudi::Property< float > m_FracClusterSize2_BIS78
Amg::Transform3D fromSimHitToLayer(const MuonGM::RpcReadoutElement *readOutEle, const Identifier &layerId) const
Returns the position of the hit expressed in the gasGap coordinate system.
Gaudi::Property< int > m_CutProjectedTracks
Gaudi::Property< double > m_timeWindowUpperOffset
Gaudi::Property< std::vector< double > > m_FracClusterSize1_A
Gaudi::Property< int > m_deadTime
Gaudi::Property< bool > m_kill_deadstrips
Gaudi::Property< bool > m_ClusterSize_fromCOOL
std::array< int, 3 > TurnOnStrips(const MuonGM::RpcReadoutElement *reEle, std::array< int, 3 > &&pcs, const Identifier &id) const
Cluster simulation: second step.
Gaudi::Property< bool > m_includePileUpTruth
SG::WriteHandleKey< RPCSimHitCollection > m_simHitValidKey
SG::WriteHandleKey< RpcDigitContainer > m_outputDigitCollectionKey
StatusCode initializeRunDependentParameters()
Gaudi::Property< bool > m_ClusterSize_BIS78_fromCOOL
Gaudi::Property< bool > m_ClusterSize1_2uncorr
Gaudi::Property< bool > m_ignoreRunDepConfig
ServiceHandle< PileUpMergeSvc > m_mergeSvc
Gaudi::Property< std::vector< double > > m_MeanClusterSizeTail_A
std::vector< std::unique_ptr< RPCSimHitCollection > > m_RPCHitCollList
SG::WriteHandleKey< MuonSimDataCollection > m_outputSDO_CollectionKey
SG::ReadCondHandleKey< RpcCondDbData > m_readKey
Gaudi::Property< std::vector< double > > m_MeanClusterSizeTail_C
std::map< Identifier, std::vector< MuonSimData::Deposit > > m_sdo_tmp_map
Gaudi::Property< std::vector< float > > m_OnlyPhiEff_C
Gaudi::Property< bool > m_validationSetup
Gaudi::Property< std::vector< double > > m_FracClusterSize1_C
SG::ReadCondHandleKey< MuonGM::MuonDetectorManager > m_detMgrKey
virtual StatusCode prepareEvent(const EventContext &ctx, const unsigned int) override final
When being run from PileUpToolsAlgs, this method is called at the start of the subevts loop.
Gaudi::Property< double > m_UncorrJitter_BIS78
virtual StatusCode processBunchXing(int bunchXing, SubEventIterator bSubEvents, SubEventIterator eSubEvents) override final
When being run from PileUpToolsAlgs, this method is called for each active bunch-crossing to process ...
double PropagationTime(const MuonGM::RpcReadoutElement *reEle, const Identifier &id, const Amg::Vector3D &globPos) const
Calculates the propagation time along the strip.
Gaudi::Property< std::vector< float > > m_OnlyPhiEff_A
Gaudi::Property< float > m_FracClusterSizeTail_BIS78
Gaudi::Property< std::vector< double > > m_FracClusterSizeTail_A
Gaudi::Property< float > m_OnlyEtaEff_BIS78
Gaudi::Property< std::vector< double > > m_FracClusterSize2_C
std::unique_ptr< TimedHitCollection< RPCSimHit > > m_thpcRPC
StatusCode getNextEvent(const EventContext &ctx)
Get next event and extract collection of hit collections:
Gaudi::Property< std::vector< float > > m_OnlyEtaEff_C
Gaudi::Property< std::string > m_RPC_TimeSchema
const RpcIdHelper * m_idHelper
static void UnPackMCTruth(double theWord, float &proptime, float &tof, float &posy, float &posz)
static double timeOverThreshold(CLHEP::HepRandomEngine *rndmEngine)
Gaudi::Property< float > m_CutMaxClusterSize
SG::ReadHandleKey< RPCSimHitCollection > m_hitsContainerKey
Gaudi::Property< float > m_OnlyPhiEff_BIS78
Gaudi::Property< bool > m_turnON_clustersize
Gaudi::Property< std::vector< float > > m_PhiAndEtaEff_A
Gaudi::Property< bool > m_RPCInfoFromDb
Gaudi::Property< double > m_rpc_time_shift
Gaudi::Property< float > m_FracClusterSize1_BIS78
Gaudi::Property< double > m_CorrJitter
Gaudi::Property< std::vector< float > > m_OnlyEtaEff_A
Gaudi::Property< bool > m_Efficiency_fromCOOL
Gaudi::Property< double > m_UncorrJitter
Calculates the position of the hit wrt to the strip panel this transformation is needed since the imp...
std::string m_inputHitCollectionName
StatusCode doDigitization(const EventContext &ctx, Collections_t &collections, MuonSimDataCollection *sdoContainer)
Digitization functionality shared with RPC_PileUpTool.
virtual StatusCode processAllSubEvents(const EventContext &ctx) override final
alternative interface which uses the PileUpMergeSvc to obtain all the required SubEvents.
StatusCode retrieveCondData(const EventContext &ctx, const SG::ReadCondHandleKey< CondType > &key, const CondType *&condPtr) const
Gaudi::Property< bool > m_onlyUseContainerName
Gaudi::Property< std::vector< float > > m_PhiAndEtaEff_C
Gaudi::Property< bool > m_Efficiency_BIS78_fromCOOL
Gaudi::Property< std::vector< double > > m_FracClusterSizeTail_C
Gaudi::Property< bool > m_BOG_BOF_DoubletR2_OFF
Gaudi::Property< bool > m_muonOnlySDOs
Gaudi::Property< bool > m_sdoAreOnlyDigits
Gaudi::Property< float > m_PhiAndEtaEff_BIS78
static const RpcHitIdHelper * GetHelper(unsigned int nGasGaps=2)
const_pointer_type cptr()
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
std::string store() const
Return the name of the store holding the object we are proxying.
const std::string & name() const
Return the StoreGate ID for the referenced object.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
pointer_type ptr()
Dereference the pointer.
TimedVector::const_iterator const_iterator
a smart pointer to a hit that also provides access to the extended timing info of the host event.
Definition TimedHitPtr.h:18
unsigned short eventId() const
the index of the component event in PileUpEventInfo.
Definition TimedHitPtr.h:47
constexpr bool simData
Definition constants.h:36
Amg::Transform3D getTranslate3D(const double X, const double Y, const double Z)
: Returns a shift transformation along an arbitrary axis
std::optional< double > intersect(const AmgVector(N)&posA, const AmgVector(N)&dirA, const AmgVector(N)&posB, const AmgVector(N)&dirB)
Calculates the point B' along the line B that's closest to a second line A.
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
Amg::Transform3D getRotateY3D(double angle)
get a rotation transformation around Y-axis
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
bool ignoreTruthLink(const T &p, bool vetoPileUp)
Helper function for SDO creation in PileUpTools.
bool isGenericMultichargedParticle(const T &p)
In addition, there is a need to identify ”Q-ball” and similar very exotic (multi-charged) particles w...
Ensure that the Athena extensions are properly loaded.
Definition GeoMuonHits.h:27
USAGE: openCoraCool.exe "COOLONL_SCT/COMP200".
Definition index.py:1
std::list< value_t > type
type of the collection of timed data object
a struct encapsulating the identifier of a pile-up event
index_type index() const
the index of the component event in PileUpEventInfo
PileUpType type() const
the pileup type - minbias, cavern, beam halo, signal?
time_type time() const
bunch xing time in ns
#define THROW_EXCEPTION(MESSAGE)
Definition throwExcept.h:10
int run(int argc, char *argv[])