ATLAS Offline Software
Loading...
Searching...
No Matches
SiSmearedDigitizationTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6// SiSmearedDigitizationTool.cxx
7// Implementation file for class SiSmearedDigitizationTool
9
10// Pixel digitization includes
12
13// Det Descr
14#include "Identifier/Identifier.h"
17
23
24// Random numbers
26#include "CLHEP/Random/RandGaussZiggurat.h"
27#include "CLHEP/Random/RandFlat.h"
28#include "CLHEP/Random/RandGauss.h"
29#include "CLHEP/Random/RandLandau.h"
30#include "CLHEP/Vector/ThreeVector.h"
31
32// DataHandle
34
35// Pile-up
37
42
43// Fatras
54
56
57// Root
58#include "TTree.h"
59#include "TFile.h"
60
61#include <cmath>
62
63using namespace InDetDD;
64
65// Constructor with parameters:
66SiSmearedDigitizationTool::SiSmearedDigitizationTool(const std::string &type, const std::string &name,
67 const IInterface* parent):
68 PileUpToolBase(type, name, parent),
69 m_pixel_ID(nullptr),
70 m_sct_ID(nullptr),
71 m_randomEngineName("SiSmearedDigitization"),
72 m_pitch_X(0),
73 m_pitch_Y(0),
74 m_merge(false),
75 m_nSigma(0.),
76 m_useDiscSurface(false),
78 m_sctClusterContainer(nullptr),
79 m_mergeSvc("PileUpMergeSvc",name),
82 m_prdTruthNamePixel("PRD_MultiTruthPixel"),
83 m_prdTruthNameSCT("PRD_MultiTruthSCT"),
84 m_SmearPixel(true), //true: smear pixel --- false: smear SCT
85 m_emulateAtlas(true), // error rotation for endcap SCT
86 m_checkSmear(false),
87 m_outputFile(nullptr),
88 m_currentTree(nullptr),
89 m_x_pixel(0),
90 m_y_pixel(0),
100 m_x_SCT(0),
101 m_x_exit_SCT(0),
102 m_y_exit_SCT(0),
103 m_z_exit_SCT(0),
104 m_x_entry_SCT(0),
105 m_y_entry_SCT(0),
106 m_z_entry_SCT(0),
113 m_Err_x_pixel(0),
114 m_Err_y_pixel(0),
115 m_Err_x_SCT(0),
116 m_Err_y_SCT(0)
117{
118 declareProperty("RndmEngine", m_randomEngineName, "Random engine name");
119 declareProperty("InputObjectName", m_inputObjectName="PixelHits", "Input Object name" );
120 declareProperty("pitch_X", m_pitch_X);
121 declareProperty("pitch_Y", m_pitch_Y);
122 declareProperty("MergeClusters", m_merge);
123 declareProperty("Nsigma", m_nSigma);
124 declareProperty("SmearPixel", m_SmearPixel, "Enable Pixel or SCT Smearing");
125 declareProperty("PixelClusterContainerName", m_pixel_SiClustersName="PixelClusters");
126 declareProperty("SCT_ClusterContainerName", m_Sct_SiClustersName="SCT_Clusters");
127 declareProperty("CheckSmear", m_checkSmear);
128
129 declareProperty("HardScatterSplittingMode" , m_HardScatterSplittingMode, "Control pileup & signal splitting" );
130
131}
132
133
134// Initialize method:
136{
137
138 ATH_MSG_DEBUG ( "SiSmearedDigitizationTool::initialize()" );
139
140 //locate the AtRndmGenSvc and initialize our local ptr
141 if (!m_rndmSvc.retrieve().isSuccess())
142 {
143 ATH_MSG_ERROR ( "Could not find given RndmSvc" );
144 return StatusCode::FAILURE;
145 }
146
147 if (detStore()->retrieve(m_pixel_ID, "PixelID").isFailure()) {
148 ATH_MSG_ERROR ( "Could not get Pixel ID helper" );
149 return StatusCode::FAILURE;
150 }
151
152 if (not m_SmearPixel){ // Smear SCT
153 if (detStore()->retrieve(m_sct_ID, "SCT_ID").isFailure()) {
154 ATH_MSG_ERROR ( "Could not get SCT ID helper" );
155 return StatusCode::FAILURE;
156 }
157
158 m_inputObjectName="SCT_Hits"; // Set the input object name
159 }
160
161 // Initialize ReadCondHandleKeys
164
165 if (m_inputObjectName.empty())
166 {
167 ATH_MSG_FATAL ( "Property InputObjectName not set !" );
168 return StatusCode::FAILURE;
169 }
170 else
171 {
172 ATH_MSG_DEBUG ( "Input objects: '" << m_inputObjectName << "'" );
173 }
174
175 //locate the PileUpMergeSvc and initialize our local ptr
176 if (!m_mergeSvc.retrieve().isSuccess()) {
177 ATH_MSG_ERROR ( "Could not find PileUpMergeSvc" );
178 return StatusCode::FAILURE;
179 }
180
181 if (m_checkSmear){
182
183 // get THistSvc
184 ATH_CHECK(m_thistSvc.retrieve());
185
186 if (m_SmearPixel){
187 m_outputFile = new TFile("CheckSmearing_Pixel.root","RECREATE");
188 m_currentTree = new TTree("PixelTree","Check smearing Pixel");
189 m_currentTree->Branch("pixel_X" , &m_x_pixel , "m_x_pixel/D");
190 m_currentTree->Branch("pixel_Y" , &m_y_pixel , "m_y_pixel/D");
191 m_currentTree->Branch("pixel_X_exit" , &m_x_exit_pixel , "m_x_exit_pixel/D");
192 m_currentTree->Branch("pixel_Y_exit" , &m_y_exit_pixel , "m_y_exit_pixel/D");
193 m_currentTree->Branch("pixel_Z_exit" , &m_z_exit_pixel , "m_z_exit_pixel/D");
194 m_currentTree->Branch("pixel_X_entry" , &m_x_entry_pixel , "m_x_entry_pixel/D");
195 m_currentTree->Branch("pixel_Y_entry" , &m_y_entry_pixel , "m_y_entry_pixel/D");
196 m_currentTree->Branch("pixel_Z_entry" , &m_z_entry_pixel , "m_z_entry_pixel/D");
197 m_currentTree->Branch("pixel_X_global" , &m_x_pixel_global , "m_x_pixel_global/D");
198 m_currentTree->Branch("pixel_Y_global" , &m_y_pixel_global , "m_y_pixel_global/D");
199 m_currentTree->Branch("pixel_Z_global" , &m_z_pixel_global , "m_z_pixel_global/D");
200 m_currentTree->Branch("pixel_X_smear" , &m_x_pixel_smeared , "m_x_pixel_smeared/D");
201 m_currentTree->Branch("pixel_Y_smear" , &m_y_pixel_smeared , "m_y_pixel_smeared/D");
202 m_currentTree->Branch("pixel_Err_X" , &m_Err_x_pixel , "m_Err_x_pixel/D");
203 m_currentTree->Branch("pixel_Err_Y" , &m_Err_y_pixel , "m_Err_y_pixel/D");
204
205 if( m_thistSvc->regTree("PixelTree",m_currentTree).isFailure()) {
206 ATH_MSG_ERROR("Cannot register Ttree");
207 return StatusCode::FAILURE;
208 }else{
209 ATH_MSG_DEBUG ( "Ttree registered" );
210 }
211 }
212 else{
213 m_outputFile = new TFile("CheckSmearing_SCT.root","RECREATE");
214 m_currentTree = new TTree("SCT_Tree","Check smearing SCT");
215 m_currentTree->Branch("SCT_X" , &m_x_SCT , "m_x_SCT/D");
216 m_currentTree->Branch("SCT_exit_X" , &m_x_exit_SCT , "m_x_exit_SCT/D");
217 m_currentTree->Branch("SCT_exit_Y" , &m_y_exit_SCT , "m_y_exit_SCT/D");
218 m_currentTree->Branch("SCT_exit_Z" , &m_z_exit_SCT , "m_z_exit_SCT/D");
219 m_currentTree->Branch("SCT_entry_X" , &m_x_entry_SCT , "m_x_entry_SCT/D");
220 m_currentTree->Branch("SCT_entry_Y" , &m_y_entry_SCT , "m_y_entry_SCT/D");
221 m_currentTree->Branch("SCT_entry_Z" , &m_z_entry_SCT , "m_z_entry_SCT/D");
222 m_currentTree->Branch("SCT_X_global" , &m_x_SCT_global , "m_x_SCT_global/D");
223 m_currentTree->Branch("SCT_Y_global" , &m_y_SCT_global , "m_y_SCT_global/D");
224 m_currentTree->Branch("SCT_Z_global" , &m_z_SCT_global , "m_z_SCT_global/D");
225 m_currentTree->Branch("SCT_X_smear" , &m_x_SCT_smeared , "m_x_SCT_smeared/D");
226 m_currentTree->Branch("SCT_Err_X" , &m_Err_x_SCT , "m_Err_x_SCT/D");
227
228 if( m_thistSvc->regTree("SCT_Tree",m_currentTree).isFailure()) {
229 ATH_MSG_ERROR("Cannot register Ttree");
230 return StatusCode::FAILURE;
231 }else{
232 ATH_MSG_DEBUG ( "Ttree registered" );
233 }
234 }
235
236
237 }
238
239 return StatusCode::SUCCESS;
240}
241
242// Finalize method:
244{
245
246 if (m_checkSmear){
247 m_outputFile->cd();
248 m_currentTree->Write();
249 m_outputFile->Close();
250 ATH_MSG_DEBUG ( "SiSmearedDigitizationTool : Writing Tree" );
251
252 }
253
254 ATH_MSG_DEBUG ( "SiSmearedDigitizationTool : finalize()" );
255
256
257 return StatusCode::SUCCESS;
258
259}
260
261StatusCode SiSmearedDigitizationTool::prepareEvent(const EventContext& /*ctx*/, unsigned int)
262{
263
264 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: in pixel prepareEvent() ---" );
265
266 m_siHitCollList.clear();
268
269 return StatusCode::SUCCESS;
270}
271
272
274 SubEventIterator bSubEvents,
275 SubEventIterator eSubEvents)
276{
277 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: in pixel processBunchXing() ---" );
278 //decide if this event will be processed depending on HardScatterSplittingMode & bunchXing
279 if (m_HardScatterSplittingMode == 2 && !m_HardScatterSplittingSkipper ) { m_HardScatterSplittingSkipper = true; return StatusCode::SUCCESS; }
280 if (m_HardScatterSplittingMode == 1 && m_HardScatterSplittingSkipper ) { return StatusCode::SUCCESS; }
282
284 TimedHitCollList hitCollList;
285
286 if (!(m_mergeSvc->retrieveSubSetEvtData(m_inputObjectName, hitCollList, bunchXing,
287 bSubEvents, eSubEvents).isSuccess()) &&
288 hitCollList.empty()) {
289 ATH_MSG_ERROR("Could not fill TimedHitCollList");
290 return StatusCode::FAILURE;
291 } else {
292 ATH_MSG_VERBOSE(hitCollList.size() << " SiHitCollections with key " <<
293 m_inputObjectName << " found");
294 }
295
296 TimedHitCollList::iterator iColl(hitCollList.begin());
297 TimedHitCollList::iterator endColl(hitCollList.end());
298
299 for( ; iColl != endColl; ++iColl) {
300 SiHitCollection *siHitColl = new SiHitCollection(*iColl->second);
301 PileUpTimeEventIndex timeIndex(iColl->first);
302 ATH_MSG_DEBUG("SiHitCollection found with " << siHitColl->size() <<
303 " hits");
304 ATH_MSG_VERBOSE("time index info. time: " << timeIndex.time()
305 << " index: " << timeIndex.index()
306 << " type: " << timeIndex.type());
307 m_siHitCollList.push_back(siHitColl);
308 }
309
310 return StatusCode::SUCCESS;
311}
312
313
314StatusCode SiSmearedDigitizationTool::processAllSubEvents(const EventContext& ctx) {
315
316 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: in pixel processAllSubEvents() ---" );
317
318 InDet::SiClusterContainer* symSiContainer=nullptr;
319
320 if(m_SmearPixel){ // Smear Pixel
321 m_pixelClusterContainer = new InDet::PixelClusterContainer(m_pixel_ID->wafer_hash_max());
322
324 ATH_MSG_FATAL( "[ --- ] Could not create PixelClusterContainer");
325 return StatusCode::FAILURE;
326 }
327
328 // --------------------------------------
329 // PixelCluster container registration
330
331 m_pixelClusterContainer->cleanup();
332 if ((evtStore()->record(m_pixelClusterContainer, m_pixel_SiClustersName)).isFailure()) {
333 if ((evtStore()->retrieve(m_pixelClusterContainer, m_pixel_SiClustersName)).isFailure()) {
334 ATH_MSG_FATAL("[ hitproc ] Error while registering PixelCluster container");
335 return StatusCode::FAILURE;
336 }
337 }
338
339 // symlink the Pixel Container
340 // Pixel
341
342 if ((evtStore()->symLink(m_pixelClusterContainer,symSiContainer)).isFailure()) {
343 ATH_MSG_FATAL( "[ --- ] PixelClusterContainer could not be symlinked to SiClusterContainter in StoreGate !" );
344 return StatusCode::FAILURE;
345 } else {
346 ATH_MSG_INFO( "[ hitproc ] PixelClusterContainer symlinked to SiClusterContainer in StoreGate" );
347 }
348
349 }else{ // Smear SCT
350 m_sctClusterContainer = new InDet::SCT_ClusterContainer(m_sct_ID->wafer_hash_max());
351
353 ATH_MSG_FATAL( "[ --- ] Could not create SCT_ClusterContainer");
354 return StatusCode::FAILURE;
355 }
356
357 // --------------------------------------
358 // SCT_Cluster container registration
359 m_sctClusterContainer->cleanup();
360 if ((evtStore()->record(m_sctClusterContainer, m_Sct_SiClustersName)).isFailure()) {
361 ATH_MSG_FATAL("[ hitproc ] Error while registering SCT_Cluster container");
362 return StatusCode::FAILURE;
363 }
364 // symlink the SCT Container
365 // SCT
366
367 if ((evtStore()->symLink(m_sctClusterContainer,symSiContainer)).isFailure()) {
368 ATH_MSG_FATAL( "[ --- ] SCT_ClusterContainer could not be symlinked to SiClusterContainter in StoreGate !" );
369 return StatusCode::FAILURE;
370 } else {
371 ATH_MSG_DEBUG( "[ hitproc ] SCT_ClusterContainer symlinked to SiClusterContainer in StoreGate" );
372 }
373
374 }
375
376 if (retrieveTruth().isFailure()) {
377 ATH_MSG_FATAL ( "retrieveTruth() failed!" );
378 return StatusCode::FAILURE;
379 }
380
381 // get the container(s)
383
385
386 //this is a list<pair<time_t, DataLink<SCTUncompressedHitCollection> > >
387 TimedHitCollList hitCollList;
388 unsigned int numberOfSimHits(0);
389 if ( !(m_mergeSvc->retrieveSubEvtsData(m_inputObjectName, hitCollList, numberOfSimHits).isSuccess()) && hitCollList.empty() ) {
390 ATH_MSG_ERROR ( "Could not fill TimedHitCollList" );
391 return StatusCode::FAILURE;
392 } else {
393 ATH_MSG_DEBUG ( hitCollList.size() << " SiHitCollections with key " << m_inputObjectName << " found" );
394 }
395
396 // Define Hit Collection
397 TimedHitCollection<SiHit> thpcsi(numberOfSimHits);
398
399 //now merge all collections into one
400 TimedHitCollList::iterator iColl(hitCollList.begin());
401 TimedHitCollList::iterator endColl(hitCollList.end() );
402
404 // loop on the hit collections
405 while ( iColl != endColl ) {
407 if (m_HardScatterSplittingMode == 1 && m_HardScatterSplittingSkipper ) { ++iColl; continue; }
409 const SiHitCollection* p_collection(iColl->second);
410 thpcsi.insert(iColl->first, p_collection);
411 ATH_MSG_DEBUG ( "SiHitCollection found with " << p_collection->size() << " hits" );
412 ++iColl;
413 }
414
415 // Process the Hits straw by straw: get the iterator pairs for given straw
416 if(this->digitize(ctx, thpcsi).isFailure()) {
417 ATH_MSG_FATAL ( "digitize method failed!" );
418 return StatusCode::FAILURE;
419 }
420
421 if (m_merge)
422 if(this->mergeEvent(ctx).isFailure()) {
423 ATH_MSG_FATAL ( "merge method failed!" );
424 return StatusCode::FAILURE;
425 }
426
427 if (createAndStoreRIOs(ctx).isFailure()) {
428 ATH_MSG_FATAL ( "createAndStoreRIOs() failed!" );
429 return StatusCode::FAILURE;
430 }
431 else {
432 ATH_MSG_DEBUG ( "createAndStoreRIOs() succeeded" );
433 }
434
435 return StatusCode::SUCCESS;
436
437}
438
440
441
442 if(m_SmearPixel){ // Smear Pixel
443
445
447 if((evtStore()->retrieve(m_pixelPrdTruth, m_prdTruthNamePixel)).isFailure()){
448 ATH_MSG_FATAL("Could not retrieve collection " << m_prdTruthNamePixel);
449 return StatusCode::FAILURE;
450 }
451 }else{
452 if((evtStore()->record(m_pixelPrdTruth, m_prdTruthNamePixel)).isFailure()){
453 ATH_MSG_FATAL("Could not record collection " << m_prdTruthNamePixel);
454 return StatusCode::FAILURE;
455 }
456 }
457 }else{ // Smear SCT
459
461 if((evtStore()->retrieve(m_SCTPrdTruth, m_prdTruthNameSCT)).isFailure()){
462 ATH_MSG_FATAL("Could not retrieve collection " << m_prdTruthNameSCT);
463 return StatusCode::FAILURE;
464 }
465 }else{
466 if((evtStore()->record(m_SCTPrdTruth, m_prdTruthNameSCT)).isFailure()){
467 ATH_MSG_FATAL("Could not record collection " << m_prdTruthNameSCT);
468 return StatusCode::FAILURE;
469 }
470 }
471 }
472
473
474 return StatusCode::SUCCESS;
475
476}
477
478template<typename CLUSTER>
480
481 ATH_MSG_DEBUG("Truth map filling with cluster " << *cluster << " and link = " << hit->particleLink());
482 if (hit->particleLink().isValid()){
483 if (!HepMC::ignoreTruthLink(hit->particleLink(), m_vetoPileUpTruthLinks)) {
484 map->insert(std::make_pair(cluster->identify(), hit->particleLink()));
485 ATH_MSG_DEBUG("Truth map filled with cluster " << *cluster << " and link = " << hit->particleLink());
486 }
487 }else{
488 ATH_MSG_DEBUG("Particle link NOT valid!! Truth map NOT filled with cluster" << cluster << " and link = " << hit->particleLink());
489 }
490
491 return StatusCode::SUCCESS;
492}
493
494StatusCode SiSmearedDigitizationTool::mergeEvent(const EventContext& /*ctx*/){
495
496 if (m_SmearPixel)
498 else
500}
501
502template<typename CLUSTER>
504
505 // take needed information on the first clusters
506 Amg::Vector2D intersection_a = clusterA->localPosition();
507
508 // take needed information on the second clusters
509 Amg::Vector2D intersection_b = clusterB->localPosition();
510
511 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: intersection_a = " << intersection_a);
512 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: intersection_b = " << intersection_b);
513
514 double distX = intersection_a.x() - intersection_b.x();
515 double distY = intersection_a.y() - intersection_b.y();
516
517 return sqrt(distX*distX + distY*distY);
518}
519
520template<typename CLUSTER>
522 // take needed information on the first cluster
523 const Amg::MatrixX& clusterErr_a = clusterA->localCovariance();
524
525 // take needed information on the second clusters
526 const Amg::MatrixX& clusterErr_b = clusterB->localCovariance();
527
528 double sigmaX = sqrt(Amg::error(clusterErr_a,Trk::locX) * Amg::error(clusterErr_a,Trk::locX) +
529 Amg::error(clusterErr_b,Trk::locX) * Amg::error(clusterErr_b,Trk::locX));
530
531 double sigmaY = sqrt(Amg::error(clusterErr_a,Trk::locY) * Amg::error(clusterErr_a,Trk::locY) +
532 Amg::error(clusterErr_b,Trk::locY) * Amg::error(clusterErr_b,Trk::locY));
533
534 return sqrt(sigmaX*sigmaX + sigmaY*sigmaY);
535}
536
537template<typename CLUSTER>
539 // take needed information on the first clusters
540 Amg::Vector2D intersection_a = clusterA->localPosition();
541 const Amg::MatrixX& clusterErr_a = clusterA->localCovariance();
542
543 // take needed information on the second clusters
544 Amg::Vector2D intersection_b = clusterB->localPosition();
545 const Amg::MatrixX& clusterErr_b = clusterB->localCovariance();
546
547 double sigmaX = sqrt(Amg::error(clusterErr_a,Trk::locX) * Amg::error(clusterErr_a,Trk::locX) +
548 Amg::error(clusterErr_b,Trk::locX) * Amg::error(clusterErr_b,Trk::locX));
549
550 double sigmaY = sqrt(Amg::error(clusterErr_a,Trk::locY) * Amg::error(clusterErr_a,Trk::locY) +
551 Amg::error(clusterErr_b,Trk::locY) * Amg::error(clusterErr_b,Trk::locY));
552
553 double interX = 0.5*(intersection_a.x()+intersection_b.x());
554 double interY = 0.5*(intersection_a.y()+intersection_b.y());
555
556 Amg::Vector2D intersection(interX, interY);
557
558 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: intersection = " << intersection);
559
560 const InDet::SiWidth& siWidth_a = clusterA->width();
561 const InDet::SiWidth& siWidth_b = clusterB->width();
562
563 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: siWidth_a = " << siWidth_a);
564 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: siWidth_b = " << siWidth_b);
565
566 Amg::Vector2D colRow = siWidth_a.colRow() + siWidth_b.colRow();
567 Amg::Vector2D phiRz = siWidth_a.widthPhiRZ() + siWidth_b.widthPhiRZ();
568
569 InDet::SiWidth siWidth(colRow, phiRz);
570
571 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: siWidth = " << siWidth);
572
573 AmgSymMatrix(2) covariance;
574 covariance.setIdentity();
575 covariance(Trk::locX,Trk::locX) = sigmaX*sigmaX;
576 covariance(Trk::locY,Trk::locY) = sigmaY*sigmaY;
577 Amg::MatrixX clusterErr = covariance;
578
579 return ClusterInfo( intersection, siWidth, clusterErr );
580
581}
582
584{
585 // The idea is first to check how many cluster we have to merge and then merge them.
586 // The loop is done until there aren't any clusters to merge
587
588 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: in mergeClusters() using PixelClusters --- ");
589
590 Pixel_detElement_RIO_map::iterator i = cluster_map->begin();
591 Pixel_detElement_RIO_map::iterator e = cluster_map->end();
592
593 for (; i != e; i = cluster_map->upper_bound(i->first)){
594 IdentifierHash current_id = i->first;
595 // Check if clusters with current_id have been already considered
596
597 bool NewMerge = true;
598
599 while (NewMerge) {
600 NewMerge = false;
601 std::pair <Pixel_detElement_RIO_map::iterator, Pixel_detElement_RIO_map::iterator> range = cluster_map->equal_range(current_id);
602
603 for ( Pixel_detElement_RIO_map::iterator iter = range.first; iter != range.second; ++iter){
604 for (Pixel_detElement_RIO_map::iterator inner_iter = std::next(iter); inner_iter != range.second; ++inner_iter){
605
606 double dist = calculateDistance((*iter).second,(*inner_iter).second);
607 double sigma = calculateSigma((*iter).second,(*inner_iter).second);
608
609 if( dist <= m_nSigma * sigma) {
610
611 std::vector<Identifier> rdoList;
612
614 InDet::SiWidth siWidth;
615 Amg::MatrixX clusterErr;
616 std::tie( intersection, siWidth, clusterErr ) = calculateNewCluster( iter->second, inner_iter->second );
617
618 const InDetDD::SiDetectorElement* hitSiDetElement = (((*inner_iter).second)->detectorElement());
619 Identifier intersectionId = hitSiDetElement->identifierOfPosition(intersection);
620
621 rdoList.push_back(intersectionId);
622
623 InDetDD::SiCellId currentCellId = hitSiDetElement->cellIdFromIdentifier(intersectionId);
624
625 if ( !currentCellId.isValid() ) {
626 continue;
627 }
628
629 InDet::PixelCluster* pixelCluster = new InDet::PixelCluster(intersectionId,
631 std::move(rdoList),
632 siWidth,
633 hitSiDetElement,
634 std::move(clusterErr));
635 ((*inner_iter).second) = pixelCluster;
636
637 cluster_map->erase(iter);
638 NewMerge = true;
639 goto REPEAT_LOOP;
640 }
641 }
642 }
643 REPEAT_LOOP: ;
644 }
645 }
646
647 return StatusCode::SUCCESS;
648}
649
651{
652 // The idea is first to check how many cluster we have to merge and then merge them.
653 // The loop is done until there aren't any clusters to merge
654
655 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: in mergeClusters() using SCT_Clusters --- ");
656
657 SCT_detElement_RIO_map::iterator i = cluster_map->begin();
658 SCT_detElement_RIO_map::iterator e = cluster_map->end();
659
660 for (; i != e; i = cluster_map->upper_bound(i->first)){
661 IdentifierHash current_id = i->first;
662 // Check if clusters with current_id have been already considered
663 bool NewMerge = true;
664
665 while (NewMerge) {
666 NewMerge = false;
667
668 std::pair <SCT_detElement_RIO_map::iterator, SCT_detElement_RIO_map::iterator> range = cluster_map->equal_range(current_id);
669 for ( SCT_detElement_RIO_map::iterator iter = range.first; iter != range.second; ++iter){
670 for ( SCT_detElement_RIO_map::iterator inner_iter = std::next(iter); inner_iter != range.second; ++inner_iter){
671
672 double dist = calculateDistance((*iter).second,(*inner_iter).second);
673
674 double sigma = calculateSigma((*iter).second,(*inner_iter).second);
675
676 if( dist <= m_nSigma * sigma) {
677
678 std::vector<Identifier> rdoList;
679
681 InDet::SiWidth siWidth;
682 Amg::MatrixX clusterErr;
683 std::tie( intersection, siWidth, clusterErr ) = calculateNewCluster( iter->second, inner_iter->second );
684
685 const InDetDD::SiDetectorElement* hitSiDetElement = (((*inner_iter).second)->detectorElement());
686 Identifier intersectionId = hitSiDetElement->identifierOfPosition(intersection);
687
688 rdoList.push_back(intersectionId);
689
690 InDetDD::SiCellId currentCellId = hitSiDetElement->cellIdFromIdentifier(intersectionId);
691
692 if ( !currentCellId.isValid() ) {
693 continue;
694 }
695
696 InDet::SCT_Cluster* sctCluster = new InDet::SCT_Cluster(intersectionId,
698 std::vector<Identifier>(rdoList),
699 siWidth,
700 hitSiDetElement,
701 Amg::MatrixX(clusterErr));
702 ((*inner_iter).second) = sctCluster;
703
704 cluster_map->erase(iter);
705 NewMerge = true;
706 goto REPEAT_LOOP;
707 }
708 }
709 }
710 REPEAT_LOOP: ;
711 }
712 }
713
714 return StatusCode::SUCCESS;
715}
716
717StatusCode SiSmearedDigitizationTool::digitize(const EventContext& ctx,
719{
720 // Set the RNG to use for this event.
721 ATHRNG::RNGWrapper* rngWrapper = m_rndmSvc->getEngine(this, m_randomEngineName);
722 const std::string rngName = name()+m_randomEngineName;
723 rngWrapper->setSeed( rngName, ctx );
724 CLHEP::HepRandomEngine *rndmEngine = rngWrapper->getEngine(ctx);
725
726 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: in SiSmearedDigizationTool::digitize() ---" );
727
728 // Get PixelDetectorElementCollection
729 const InDetDD::SiDetectorElementCollection* elementsPixel = nullptr;
730 if (m_SmearPixel) {
732 elementsPixel = pixelDetEle.retrieve();
733 if (elementsPixel==nullptr) {
734 ATH_MSG_FATAL(m_pixelDetEleCollKey.fullKey() << " could not be retrieved");
735 return StatusCode::FAILURE;
736 }
737 }
738 // Get SCT_DetectorElementCollection
739 const InDetDD::SiDetectorElementCollection* elementsSCT = nullptr;
740 if (not m_SmearPixel) {
742 elementsSCT = sctDetEle.retrieve();
743 if (elementsSCT==nullptr) {
744 ATH_MSG_FATAL(m_SCTDetEleCollKey.fullKey() << " could not be retrieved");
745 return StatusCode::FAILURE;
746 }
747 }
748
750
751 if(m_SmearPixel) { // Smear Pixel
753 } else {// Smear SCT
755 }
756
757 while (thpcsi.nextDetectorElement(i, e)) {
758
759 while (i != e) {
760 m_useDiscSurface = false;
761
762 const TimedHitPtr<SiHit>& hit(*i++);
763 int barrelEC = hit->getBarrelEndcap();
764 int layerDisk = hit->getLayerDisk();
765 int phiModule = hit->getPhiModule();
766 int etaModule = hit->getEtaModule();
767 int side = 0;
768
769 const InDetDD::SiDetectorElement* hitSiDetElement = nullptr;
770
771 if(m_SmearPixel) { // Smear Pixel
772 Identifier wafer_id = m_pixel_ID->wafer_id(barrelEC,layerDisk,phiModule,etaModule);
773 IdentifierHash wafer_hash = m_pixel_ID->wafer_hash(wafer_id);
774 const InDetDD::SiDetectorElement* hitSiDetElement_temp = elementsPixel->getDetectorElement(wafer_hash);
775 ATH_MSG_DEBUG("Pixel SiDetectorElement --> barrel_ec " << barrelEC << ", layer_disk " << layerDisk << ", phi_module " << phiModule << ", eta_module " << etaModule );
776 hitSiDetElement = hitSiDetElement_temp;
777 } else { // Smear SCT
778 side = hit->getSide();
779 Identifier idwafer = m_sct_ID->wafer_id(barrelEC,layerDisk,phiModule,etaModule,side);
780 IdentifierHash idhash = m_sct_ID->wafer_hash(m_sct_ID->wafer_id(idwafer));
781 const InDetDD::SiDetectorElement* hitSiDetElement_temp = elementsSCT->getDetectorElement(idhash);
782 ATH_MSG_DEBUG("SCT SiDetectorElement --> barrel_ec " << barrelEC << ", layer_disk " << layerDisk << ", phi_module " << phiModule << ", eta_module " << etaModule << ", side " << side);
783 hitSiDetElement = hitSiDetElement_temp;
784 }
785
786 // untangling the logic: if not using custom geometry, hitSiDetElement
787 // gets dereferenced (and should be checked)
788 //
789 // if using custom geometry, hitPlanarDetElement gets dereferenced
790 // (and should be checked)
791 if (not hitSiDetElement) {
792 ATH_MSG_FATAL("hitSiDetElement is null in SiSmearedDigitizationTool:"<<__LINE__);
793 throw std::runtime_error(std::string("hitSiDetElement is null in SiSmearedDigitizationTool::digitize() "));
794 }
795
796 if (m_SmearPixel && !(hitSiDetElement->isPixel())) continue;
797 if (!m_SmearPixel && !(hitSiDetElement->isSCT())) continue;
798
799 IdentifierHash waferID;
800
801 if(m_SmearPixel) { // Smear Pixel
802 waferID = m_pixel_ID->wafer_hash(hitSiDetElement->identify());
803 } else { // Smear SCT
804 waferID = m_sct_ID->wafer_hash(hitSiDetElement->identify());
805 }
806
807 HepGeom::Point3D<double> pix_localStartPosition = hit->localStartPosition();
808 HepGeom::Point3D<double> pix_localEndPosition = hit->localEndPosition();
809
810 pix_localStartPosition = hitSiDetElement->hitLocalToLocal3D(pix_localStartPosition);
811 pix_localEndPosition = hitSiDetElement->hitLocalToLocal3D(pix_localEndPosition);
812
813 double localEntryX = pix_localStartPosition.x();
814 double localEntryY = pix_localStartPosition.y();
815 double localEntryZ = pix_localStartPosition.z();
816 double localExitX = pix_localEndPosition.x();
817 double localExitY = pix_localEndPosition.y();
818 double localExitZ = pix_localEndPosition.z();
819
820 double thickness = 0.0;
821 thickness = hitSiDetElement->thickness();
822
823 // Transform to reconstruction local coordinates (different x,y,z ordering and sign conventions compared to simulation coords)
824 if (!m_SmearPixel) { // Smear SCT
825 HepGeom::Point3D<double> sct_localStartPosition = hit->localStartPosition();
826 HepGeom::Point3D<double> sct_localEndPosition = hit->localEndPosition();
827
828 sct_localStartPosition = hitSiDetElement->hitLocalToLocal3D(sct_localStartPosition);
829 sct_localEndPosition = hitSiDetElement->hitLocalToLocal3D(sct_localEndPosition);
830
831 localEntryX = sct_localStartPosition.x();
832 localEntryY = sct_localStartPosition.y();
833 localEntryZ = sct_localStartPosition.z();
834 localExitX = sct_localEndPosition.x();
835 localExitY = sct_localEndPosition.y();
836 localExitZ = sct_localEndPosition.z();
837 }
838
839 double distX = std::abs(std::abs(localExitX)-std::abs(localEntryX));
840 double distY = std::abs(std::abs(localExitY)-std::abs(localEntryY));
841
842 if(m_SmearPixel) { // Smear Pixel
843 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: pixel start position --- " << localEntryX << ", " << localEntryY << ", " << localEntryZ );
844 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: pixel exit position --- " << localExitX << ", " << localExitY << ", " << localExitZ );
845 m_x_entry_pixel = localEntryX;
846 m_y_entry_pixel = localEntryY;
847 m_z_entry_pixel = localEntryZ;
848 m_x_exit_pixel = localExitX;
849 m_y_exit_pixel = localExitY;
850 m_z_exit_pixel = localExitZ;
851 } else { // Smear SCT
852 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: SCT start position --- " << localEntryX << ", " << localEntryY << ", " << localEntryZ );
853 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: SCT exit position --- " << localExitX << ", " << localExitY << ", " << localExitZ );
854 m_x_entry_SCT = localEntryX;
855 m_y_entry_SCT = localEntryY;
856 m_z_entry_SCT = localEntryZ;
857 m_x_exit_SCT = localExitX;
858 m_y_exit_SCT = localExitY;
859 m_z_exit_SCT = localExitZ;
860 }
861
862 Amg::Vector2D localEntry(localEntryX,localEntryY);
863 Amg::Vector2D localExit(localExitX,localExitY);
864
865 // the pixel positions and other needed stuff for the geometrical clustering
866 std::vector<Identifier> rdoList;
867
868 Amg::Vector3D localDirection(localExitX-localEntryX, localExitY-localEntryY, localExitZ-localEntryZ);
869
870 InDetDD::SiCellId entryCellId;
871 InDetDD::SiCellId exitCellId;
872
873 // get the identifier of the entry and the exit
874 Identifier entryId = hitSiDetElement->identifierOfPosition(localEntry);
875 Identifier exitId = hitSiDetElement->identifierOfPosition(localExit);
876
877 // now get the cellIds and check whether they're valid
878 entryCellId = hitSiDetElement->cellIdFromIdentifier(entryId);
879 exitCellId = hitSiDetElement->cellIdFromIdentifier(exitId);
880
881 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: entryId " << entryId << " --- exitId " << exitId );
882 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: entryCellId " << entryCellId << " --- exitCellId " << exitCellId );
883
884 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: surface " << hitSiDetElement->surface());
885
886 // entry / exit validity
887 bool entryValid = entryCellId.isValid();
888 bool exitValid = exitCellId.isValid();
889
890 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: entryValid? " << entryValid << " --- exitValid? " << exitValid );
891
892 if (!entryValid && !exitValid) continue;
893
894 // the intersecetion id and cellId of it
895 double interX = 0.5*(localEntryX+localExitX);
896 double interY = 0.5*(localEntryY+localExitY);
897
898 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: " << (m_SmearPixel ? "pixel" : "SCT") << " inter X --- " << interX );
899 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: " << (m_SmearPixel ? "pixel" : "SCT") << " inter Y --- " << interY );
900 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: " << (m_SmearPixel ? "pixel" : "SCT") << " dist X --- " << distX );
901 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: " << (m_SmearPixel ? "pixel" : "SCT") << " dist Y --- " << distY );
902
903 //the particle crosses at least n pixels (in x direction you have timesX, in y direction you have timesY)
904 double timesX = (m_pitch_X) ? floor(distX/m_pitch_X) : 0;
905 double timesY = (m_pitch_Y) ? floor(distY/m_pitch_Y) : 0;
906
907 double newdistX = distX - (timesX*m_pitch_X);
908 double newdistY = distY - (timesY*m_pitch_Y);
909
910 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: times X --- " << timesX );
911 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: times Y --- " << timesY );
912 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: new dist X --- " << newdistX );
913 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: new dist Y --- " << newdistY );
914 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: thickness --- " << thickness );
915
916 //Probability
917
918 double ProbY = 2*newdistY/(m_pitch_Y+newdistY);
919 double ProbX = 2*newdistX/(m_pitch_X+newdistX);
920
921 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: ProbX --- " << ProbX );
922 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: ProbY --- " << ProbY );
923
924 // create the errors
925 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: pitch X --- " << m_pitch_X );
926 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: pitch Y --- " << m_pitch_Y );
927
928 double sigmaX = m_pitch_X/sqrt(12.);
929 double sigmaY = m_pitch_Y/sqrt(12.);
930
931 int elementX = timesX+1;
932 int elementY = timesY+1;
933
934 if(m_SmearPixel) {
935 if (CLHEP::RandFlat::shoot(rndmEngine, 0.0, 1.0) < ProbY) { // number of crossed pixel is (timesY+1)+1
936 sigmaY = (double)(timesY+2)*m_pitch_Y/sqrt(12.);
937 elementY++;
938 } else // number of crossed pixel is (timesY+1)
939 sigmaY = (double)(timesY+1)*m_pitch_Y/sqrt(12.);
940
941 if (CLHEP::RandFlat::shoot(rndmEngine, 0.0, 1.0) < ProbX) { // number of crossed pixel is (timesY+1)+1
942 sigmaX = (double)(timesX+2)*m_pitch_X/sqrt(12.);
943 elementX++;
944 } else // number of crossed pixel is (timesY+1)
945 sigmaX = (double)(timesX+1)*m_pitch_X/sqrt(12.);
946 }
947
948 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool " << (m_SmearPixel ? "pixel" : "SCT") << " sigma X --- " << sigmaX);
949 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool " << (m_SmearPixel ? "pixel" : "SCT") << " sigma Y --- " << sigmaY);
950
951
952 double temp_X = interX;
953 double temp_Y = interY;
954
955 Amg::Vector2D intersection(interX, interY);
956
957 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: Intersection after smearing: " << intersection);
958
959 Identifier intersectionId;
960 intersectionId = hitSiDetElement->identifierOfPosition(intersection);
961
962 rdoList.push_back(intersectionId);
963 InDetDD::SiCellId currentCellId = hitSiDetElement->cellIdFromIdentifier(intersectionId);
964
965 if (!currentCellId.isValid()) continue;
966
967 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: Intersection Id = " << intersectionId << " --- currentCellId = " << currentCellId );
968
969 if(m_SmearPixel) { // Smear Pixel --> Create Pixel Cluster
970
971 double lengthX = (m_pitch_X) ? elementX*m_pitch_X : 1.;
972 double lengthY = (m_pitch_Y) ? elementY*m_pitch_Y : 1.;
973
974 if (m_pitch_X == 0. || m_pitch_Y == 0.)
975 ATH_MSG_WARNING( "--- SiSmearedDigitizationTool: pitchX and/or pitchY are 0. Cluster length is forced to be 1. mm");
976
977 InDet::SiWidth siWidth(Amg::Vector2D(elementX,elementY), Amg::Vector2D(lengthX, lengthY));
978
979 InDet::PixelCluster* pixelCluster = nullptr;
980
981 AmgSymMatrix(2) covariance;
982 covariance.setIdentity();
983 covariance(Trk::locX,Trk::locX) = sigmaX*sigmaX;
984 covariance(Trk::locY,Trk::locY) = sigmaY*sigmaY;
985
986 // create the cluster
987 pixelCluster = new InDet::PixelCluster(intersectionId,
989 std::vector<Identifier>(rdoList),
990 siWidth,
991 hitSiDetElement,
992 Amg::MatrixX(covariance));
993 m_pixelClusterMap->insert(std::pair<IdentifierHash, InDet::PixelCluster* >(waferID, pixelCluster));
994
995 if (FillTruthMap(m_pixelPrdTruth, pixelCluster, hit).isFailure()) {
996 ATH_MSG_FATAL ( "FillTruthMap() for pixel failed!" );
997 return StatusCode::FAILURE;
998 }
999
1000 if (m_checkSmear) {
1001
1002 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: pixelCluster --> " << *pixelCluster);
1003 //Take info to store in the tree
1004 m_x_pixel = temp_X;
1005 m_y_pixel = temp_Y;
1006
1007 m_x_pixel_smeared = (pixelCluster->localPosition()).x();
1008 m_y_pixel_smeared = (pixelCluster->localPosition()).y();
1009
1010 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: BEFORE SMEARING LocalPosition --> X = " << m_x_pixel << " Y = " << m_y_pixel);
1011 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: LocalPosition --> X = " << m_x_pixel_smeared << " Y = " << m_y_pixel_smeared);
1012
1013 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: GlobalPosition --> X = " << (pixelCluster->globalPosition()).x() << " Y = " << (pixelCluster->globalPosition()).y());
1014 m_x_pixel_global = (pixelCluster->globalPosition()).x();
1015 m_y_pixel_global = (pixelCluster->globalPosition()).y();
1016 m_z_pixel_global = (pixelCluster->globalPosition()).z();
1017
1018 m_Err_x_pixel = Amg::error(pixelCluster->localCovariance(), Trk::locX);
1019 m_Err_y_pixel = Amg::error(pixelCluster->localCovariance(), Trk::locY);
1020
1021 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: Error --> X = " << m_Err_x_pixel << " Y = " << m_Err_y_pixel);
1022
1023 m_currentTree -> Fill();
1024 } // End of Pix smear check
1025
1026 } else { // Smear SCT --> Create SCT Cluster
1027
1028 // prepare the clusters
1029 InDet::SCT_Cluster * sctCluster = nullptr;
1030
1031 // Pixel Design needed -------------------------------------------------------------
1032 const InDetDD::SCT_ModuleSideDesign* design_sct;
1033
1034 design_sct = dynamic_cast<const InDetDD::SCT_ModuleSideDesign*>(&hitSiDetElement->design());
1035
1036 if (!design_sct) {
1037 ATH_MSG_INFO ( "Could not get design"<< design_sct) ;
1038 continue;
1039 }
1040
1041 // Find length of strip at centre
1042 double clusterWidth = rdoList.size()*hitSiDetElement->phiPitch(intersection);
1043 const std::pair<InDetDD::SiLocalPosition, InDetDD::SiLocalPosition> ends(design_sct->endsOfStrip(intersection));
1044 double stripLength = std::abs(ends.first.xEta()-ends.second.xEta());
1045
1046 InDet::SiWidth siWidth(Amg::Vector2D(int(rdoList.size()),1),
1047 Amg::Vector2D(clusterWidth,stripLength) );
1048
1049 const Amg::Vector2D& colRow = siWidth.colRow();
1050
1051 AmgSymMatrix(2) mat;
1052 mat.setIdentity();
1053 mat(Trk::locX,Trk::locX) = sigmaX*sigmaX;
1054 mat(Trk::locY,Trk::locY) = hitSiDetElement->length()*hitSiDetElement->length()/12.;
1055
1056 // single strip - resolution close to pitch/sqrt(12)
1057 // two-strip hits: better resolution, approx. 40% lower
1058
1059 InDetDD::DetectorShape elShape = hitSiDetElement->design().shape();
1060 if(m_emulateAtlas && elShape == InDetDD::Trapezoid)
1061 { // rotation for endcap SCT
1062
1063 if(colRow.x() == 1) {
1064 mat(Trk::locX,Trk::locX) = pow(siWidth.phiR(),2)/12;
1065 }
1066 else if(colRow.x() == 2) {
1067 mat(Trk::locX,Trk::locX) = pow(0.27*siWidth.phiR(),2)/12;
1068 }
1069 else {
1070 mat(Trk::locX,Trk::locX) = pow(siWidth.phiR(),2)/12;
1071 }
1072
1073 mat(Trk::locY,Trk::locY) = pow(siWidth.z()/colRow.y(),2)/12;
1074 double sn = hitSiDetElement->sinStereoLocal(intersection);
1075 double sn2 = sn*sn;
1076 double cs2 = 1.-sn2;
1077 double w = hitSiDetElement->phiPitch(intersection)/hitSiDetElement->phiPitch();
1078 double v0 = mat(Trk::locX, Trk::locX)*w*w;
1079 double v1 = mat(Trk::locY, Trk::locY);
1080 mat(Trk::locX, Trk::locX) = (cs2*v0+sn2*v1);
1081 mat(Trk::locY, Trk::locX) = (sn*sqrt(cs2)*(v0-v1));
1082 mat(Trk::locY, Trk::locY) = (sn2*v0+cs2*v1);
1083 } // End of rotation endcap SCT
1084
1085
1086 sctCluster = new InDet::SCT_Cluster(intersectionId,
1088 std::vector<Identifier>(rdoList),
1089 siWidth,
1090 hitSiDetElement,
1091 Amg::MatrixX(mat));
1092
1093 m_sctClusterMap->insert(std::pair<IdentifierHash, InDet::SCT_Cluster* >(waferID, sctCluster));
1094
1095 if (FillTruthMap(m_SCTPrdTruth, sctCluster, hit).isFailure()) {
1096 ATH_MSG_FATAL ( "FillTruthMap() for SCT failed!" );
1097 return StatusCode::FAILURE;
1098 }
1099
1100 if (m_checkSmear) {
1101
1102 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: SCT_Cluster --> " << *sctCluster);
1103
1104 //Take info to store in the tree
1105 m_x_SCT = m_useDiscSurface ? temp_Y : temp_X;
1106 m_x_SCT_smeared = (sctCluster->localPosition()).x();
1107
1108 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: BEFORE SMEARING LocalPosition --> X = " << m_x_SCT );
1109 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: LocalPosition --> X = " << m_x_SCT_smeared );
1110
1111 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: GlobalPosition --> X = " << (sctCluster->globalPosition()).x() << " Y = " << (sctCluster->globalPosition()).y());
1112 m_x_SCT_global = (sctCluster->globalPosition()).x();
1113 m_y_SCT_global = (sctCluster->globalPosition()).y();
1114 m_z_SCT_global = (sctCluster->globalPosition()).z();
1115
1116 m_currentTree -> Fill();
1117 } // End smear checking
1118 } // End SCT
1119 } // while
1120 } //while
1121 return StatusCode::SUCCESS;
1122}
1123
1124
1125StatusCode SiSmearedDigitizationTool::createAndStoreRIOs(const EventContext& ctx)
1126{
1127 // Get PixelDetectorElementCollection
1128 const InDetDD::SiDetectorElementCollection* elementsPixel = nullptr;
1129 if (m_SmearPixel) {
1131 elementsPixel = pixelDetEle.retrieve();
1132 if (elementsPixel==nullptr) {
1133 ATH_MSG_FATAL(m_pixelDetEleCollKey.fullKey() << " could not be retrieved");
1134 return StatusCode::FAILURE;
1135 }
1136 }
1137 // Get SCT_DetectorElementCollection
1138 const InDetDD::SiDetectorElementCollection* elementsSCT = nullptr;
1139 if (not m_SmearPixel) {
1141 elementsSCT = sctDetEle.retrieve();
1142 if (elementsSCT==nullptr) {
1143 ATH_MSG_FATAL(m_SCTDetEleCollKey.fullKey() << " could not be retrieved");
1144 return StatusCode::FAILURE;
1145 }
1146 }
1147
1148 if ( m_SmearPixel ) { // Store Pixel RIOs
1149
1150 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: in pixel createAndStoreRIOs() ---" );
1151
1152 Pixel_detElement_RIO_map::iterator i = m_pixelClusterMap->begin();
1153 Pixel_detElement_RIO_map::iterator e = m_pixelClusterMap->end();
1154
1155 for ( ; i != e; i = m_pixelClusterMap->upper_bound(i->first) ) {
1156
1157 std::pair <Pixel_detElement_RIO_map::iterator, Pixel_detElement_RIO_map::iterator> range;
1158 range = m_pixelClusterMap->equal_range(i->first);
1159
1160 Pixel_detElement_RIO_map::iterator firstDetElem;
1161 firstDetElem = range.first;
1162
1163 IdentifierHash waferID;
1164 waferID = firstDetElem->first;
1165
1166 const InDetDD::SiDetectorElement* detElement = elementsPixel->getDetectorElement(waferID);
1167
1168 InDet::PixelClusterCollection *clusterCollection = new InDet::PixelClusterCollection(waferID);
1169 clusterCollection->setIdentifier(detElement->identify());
1170
1171 for ( Pixel_detElement_RIO_map::iterator iter = range.first; iter != range.second; ++iter ) {
1172
1173 InDet::PixelCluster* pixelCluster = (*iter).second;
1174 pixelCluster->setHashAndIndex(clusterCollection->identifyHash(),clusterCollection->size());
1175 clusterCollection->push_back(pixelCluster);
1176 }
1177
1178 if ( m_pixelClusterContainer->addCollection( clusterCollection, waferID ).isFailure() ) {
1179 ATH_MSG_WARNING( "Could not add collection to Identifyable container !" );
1180 }
1181 } // end for
1182
1183 m_pixelClusterMap->clear();
1184
1185 }
1186 else { // Store SCT RIOs
1187
1188 ATH_MSG_DEBUG( "--- SiSmearedDigitizationTool: in SCT createAndStoreRIOs() ---" );
1189
1190 SCT_detElement_RIO_map::iterator i = m_sctClusterMap->begin();
1191 SCT_detElement_RIO_map::iterator e = m_sctClusterMap->end();
1192
1193 for ( ; i != e; i = m_sctClusterMap->upper_bound(i->first) ) {
1194 std::pair <SCT_detElement_RIO_map::iterator, SCT_detElement_RIO_map::iterator> range;
1195 range = m_sctClusterMap->equal_range(i->first);
1196
1197 SCT_detElement_RIO_map::iterator firstDetElem;
1198 firstDetElem = range.first;
1199
1200 IdentifierHash waferID;
1201 waferID = firstDetElem->first;
1202 const InDetDD::SiDetectorElement* detElement = elementsSCT->getDetectorElement(waferID);
1203
1204 InDet::SCT_ClusterCollection *clusterCollection = new InDet::SCT_ClusterCollection(waferID);
1205 clusterCollection->setIdentifier(detElement->identify());
1206
1207
1208 for ( SCT_detElement_RIO_map::iterator iter = range.first; iter != range.second; ++iter ) {
1209 InDet::SCT_Cluster* sctCluster = (*iter).second;
1210 sctCluster->setHashAndIndex(clusterCollection->identifyHash(),clusterCollection->size());
1211 clusterCollection->push_back(sctCluster);
1212 }
1213
1214 if ( m_sctClusterContainer->addCollection( clusterCollection, clusterCollection->identifyHash() ).isFailure() ) {
1215 ATH_MSG_WARNING( "Could not add collection to Identifyable container !" );
1216 }
1217
1218 } // end for
1219
1220 m_sctClusterMap->clear();
1221 }
1222
1223 return StatusCode::SUCCESS;
1224}
#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)
xAOD::CaloCluster CLUSTER
defines an "iterator" over instances of a given type in StoreGateSvc
#define AmgSymMatrix(dim)
std::vector< xAOD::EventInfo::SubEvent >::const_iterator SubEventIterator
Definition IPileUpTool.h:22
the preferred mechanism to access information from the different event stores in a pileup job.
This is an Identifier helper class for the Pixel subdetector.
This is an Identifier helper class for the SCT subdetector.
AtlasHitsVector< SiHit > SiHitCollection
std::tuple< Amg::Vector2D, InDet::SiWidth, Amg::MatrixX > ClusterInfo
constexpr int pow(int base, int exp) noexcept
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
This is a "hash" representation of an Identifier.
virtual DetectorShape shape() const
Shape of element.
Base class for the SCT module side design, extended by the Forward and Barrel module design.
virtual std::pair< SiLocalPosition, SiLocalPosition > endsOfStrip(const SiLocalPosition &position) const override=0
give the ends of strips
Identifier for the strip or pixel cell.
Definition SiCellId.h:29
bool isValid() const
Test if its in a valid state.
Definition SiCellId.h:136
Class to hold the SiDetectorElement objects to be put in the detector store.
const SiDetectorElement * getDetectorElement(const IdentifierHash &hash) const
Class to hold geometrical description of a silicon detector element.
virtual SiCellId cellIdFromIdentifier(const Identifier &identifier) const override final
SiCellId from Identifier.
virtual const SiDetectorDesign & design() const override final
access to the local description (inline):
double phiPitch() const
Pitch (inline methods)
double sinStereoLocal(const Amg::Vector2D &localPos) const
Angle of strip in local frame with respect to the etaAxis.
double length() const
Length in eta direction (z - barrel, r - endcap)
HepGeom::Point3D< double > hitLocalToLocal3D(const HepGeom::Point3D< double > &hitPosition) const
Same as previuos method but 3D.
virtual Identifier identify() const override final
identifier of this detector element (inline)
Identifier identifierOfPosition(const Amg::Vector2D &localPos) const
Full identifier of the cell for a given position: assumes a raw local position (no Lorentz shift)
Trk::Surface & surface()
Element Surface.
const Amg::Vector3D & globalPosition() const
return global position reference
double z() const
Definition SiWidth.h:131
double phiR() const
Definition SiWidth.h:126
const Amg::Vector2D & widthPhiRZ() const
Definition SiWidth.h:121
const Amg::Vector2D & colRow() const
Definition SiWidth.h:115
A PRD is mapped onto all contributing particles.
Gaudi::Property< int > m_vetoPileUpTruthLinks
PileUpToolBase(const std::string &type, const std::string &name, const IInterface *parent)
const_pointer_type retrieve()
std::string m_randomEngineName
Name of the random number stream.
const SCT_ID * m_sct_ID
Handle to the ID helper.
InDet::SCT_ClusterContainer * m_sctClusterContainer
the SCT_ClusterContainer
double calculateDistance(CLUSTER *clusterA, CLUSTER *clusterB)
StatusCode mergeClusters(Pixel_detElement_RIO_map *cluster_map)
StatusCode createAndStoreRIOs(const EventContext &ctx)
TTree * m_currentTree
the tree to store information from pixel and SCT (before and after smearing)
SG::ReadCondHandleKey< InDetDD::SiDetectorElementCollection > m_pixelDetEleCollKey
SG::ReadCondHandleKey< InDetDD::SiDetectorElementCollection > m_SCTDetEleCollKey
SCT_detElement_RIO_map * m_sctClusterMap
std::multimap< IdentifierHash, InDet::SCT_Cluster * > SCT_detElement_RIO_map
const PixelID * m_pixel_ID
Handle to the ID helper.
StatusCode digitize(const EventContext &ctx, TimedHitCollection< SiHit > &thpcsi)
PRD_MultiTruthCollection * m_pixelPrdTruth
StatusCode processAllSubEvents(const EventContext &ctx)
ServiceHandle< PileUpMergeSvc > m_mergeSvc
PileUp Merge service.
std::multimap< IdentifierHash, InDet::PixelCluster * > Pixel_detElement_RIO_map
StatusCode processBunchXing(int bunchXing, SubEventIterator bSubEvents, SubEventIterator eSubEvents)
StatusCode FillTruthMap(PRD_MultiTruthCollection *, CLUSTER *, const TimedHitPtr< SiHit > &)
ServiceHandle< IAthRNGSvc > m_rndmSvc
Random number service.
ServiceHandle< ITHistSvc > m_thistSvc
ClusterInfo calculateNewCluster(CLUSTER *clusterA, CLUSTER *clusterB)
PRD_MultiTruthCollection * m_SCTPrdTruth
InDet::PixelClusterContainer * m_pixelClusterContainer
the PixelClusterContainer
std::vector< SiHitCollection * > m_siHitCollList
name of the sub event hit collections.
int m_HardScatterSplittingMode
Process all SiHit or just those from signal or background events.
StatusCode prepareEvent(const EventContext &ctx, unsigned int)
Pixel_detElement_RIO_map * m_pixelClusterMap
StatusCode mergeEvent(const EventContext &ctx)
double calculateSigma(CLUSTER *clusterA, CLUSTER *clusterB)
bool nextDetectorElement(const_iterator &b, const_iterator &e)
sets an iterator range with the hits of current detector element returns a bool when done
TimedVector::const_iterator const_iterator
void insert(const PileUpTimeEventIndex &timeEventIndex, const AtlasHitsVector< HIT > *inputCollection)
a smart pointer to a hit that also provides access to the extended timing info of the host event.
Definition TimedHitPtr.h:18
const Amg::Vector2D & localPosition() const
return the local position reference
void setHashAndIndex(unsigned short collHash, unsigned short objIndex)
TEMP for testing: might make some classes friends later ...
STL class.
std::vector< std::string > intersection(std::vector< std::string > &v1, std::vector< std::string > &v2)
bool contains(const std::string &s, const std::string &regx)
does a string contain the substring
Definition hcg.cxx:114
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
double error(const Amg::MatrixX &mat, int index)
return diagonal error of the matrix caller should ensure the matrix is symmetric and the index is in ...
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D
bool ignoreTruthLink(const T &p, bool vetoPileUp)
Helper function for SDO creation in PileUpTools.
Message Stream Member.
@ locY
local cartesian
Definition ParamDefs.h:38
@ locX
Definition ParamDefs.h:37
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