ATLAS Offline Software
MaterialMapping.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // MaterialMapping.cxx, (c) ATLAS Detector software
8 
9 // Gaudi Units
10 #include "GaudiKernel/SystemOfUnits.h"
11 //TrkDetDescr Algs, Interfaces, Utils
12 #include "MaterialMapping.h"
16 // TrkGeometry
27 // TrkEvent
29 // Amg
31 
32 #ifdef TRKDETDESCR_MEMUSAGE
33 #include <unistd.h>
34 #endif
35 
36 Trk::MaterialMapping::MaterialMapping(const std::string& name, ISvcLocator* pSvcLocator)
37 : AthAlgorithm(name,pSvcLocator),
38  m_checkForEmptyHits(true),
39  m_mappingVolumeName("Atlas"),
40  m_mappingVolume(nullptr),
41  m_inputMaterialStepCollection("MaterialStepRecords"),
42  m_etaCutOff(6.0),
43  m_etaSide(0),
44  m_useLayerThickness(false),
45  m_associationType(1),
46  m_mapMaterial(true),
47  m_mapComposition(false),
48  m_minCompositionFraction(0.005),
49  m_elementTable(nullptr),
50  m_inputEventElementTable("ElementTable"),
51  m_accumulatedMaterialXX0(0.),
52  m_accumulatedRhoS(0.),
53  m_mapped(0),
54  m_unmapped(0),
55  m_skippedOutside(0),
56  m_layerMaterialScreenOutput(0)
57 #ifdef TRKDETDESCR_MEMUSAGE
58  ,m_memoryLogger()
59 #endif
60 {
61  // the name of the volume to map
62  declareProperty("MappingVolumeName" , m_mappingVolumeName);
63  // the extrapolation engine
64  declareProperty("CheckForEmptyHits" , m_checkForEmptyHits);
65  // general steering
66  declareProperty("EtaCutOff" , m_etaCutOff);
67  declareProperty("EtaSide" , m_etaSide);
68  // the toolhandle of the MaterialMapper to be used
69  declareProperty("MapMaterial" , m_mapMaterial);
70  // Composition related parameters
71  declareProperty("MapComposition" , m_mapComposition);
72  declareProperty("MinCompositionFraction" , m_minCompositionFraction);
73  // Steer the layer thickness
74  declareProperty("UseActualLayerThicknesss" , m_useLayerThickness);
75  // some job setup
76  declareProperty("MaterialAssociationType" , m_associationType);
77  declareProperty("InputMaterialStepCollection" , m_inputMaterialStepCollection);
78  declareProperty("InputElementTable" , m_inputEventElementTable);
79  // Output screen stuff
80  declareProperty("MaterialScreenOutputLevel" , m_layerMaterialScreenOutput);
81 
82 }
83 
85 = default;
86 
87 
89 {
90 
91  ATH_MSG_INFO("initialize()");
92 
93  ATH_CHECK( m_trackingGeometryReadKey.initialize(!m_trackingGeometryReadKey.key().empty()) );
94 
95  ATH_CHECK(m_extrapolationEngine.retrieve());
96 
97  if ( !m_materialMapper.empty() )
98  ATH_CHECK( m_materialMapper.retrieve() );
99 
100  if ( !m_layerMaterialRecordAnalyser.empty() )
101  ATH_CHECK( m_layerMaterialRecordAnalyser.retrieve() );
102 
103  if ( !m_layerMaterialCreators.empty() )
104  ATH_CHECK( m_layerMaterialCreators.retrieve() );
105 
106  if ( !m_layerMaterialAnalysers.empty() )
107  ATH_CHECK( m_layerMaterialAnalysers.retrieve() );
108 
109  ATH_CHECK( m_inputMaterialStepCollection.initialize() );
110  ATH_CHECK( m_inputEventElementTable.initialize() );
111 
112  return StatusCode::SUCCESS;
113 }
114 
115 
117 {
118  ATH_MSG_VERBOSE("MaterialMapping execute() start");
119 
120  // ------------------------------- get the trackingGeometry at first place
121  if (!m_mappingVolume) {
122  StatusCode retrieveCode = handleTrackingGeometry();
123  if (retrieveCode.isFailure()){
124  ATH_MSG_INFO("Could not retrieve mapping volume from tracking geometry. Exiting.");
125  return retrieveCode;
126  }
127  }
128  if (m_mappingVolume)
129  ATH_MSG_VERBOSE("Mapping volume correctly retrieved from tracking geometry");
130 
131 
132  SG::ReadHandle<MaterialStepCollection> materialStepCollection(m_inputMaterialStepCollection);
133 
134  // --------- prepare the element table ---------------------------------------------------
135 
136  if (m_mapComposition) {
137  SG::ReadHandle<Trk::ElementTable> eTableEvent(m_inputEventElementTable);
138  (*m_elementTable) += (*eTableEvent); // accummulate the table
139  }
140 
141 
142  // event parameters - associated asteps, and layers hit per event
143  int associatedSteps = 0;
144  m_accumulatedMaterialXX0 = 0.;
145  m_accumulatedRhoS = 0.;
146  m_layersRecordedPerEvent.clear();
147  // clearing the recorded layers per event
148  if (materialStepCollection.isValid() && !materialStepCollection->empty()){
149 
150  // get the number of material steps
151  size_t materialSteps = materialStepCollection->size();
152  ATH_MSG_DEBUG("[+] Successfully read "<< materialSteps << " geantino steps");
153 
154  // create a direction out of the last material step
155  double dirx = (*materialStepCollection)[materialSteps-1]->hitX();
156  double diry = (*materialStepCollection)[materialSteps-1]->hitY();
157  double dirz = (*materialStepCollection)[materialSteps-1]->hitZ();
158  Amg::Vector3D direction = Amg::Vector3D(dirx,diry,dirz).unit();
159 
160  double eta = direction.eta();
161  // skip the event if the eta cut is not met
162  if ( fabs(eta) > m_etaCutOff || (m_etaSide && m_etaSide*eta < 0.) ) {
163  ATH_MSG_VERBOSE("[-] Event is outside eta acceptance of " << m_etaCutOff << ". Skipping it.");
164  return StatusCode::SUCCESS;
165  }
166 
167  // now propagate through the full detector and collect the layers
168  Trk::NeutralCurvilinearParameters ncP(Amg::Vector3D(0.,0.,0.), direction, 0.);
169  // create a neutral extrapolation cell
171  ecc.navigationCurvilinear = false;
175 
176  // let's extrapolate through the detector and remember which layers (with material) should have been hit
177  std::vector< std::pair<const Trk::Layer*, Amg::Vector3D> > layersAndHits;
178  // call the extrapolation engine
179  Trk::ExtrapolationCode eCode = m_extrapolationEngine->extrapolate(ecc);
180  // end the parameters if there
181  if (eCode.isSuccess()){
182  // name of passive surfaces found
183  size_t nLayersHit = ecc.extrapolationSteps.size();
184  ATH_MSG_VERBOSE("[+] Extrapolation to layers did succeed and found " << nLayersHit << " layers.");
185  // reserve the size of the vectors
186  layersAndHits.reserve(nLayersHit);
187  // for screen output
188  size_t ilayer = 0;
189  // find all the intersected material - remember the last parameters
190  const Trk::NeutralParameters* parameters = nullptr;
191  // loop over the collected information
192  for (auto& es : ecc.extrapolationSteps){
193  // continue if we have parameters
194  parameters = es.parameters;
195  if (parameters){
196  const Trk::Surface& pSurface = parameters->associatedSurface();
197  // get the surface with associated layer (that has material)
198  ATH_MSG_VERBOSE("[L] Testing layer with associatedLayer() " << pSurface.associatedLayer() << " and materialLayer() " << pSurface.materialLayer() );
199  //
200  if ( ( pSurface.associatedLayer() && pSurface.associatedLayer()->layerMaterialProperties() ) || pSurface.materialLayer() ){
201  // material layer
202  const Trk::Layer* mLayer = pSurface.materialLayer() ? pSurface.materialLayer() : pSurface.associatedLayer();
203  // record that one
204  std::pair<const Trk::Layer*, Amg::Vector3D> layerHitPair(mLayer, parameters->position());
205  ATH_MSG_VERBOSE("[L] Layer " << ++ilayer << " with index " << mLayer->layerIndex().value() << " hit at " << Amg::toString(parameters->position()));
206  layersAndHits.push_back(layerHitPair);
207  }
208  delete parameters;
209  }
210  }
211  // cleanup of the final hits
212  if (ecc.endParameters != parameters) delete ecc.endParameters;
213 
214  // we have no layers and Hits
215  if (layersAndHits.empty()){
216  ATH_MSG_VERBOSE("[!] No Layer was intersected - skipping.");
217  return StatusCode::SUCCESS;
218  }
219 
220  // layers are ordered, hence you can move the starting point along
221  size_t currentLayer = 0;
222  // loop through hits and find the closest layer, the start point moves outwards as we go
223  for ( const Trk::MaterialStep* step : *materialStepCollection ) {
224  // verbose output
225  ATH_MSG_VERBOSE("[L] starting from layer " << currentLayer << " from layer collection for this step.");
226  // step length and position
227  double t = step->steplength();
228  Amg::Vector3D pos(step->hitX(), step->hitY(), step->hitZ());
229  // skip if :
230  // -- 0) no mapping volume exists
231  // -- 1) outside the mapping volume
232  // -- 2) outside the eta acceptance
233  if (!m_mappingVolume || !(m_mappingVolume->inside(pos)) || fabs(pos.eta()) > m_etaCutOff ){
234  ++m_skippedOutside;
235  continue;
236  }
237  // now find the closest layer
238  // (a) if the currentLayer is the last layer and the hit is still inside -> assign
239  if (currentLayer < nLayersHit-1) {
240  // search through the layers - this is the reference distance for projection
241  double currentDistance = (pos-layersAndHits[currentLayer].second).mag();
242  ATH_MSG_VERBOSE("- current distance is " << currentDistance << " from " << Amg::toString(pos) << " and " << Amg::toString(layersAndHits[currentLayer].second) );
243  for (size_t testLayer = (currentLayer+1); testLayer < nLayersHit; ++testLayer){
244  // calculate teh distance to the testLayer
245  double testDistance = (pos-layersAndHits[testLayer].second).mag();
246  ATH_MSG_VERBOSE("[L] Testing layer " << testLayer << " from layer collection for this step.");
247  ATH_MSG_VERBOSE("- test distance is " << testDistance << " from " << Amg::toString(pos) << " and " << Amg::toString(layersAndHits[testLayer].second) );
248  if ( testDistance < currentDistance ){
249  // screen output
250  ATH_MSG_VERBOSE("[L] Skipping over to current layer " << testLayer << " because " << testDistance << " < " << currentDistance);
251  // the test distance did shrink - update currentLayer
252  currentLayer = testLayer;
253  currentDistance = testDistance;
254  } else {
255  // stick to the layer you have
256  break;
257  }
258  }
259  }
260  // the currentLayer *should* be correct now
261  const Trk::Layer* assignedLayer = layersAndHits[currentLayer].first;
262  Amg::Vector3D assignedPosition = layersAndHits[currentLayer].second;
263  // associate the hit
264  // (1) count it
265  ++associatedSteps;
266  // (2) associate it
267  associateHit(*assignedLayer, pos, assignedPosition, t, step->fullMaterial());
268  } // loop over material Steps
269 
270  // check for the empty hits - they need to be taken into account
271  ATH_MSG_VERBOSE("Found " << layersAndHits.size() << " intersected layers - while having " << m_layersRecordedPerEvent.size() << " recorded ones.");
272 
273  // now - cross-chek if you have additional layers
274  for ( auto& lhp : layersAndHits){
275  // check if you find the layer int he already done record-map : not found - we need to do an empty hit scaling
276  if (m_layersRecordedPerEvent.find(lhp.first) == m_layersRecordedPerEvent.end()){
277  // try to find the layer material record
278  auto clIter = m_layerRecords.find(lhp.first);
279  if (clIter != m_layerRecords.end() ){
280  (*clIter).second.associateEmptyHit(lhp.second);
281  ATH_MSG_VERBOSE("- to layer with index "<< lhp.first->layerIndex().value() << " with empty hit detected.");
282  } else
283  ATH_MSG_WARNING("- no Layer found in the associated map! Should not happen.");
284  }
285  }
286 
287  // check whether the event was good for at least one hit
288  if (associatedSteps) {
289  ATH_MSG_VERBOSE("There are associated steps, need to call finalizeEvent() & record to the MaterialMapper.");
290  // finalize the event --------------------- Layers ---------------------------------------------
291  for (auto& lRecord : m_layerRecords ) {
292  // associated material
293  Trk::AssociatedMaterial* assMatHit = lRecord.second.finalizeEvent((*lRecord.first));
294  // record the full layer hit
295  if (assMatHit && !m_materialMapper.empty()) m_materialMapper->recordLayerHit(*assMatHit, true);
296  delete assMatHit;
297  // call the material mapper finalize method
298  ATH_MSG_VERBOSE("Calling finalizeEvent on the MaterialMapper ...");
299  }
300  } // the event had at least one associated hit
301 
302  } // end of eCode.success : needed for new mapping schema
303 
304  } // material steps existed
305 
306 
307  return StatusCode::SUCCESS;
308 }
309 
310 
311 bool Trk::MaterialMapping::associateHit( const Trk::Layer& associatedLayer,
312  const Amg::Vector3D& pos,
313  const Amg::Vector3D& positionOnLayer,
314  double stepl,
315  const Trk::Material& mat)
316 {
317  // get the intersection with the layer
318  ++m_mapped;
319  // get the layer
320  const Trk::Layer* layer = &associatedLayer;
321 
322  // get the associated volume
323  const Trk::TrackingVolume* associatedVolume = trackingGeometry().lowestTrackingVolume(pos);
324 
325  // try to find the layer material record
326  auto clIter = m_layerRecords.find(layer);
327  if (clIter != m_layerRecords.end() ){
328  // remember that you actually hit this layer
329  m_layersRecordedPerEvent[layer] = true;
330  // LayerMaterialRecord found, add the hit
331  (*clIter).second.associateGeantinoHit(positionOnLayer, stepl, mat);
332  ATH_MSG_VERBOSE("- associate Geantino Information at intersection [r/z] = " << positionOnLayer.perp() << "/"<< positionOnLayer.z() );
333  ATH_MSG_VERBOSE(" mapping distance = " << (pos-positionOnLayer).mag() );
334  ATH_MSG_VERBOSE("- associate Geantino Information ( s, s/x0 , x0 , l0, a, z, rho ) = "
335  << stepl << ", "<< stepl/mat.X0 << ", "<< mat.X0 << ", "<< mat.L0 << ", "<< mat.A << ", "<< mat.Z << ", "<< mat.rho );
336  m_accumulatedMaterialXX0 += stepl/mat.X0;
337  m_accumulatedRhoS += stepl*mat.rho;
338  ATH_MSG_VERBOSE("- accumulated material information X/x0 = " << m_accumulatedMaterialXX0);
339  ATH_MSG_VERBOSE("- accumulated effective densitity rho*s = " << m_accumulatedRhoS);
340  ATH_MSG_VERBOSE("- to layer with index "<< layer->layerIndex().value() << " from '"<< associatedVolume->volumeName() << "'.");
341  // record the plain information w/o correction & intersection
342  if (m_mapMaterial && !m_materialMapper.empty()){
343  Trk::AssociatedMaterial am(pos, stepl, mat, 1., associatedVolume, layer);
344  m_materialMapper->recordMaterialHit(am, positionOnLayer);
345  ATH_MSG_VERBOSE(" - associated material produced as : " << am);
346  }
347 
348  } else if (layer) {
349  ATH_MSG_WARNING("- associate hit - the layer with index " << layer->layerIndex().value() << " was not found - should not happen!");
350  }
351  // return
352  return true;
353 }
354 
355 
357  Trk::LayerMaterialMap* propSet)
358 {
359 
360  if (!propSet) return;
361 
362  ATH_MSG_INFO("Processing TrackingVolume: "<< tvol.volumeName() );
363 
364  // ----------------------------------- loop over confined layers ------------------------------------------
365  Trk::BinnedArray< Trk::Layer >* confinedLayers = tvol.confinedLayers();
366  if (confinedLayers) {
367  // get the objects in a vector-like format
369  ATH_MSG_INFO("--> found : "<< layers.size() << "confined Layers");
370  // the iterator over the vector
371  // loop over layers
372  for (Trk::Layer* layer : layers) {
373  // assign the material and output
374  if (layer && (*layer).layerIndex().value() ) {
375  ATH_MSG_INFO(" > LayerIndex: "<< (*layer).layerIndex() );
376  // set the material!
377  if (propSet) {
378  // find the pair
379  auto curIt = propSet->find((*layer).layerIndex());
380  if (curIt != propSet->end()) {
381  ATH_MSG_INFO("LayerMaterial assigned for Layer with index: "<< (*layer).layerIndex() );
382  // set it to the layer
383  layer->assignMaterialProperties(*((*curIt).second), 1.);
384  }
385  }
386  }
387  }
388  }
389 
390  // ----------------------------------- loop over confined volumes -----------------------------
392  if (confinedVolumes) {
393  // get the objects in a vector-like format
395  ATH_MSG_INFO("--> found : "<< volumes.size() << "confined TrackingVolumes");
396  // loop over volumes
397  for (const auto & volume : volumes) {
398  // assing the material and output
399  if (volume)
400  assignLayerMaterialProperties(*volume, propSet); // call itself recursively
401  }
402  }
403 }
404 
405 
407 {
408 
409  ATH_MSG_INFO("========================================================================================= ");
410  ATH_MSG_INFO("finalize() starts ...");
411 
412 #ifdef TRKDETDESCR_MEMUSAGE
413  m_memoryLogger.refresh(getpid());
414  ATH_MSG_INFO("[ memory usage ] Start building of material maps: " );
415  ATH_MSG_INFO( m_memoryLogger );
416 #endif
417 
418  // create a dedicated LayerMaterialMap by layerMaterialCreator;
419  std::map< std::string, Trk::LayerMaterialMap* > layerMaterialMaps;
420  for ( auto& lmcIter : m_layerMaterialCreators ){
421  ATH_MSG_INFO("-> Creating material map '"<< lmcIter->layerMaterialName() << "' from creator "<< lmcIter.typeAndName() );
422  layerMaterialMaps[lmcIter->layerMaterialName()] = new Trk::LayerMaterialMap();
423  }
424 
425  ATH_MSG_INFO( m_layerRecords.size() << " LayerMaterialRecords to be finalized for this material mapping run.");
426 
427  // loop over the layers and output the stuff --- fill the associatedLayerMaterialProperties
428  for ( auto& lIter : m_layerRecords ) {
429  // Get the key map, the layer & the volume name
430  const Trk::Layer* layer = lIter.first;
431  Trk::LayerIndex layerKey = layer->layerIndex();
432  // get the enclosing tracking volume
433  const Trk::TrackingVolume* eVolume = layer->enclosingTrackingVolume();
434  // assign the string
435  std::string vName = eVolume ? (eVolume->volumeName()) : " BoundaryCollection ";
436  ATH_MSG_INFO("Finalize MaterialAssociation for Layer "<< layerKey.value() << " in " << vName );
437  // normalize - use m_finalizeRunDebug
438  (lIter.second).finalizeRun(m_mapComposition);
439  // output the material to the analyser if registered
440  if (!m_layerMaterialRecordAnalyser.empty() && m_layerMaterialRecordAnalyser->analyseLayerMaterial(*layer, lIter.second).isFailure() )
441  ATH_MSG_WARNING("Could not analyse the LayerMaterialRecord for layer "<< layerKey.value() );
442  // check if we have analysers per creator
443  bool analyse = (m_layerMaterialCreators.size() == m_layerMaterialAnalysers.size());
444  // and now use the creators to make the maps out of the LayerMaterialRecord
445  size_t ilmc = 0;
446  for ( auto& lmcIter : m_layerMaterialCreators ){
447  // call the creator and register in the according map
448 #ifdef TRKDETDESCR_MEMUSAGE
449  m_memoryLogger.refresh(getpid());
450  ATH_MSG_INFO("[ memory usage ] Before building the map for Layer "<< layerKey.value() );
451  ATH_MSG_INFO( m_memoryLogger );
452 #endif
453  const Trk::LayerMaterialProperties* lMaterial = lmcIter->createLayerMaterial(lIter.second);
454 #ifdef TRKDETDESCR_MEMUSAGE
455  m_memoryLogger.refresh(getpid());
456  ATH_MSG_INFO("[ memory usage ] After building the map for Layer "<< layerKey.value() );
457  ATH_MSG_INFO( m_memoryLogger );
458 #endif
459  if (lMaterial)
460  ATH_MSG_VERBOSE("LayerMaterial map created as "<< *lMaterial );
461  // insert the created map for the given layer
462  (*layerMaterialMaps[lmcIter->layerMaterialName()])[layerKey.value()] = lMaterial;
463  // analyse the it if configured
464  if (analyse && lMaterial && (m_layerMaterialAnalysers[ilmc]->analyseLayerMaterial(*layer, *lMaterial)).isFailure() )
465  ATH_MSG_WARNING("Could not analyse created LayerMaterialProperties for layer "<< layerKey.value() );
466  ++ilmc;
467  }
468  }
469 
470  ATH_MSG_INFO("Finalize map synchronization and write the maps to the DetectorStore.");
471 
472  for (auto& ilmIter : layerMaterialMaps ){
473  // elementTable handling - if existent
474  if (m_mapComposition){
475  Trk::SharedObject<const Trk::ElementTable> tElementTable(new Trk::ElementTable(*m_elementTable));
476  ilmIter.second->updateElementTable(tElementTable);
477  if (ilmIter.second->elementTable()){
478  ATH_MSG_INFO("ElementTable for LayerMaterialMap '" << ilmIter.first << "' found and syncrhonized." );
479  ATH_MSG_INFO( *(ilmIter.second->elementTable()) );
480  }
481  }
482  // detector store writing
483  if ( (detStore()->record(ilmIter.second, ilmIter.first, false)).isFailure()){
484  ATH_MSG_ERROR( "Writing of LayerMaterialMap with name '" << ilmIter.first << "' was not successful." );
485  delete ilmIter.second;
486  } else ATH_MSG_INFO( "LayerMaterialMap: " << ilmIter.first << " written to the DetectorStore!" );
487  }
488  delete m_elementTable;
489 
490 #ifdef TRKDETDESCR_MEMUSAGE
491  m_memoryLogger.refresh(getpid());
492  ATH_MSG_INFO( "[ memory usage ] At the end of the material map creation.");
493  ATH_MSG_INFO( m_memoryLogger );
494 #endif
495 
496  ATH_MSG_INFO( "========================================================================================= " );
497  ATH_MSG_INFO( " -> Total mapped hits : " << m_mapped );
498  double unmapped = (m_unmapped+m_mapped) ? double(m_unmapped)/double(m_unmapped+m_mapped) : 0.;
499  ATH_MSG_INFO( " -> Total (rel.) unmapped hits : " << m_unmapped << " (" << unmapped << ")" );
500  ATH_MSG_INFO( " -> Skipped (outisde) : " << m_skippedOutside );
501  ATH_MSG_INFO( "========================================================================================= " );
502  ATH_MSG_INFO( "finalize() successful");
503  return StatusCode::SUCCESS;
504 }
505 
506 
508 {
509  // either get a string volume or the highest one
510  const Trk::TrackingVolume* trackingVolume = trackingGeometry().highestTrackingVolume();
511 
512  // prepare the mapping volume
513  m_mappingVolume = trackingGeometry().trackingVolume(m_mappingVolumeName);
514 
515  // register the confined layers from the TrackingVolume
516  registerVolume(*trackingVolume, 0);
517 
518  ATH_MSG_INFO("Add "<< m_layerRecords.size() << " confined volume layers to mapping setup.");
519  ATH_MSG_INFO("Add "<< trackingGeometry().numBoundaryLayers() << " boundary layers to mapping setup.");
520 
521  // register the layers from boundary surfaces
522  for (const auto bLayerIter : trackingGeometry().boundaryLayers())
523  insertLayerMaterialRecord(*(bLayerIter.first));
524 
525  ATH_MSG_INFO("Map for "<< m_layerRecords.size() << " layers booked & prepared for mapping procedure");
526 
527  return StatusCode::SUCCESS;
528 
529 }
530 
531 
533 {
534  int sublevel = lvl+1;
535 
536  for (int indent=0; indent<sublevel; ++indent)
537  std::cout << " ";
538  std::cout << "TrackingVolume name: "<< tvol.volumeName() << std::endl;
539 
540  // all those to be processed
541  std::vector<const Trk::Layer*> volumeLayers;
542 
543  // collect all material layers that have layerMaterial
544  const Trk::BinnedArray< Trk::Layer >* confinedLayers = tvol.confinedLayers();
545  if (confinedLayers) {
546  // this go ahead with the layers
548  for (int indent=0; indent<sublevel; ++indent)
549  std::cout << " ";
550  std::cout << "- found : "<< layers.size() << "confined Layers"<< std::endl;
551  // loop over and fill them
552  auto clIter = layers.begin();
553  auto clIterE = layers.end();
554  for ( ; clIter != clIterE; ++clIter ) {
555  // only take layers with MaterialProperties defined and which are within the mapping volume
556  const Amg::Vector3D& sReferencePoint = (*clIter)->surfaceRepresentation().globalReferencePoint();
557  bool insideMappingVolume = m_mappingVolume ? m_mappingVolume->inside(sReferencePoint) : true;
558  if ((*clIter)->layerMaterialProperties() && insideMappingVolume)
559  volumeLayers.push_back((*clIter));
560  }
561  }
562 
563  // now create LayerMaterialRecords for all
564  for ( auto& lIter : volumeLayers )
565  insertLayerMaterialRecord(*lIter);
566 
567  // step dopwn the navigation tree to reach the confined volumes
568  const Trk::BinnedArray<Trk::TrackingVolume >* confinedVolumes = tvol.confinedVolumes();
569  if (confinedVolumes) {
571 
572  for (int indent=0; indent<sublevel; ++indent)
573  std::cout << " ";
574  std::cout << "- found : "<< volumes.size() << "confined TrackingVolumes"<< std::endl;
575  // loop over the confined volumes
576  auto volumesIter = volumes.begin();
577  for (; volumesIter != volumes.end(); ++volumesIter)
578  if (*volumesIter) {
579  registerVolume(**volumesIter, sublevel);
580  }
581  }
582 
583 }
584 
586  // first occurrance, create a new LayerMaterialRecord
587  // - get the bin utility for the binned material (if necessary)
588  // - get the material first
589  const Trk::LayerMaterialProperties* layerMaterialProperties = lay.layerMaterialProperties();
590  // - dynamic cast to the BinnedLayerMaterial
591  const Trk::BinnedLayerMaterial* layerBinnedMaterial
592  = dynamic_cast<const Trk::BinnedLayerMaterial*>(layerMaterialProperties);
593  // get the binned array
594  const Trk::BinUtility* layerMaterialBinUtility = (layerBinnedMaterial) ? layerBinnedMaterial->binUtility() : nullptr;
595  // now fill the layer material record
596  if (layerMaterialBinUtility){
597  // create a new Layer Material record in the map
598  Trk::LayerMaterialRecord lmr((m_useLayerThickness ? lay.thickness() : 1.),
599  layerMaterialBinUtility,
600  (Trk::MaterialAssociationType)m_associationType);
601  // and fill it into the map
602  m_layerRecords[&lay] = lmr;
603  }
604 }
605 
607  std::stringstream msg;
608  msg << "Failed to get conditions data " << m_trackingGeometryReadKey.key() << ".";
609  throw std::runtime_error(msg.str());
610 }
611 
612 
python.SystemOfUnits.second
int second
Definition: SystemOfUnits.py:120
Trk::ExtrapolationMode::CollectPassive
@ CollectPassive
Definition: ExtrapolationCell.h:57
MaterialMapping.h
Trk::ExtrapolationCell::extrapolationSteps
std::vector< ExtrapolationStep< T > > extrapolationSteps
parameters on sensitive detector elements
Definition: ExtrapolationCell.h:292
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
Trk::ExtrapolationCell::addConfigurationMode
void addConfigurationMode(ExtrapolationMode::eMode em)
add a configuration mode
Definition: ExtrapolationCell.h:344
MaterialProperties.h
Trk::ExtrapolationMode::CollectBoundary
@ CollectBoundary
Definition: ExtrapolationCell.h:58
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:79
SG::ReadHandle
Definition: StoreGate/StoreGate/ReadHandle.h:70
Trk::MaterialMapping::m_inputEventElementTable
SG::ReadHandleKey< Trk::ElementTable > m_inputEventElementTable
input event table
Definition: MaterialMapping.h:167
AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T > &t)
Definition: AthCommonDataStore.h:145
mat
GeoMaterial * mat
Definition: LArDetectorConstructionTBEC.cxx:53
Trk::MaterialMapping::registerVolume
void registerVolume(const Trk::TrackingVolume &tvol, int lvl)
Output information with Level.
Definition: MaterialMapping.cxx:532
BinUtility.h
module_driven_slicing.layers
layers
Definition: module_driven_slicing.py:114
Trk::MaterialMapping::m_etaCutOff
double m_etaCutOff
general steering
Definition: MaterialMapping.h:131
HomogeneousLayerMaterial.h
Trk::MaterialMapping::handleTrackingGeometry
StatusCode handleTrackingGeometry()
Retrieve the TrackingGeometry and its informations.
Definition: MaterialMapping.cxx:507
NeutralParameters.h
read_hist_ntuple.t
t
Definition: read_hist_ntuple.py:5
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
python.compareNtuple.vName
vName
Definition: compareNtuple.py:23
Trk::MaterialMapping::throwFailedToGetTrackingGeometry
void throwFailedToGetTrackingGeometry() const
Definition: MaterialMapping.cxx:606
Trk::LayerMaterialMap
Definition: LayerMaterialMap.h:32
Trk::MaterialMapping::m_minCompositionFraction
double m_minCompositionFraction
minimal fraction to be accounted for the composition recording
Definition: MaterialMapping.h:162
Trk::TrackingVolume::confinedLayers
const LayerArray * confinedLayers() const
Return the subLayer array.
Trk::MaterialStep
Definition: MaterialStep.h:34
GeometryStatics.h
Trk::MaterialMapping::m_associationType
int m_associationType
Definition: MaterialMapping.h:134
Amg::toString
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
Definition: GeoPrimitivesToStringConverter.h:40
Trk::MaterialMapping::assignLayerMaterialProperties
void assignLayerMaterialProperties(Trk::TrackingVolume &tvol, Trk::LayerMaterialMap *propSet)
create the LayerMaterialRecord *‍/
Definition: MaterialMapping.cxx:356
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
geometry_dat_to_json.indent
indent
Definition: geometry_dat_to_json.py:18
Trk::MaterialMapping::m_etaSide
int m_etaSide
needed for debugging: -1 negative | 0 all | 1 positive
Definition: MaterialMapping.h:132
Trk::LayerIndex
Definition: LayerIndex.h:37
LayerMaterialRecord.h
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
Trk::MaterialMapping::m_mapComposition
bool m_mapComposition
map the composition of the material
Definition: MaterialMapping.h:161
TRT::Hit::layer
@ layer
Definition: HitInfo.h:79
Trk::ExtrapolationCode::isSuccess
bool isSuccess() const
return success
Definition: ExtrapolationCell.h:153
Trk::LayerMaterialProperties
Definition: LayerMaterialProperties.h:62
Trk::ExtrapolationCode
Definition: ExtrapolationCell.h:105
Trk::ExtrapolationCell::endParameters
T * endParameters
by pointer - are newly created and can be optionally 0
Definition: ExtrapolationCell.h:238
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
BinnedLayerMaterial.h
Trk::Surface::materialLayer
const Trk::Layer * materialLayer() const
return the material Layer
Trk::ParametersBase
Definition: ParametersBase.h:55
Trk::CurvilinearParametersT
Definition: CurvilinearParametersT.h:48
Trk::BinnedLayerMaterial
Definition: BinnedLayerMaterial.h:33
xAOD::double
double
Definition: CompositeParticle_v1.cxx:159
Trk::MaterialMapping::m_mappingVolumeName
std::string m_mappingVolumeName
Definition: MaterialMapping.h:123
Trk::LayerIndex::value
int value() const
layerIndex expressed in an integer
Definition: LayerIndex.h:71
LayerIndex.h
AthAlgorithm
Definition: AthAlgorithm.h:47
SG::ReadHandle::isValid
virtual bool isValid() override final
Can the handle be successfully dereferenced?
Trk::MaterialMapping::~MaterialMapping
~MaterialMapping()
Default Destructor.
Trk::AssociatedMaterial
Definition: AssociatedMaterial.h:33
Trk::BinUtility
Definition: BinUtility.h:39
Trk::MaterialMapping::initialize
StatusCode initialize()
standard Athena-Algorithm method
Definition: MaterialMapping.cxx:88
Trk::LayerMaterialRecord
Definition: LayerMaterialRecord.h:42
python.PyKernel.detStore
detStore
Definition: PyKernel.py:41
Trk::TrackingVolume::volumeName
const std::string & volumeName() const
Returns the VolumeName - for debug reason, might be depreciated later.
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
Trk::BinnedArray::arrayObjects
virtual BinnedArraySpan< T *const > arrayObjects()=0
Return all objects of the Array non-const we can still modify the T.
Trk::MaterialMapping::execute
StatusCode execute()
standard Athena-Algorithm method
Definition: MaterialMapping.cxx:116
Trk::Layer::thickness
double thickness() const
Return the Thickness of the Layer.
LayerMaterialMap.h
Trk::MaterialAssociationType
MaterialAssociationType
Definition: MaterialAssociationType.h:13
Amg::Vector3D
Eigen::Matrix< double, 3, 1 > Vector3D
Definition: GeoPrimitives.h:47
python.LumiBlobConversion.pos
pos
Definition: LumiBlobConversion.py:18
Trk::MaterialMapping::m_useLayerThickness
bool m_useLayerThickness
use the actual layer thickness
Definition: MaterialMapping.h:133
Trk::MaterialMapping::m_mapMaterial
bool m_mapMaterial
Mapper and Inspector.
Definition: MaterialMapping.h:156
Trk::Surface::associatedLayer
const Trk::Layer * associatedLayer() const
return the associated Layer
Trk::MaterialMapping::associateHit
bool associateHit(const Trk::Layer &tvol, const Amg::Vector3D &pos, const Amg::Vector3D &layerHitPosition, double stepl, const Trk::Material &mat)
Associate the Step to the Layer.
Definition: MaterialMapping.cxx:311
Trk::BinnedLayerMaterial::binUtility
virtual const BinUtility * binUtility() const override
Return the BinUtility.
Definition: BinnedLayerMaterial.h:121
TrackingVolume.h
Trk::ExtrapolationCell
Definition: ExtrapolationCell.h:231
Trk::MaterialMapping::m_checkForEmptyHits
bool m_checkForEmptyHits
use extrapoaltion engine to check for empty hits
Definition: MaterialMapping.h:115
Trk::SharedObject
std::shared_ptr< T > SharedObject
Definition: SharedObject.h:24
Trk::ElementTable
Definition: ElementTable.h:24
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
Trk::TrackingVolume::confinedVolumes
const TrackingVolumeArray * confinedVolumes() const
Return the subLayer array.
Trk::MaterialMapping::MaterialMapping
MaterialMapping(const std::string &name, ISvcLocator *pSvcLocator)
Standard Athena-Algorithm Constructor.
Definition: MaterialMapping.cxx:36
Trk::Layer::layerMaterialProperties
const LayerMaterialProperties * layerMaterialProperties() const
getting the LayerMaterialProperties including full/pre/post update
Trk::ExtrapolationCell::navigationCurvilinear
bool navigationCurvilinear
stay in curvilinear parameters where possible, default : true
Definition: ExtrapolationCell.h:284
GeoPrimitivesToStringConverter.h
LArCellBinning.step
step
Definition: LArCellBinning.py:158
Trk::MaterialMapping::insertLayerMaterialRecord
void insertLayerMaterialRecord(const Trk::Layer &lay)
Definition: MaterialMapping.cxx:585
AssociatedMaterial.h
physics_parameters.parameters
parameters
Definition: physics_parameters.py:144
Trk::Material
Definition: Material.h:116
Trk::MaterialMapping::m_inputMaterialStepCollection
SG::ReadHandleKey< MaterialStepCollection > m_inputMaterialStepCollection
output / input steering
Definition: MaterialMapping.h:127
Trk::MaterialMapping::finalize
StatusCode finalize()
standard Athena-Algorithm method
Definition: MaterialMapping.cxx:406
Trk::BinnedArray
Definition: BinnedArray.h:38
Trk::Surface
Definition: Tracking/TrkDetDescr/TrkSurfaces/TrkSurfaces/Surface.h:75
Trk::ExtrapolationMode::StopAtBoundary
@ StopAtBoundary
Definition: ExtrapolationCell.h:55
Trk::TrackingVolume
Definition: TrackingVolume.h:121
Trk::BinnedArraySpan
std::span< T > BinnedArraySpan
Definition: BinnedArray.h:34
python.AutoConfigFlags.msg
msg
Definition: AutoConfigFlags.py:7
mag
Scalar mag() const
mag method
Definition: AmgMatrixBasePlugin.h:25
LayerMaterialProperties.h
Trk::MaterialMapping::m_layerMaterialScreenOutput
int m_layerMaterialScreenOutput
Definition: MaterialMapping.h:183
HitType::unmapped
@ unmapped
Trk::Layer
Definition: Layer.h:73
Trk::Layer::layerIndex
const LayerIndex & layerIndex() const
get the layerIndex
MaterialStep.h
CompressedLayerMaterial.h