CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
MaterialEffects Class Reference

#include <MaterialEffects.h>

Public Member Functions

double energyLoss () const
 Return the energy loss by ionization in the current layer. More...
 
EnergyLossSimulatorenergyLossSimulator () const
 Return the Energy Loss engine. More...
 
void interact (FSimEvent &simEvent, const TrackerLayer &layer, ParticlePropagator &PP, unsigned i, RandomEngineAndDistribution const *)
 
 MaterialEffects (const edm::ParameterSet &matEff)
 Constructor. More...
 
MultipleScatteringSimulatormultipleScatteringSimulator () const
 Return the Multiple Scattering engine. More...
 
MuonBremsstrahlungSimulatormuonBremsstrahlungSimulator () const
 Return the Muon Bremsstrahlung engine. More...
 
void save ()
 Save nuclear interaction information. More...
 
double thickness () const
 Return the thickness of the current layer. More...
 
 ~MaterialEffects ()
 Default destructor. More...
 

Private Member Functions

GlobalVector normalVector (const TrackerLayer &layer, ParticlePropagator &myTrack) const
 The vector normal to the surface traversed. More...
 
double radLengths (const TrackerLayer &layer, ParticlePropagator &myTrack)
 The number of radiation lengths traversed. More...
 

Private Attributes

BremsstrahlungSimulatorBremsstrahlung
 
EnergyLossSimulatorEnergyLoss
 
MultipleScatteringSimulatorMultipleScattering
 
MuonBremsstrahlungSimulatorMuonBremsstrahlung
 
MaterialEffectsSimulatorNuclearInteraction
 
PairProductionSimulatorPairProduction
 
double pTmin
 
double theEnergyLoss
 
GlobalVector theNormalVector
 
double theTECFudgeFactor
 
double theThickness
 
bool use_hardcoded
 

Detailed Description

Definition at line 51 of file MaterialEffects.h.

Constructor & Destructor Documentation

MaterialEffects::MaterialEffects ( const edm::ParameterSet matEff)

Constructor.

Definition at line 26 of file MaterialEffects.cc.

References patCaloMETCorrections_cff::A, Bremsstrahlung, candidateVertexArbitrator_cfi::distCut, EnergyLoss, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), GeV, mps_fire::i, analyzePatCleaning_cfg::inputFile, MultipleScattering, MuonBremsstrahlung, NuclearInteraction, PairProduction, pTmin, plotFactory::ratios, AlCaHLTBitMon_QueryRunRegistry::string, theTECFudgeFactor, use_hardcoded, and DOFs::Z.

27  : PairProduction(nullptr), Bremsstrahlung(nullptr),MuonBremsstrahlung(nullptr),
28  MultipleScattering(nullptr), EnergyLoss(nullptr),
29  NuclearInteraction(nullptr),
30  pTmin(999.), use_hardcoded(true)
31 {
32  // Set the minimal photon energy for a Brem from e+/-
33 
34  use_hardcoded = matEff.getParameter<bool >("use_hardcoded_geometry");
35 
36  bool doPairProduction = matEff.getParameter<bool>("PairProduction");
37  bool doBremsstrahlung = matEff.getParameter<bool>("Bremsstrahlung");
38  bool doEnergyLoss = matEff.getParameter<bool>("EnergyLoss");
39  bool doMultipleScattering = matEff.getParameter<bool>("MultipleScattering");
40  bool doNuclearInteraction = matEff.getParameter<bool>("NuclearInteraction");
41  bool doG4NuclInteraction = matEff.getParameter<bool>("G4NuclearInteraction");
42  bool doMuonBremsstrahlung = matEff.getParameter<bool>("MuonBremsstrahlung");
43 
44  double A = matEff.getParameter<double>("A");
45  double Z = matEff.getParameter<double>("Z");
46  double density = matEff.getParameter<double>("Density");
47  double radLen = matEff.getParameter<double>("RadiationLength");
48 
49  // Set the minimal pT before giving up the dE/dx treatment
50 
51  if ( doPairProduction ) {
52 
53  double photonEnergy = matEff.getParameter<double>("photonEnergy");
54  PairProduction = new PairProductionSimulator(photonEnergy);
55  }
56 
57  if ( doBremsstrahlung ) {
58 
59  double bremEnergy = matEff.getParameter<double>("bremEnergy");
60  double bremEnergyFraction = matEff.getParameter<double>("bremEnergyFraction");
61  Bremsstrahlung = new BremsstrahlungSimulator(bremEnergy,
62  bremEnergyFraction);
63  }
64  //muon Brem+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
65  if ( doMuonBremsstrahlung ) {
66 
67  double bremEnergy = matEff.getParameter<double>("bremEnergy");
68  double bremEnergyFraction = matEff.getParameter<double>("bremEnergyFraction");
69  MuonBremsstrahlung = new MuonBremsstrahlungSimulator(A,Z,density,radLen,bremEnergy,
70  bremEnergyFraction);
71 
72  }
73 
74  //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
75 
76  if ( doEnergyLoss ) {
77 
78  pTmin = matEff.getParameter<double>("pTmin");
79  EnergyLoss = new EnergyLossSimulator(A,Z,density,radLen);
80 
81  }
82 
83  if ( doMultipleScattering ) {
84 
85  MultipleScattering = new MultipleScatteringSimulator(A,Z,density,radLen);
86 
87  }
88 
89 
90  if ( doNuclearInteraction ) {
91 
92  // The energies simulated
93  std::vector<double> hadronEnergies
94  = matEff.getUntrackedParameter<std::vector<double> >("hadronEnergies");
95 
96  // The particle types simulated
97  std::vector<int> hadronTypes
98  = matEff.getUntrackedParameter<std::vector<int> >("hadronTypes");
99 
100  // The corresponding particle names
101  std::vector<std::string> hadronNames
102  = matEff.getUntrackedParameter<std::vector<std::string> >("hadronNames");
103 
104  // The corresponding particle masses
105  std::vector<double> hadronMasses
106  = matEff.getUntrackedParameter<std::vector<double> >("hadronMasses");
107 
108  // The smallest momentum for inelastic interactions
109  std::vector<double> hadronPMin
110  = matEff.getUntrackedParameter<std::vector<double> >("hadronMinP");
111 
112  // The interaction length / radiation length ratio for each particle type
113  std::vector<double> lengthRatio
114  = matEff.getParameter<std::vector<double> >("lengthRatio");
115  // std::map<int,double> lengthRatio;
116  // for ( unsigned i=0; i<theLengthRatio.size(); ++i )
117  // lengthRatio[ hadronTypes[i] ] = theLengthRatio[i];
118 
119  // A global fudge factor for TEC layers (which apparently do not react to
120  // hadrons the same way as all other layers...
121  theTECFudgeFactor = matEff.getParameter<double>("fudgeFactor");
122 
123  // The evolution of the interaction lengths with energy
124  std::vector<double> theRatios
125  = matEff.getUntrackedParameter<std::vector<double> >("ratios");
126  //std::map<int,std::vector<double> > ratios;
127  //for ( unsigned i=0; i<hadronTypes.size(); ++i ) {
128  // for ( unsigned j=0; j<hadronEnergies.size(); ++j ) {
129  // ratios[ hadronTypes[i] ].push_back(theRatios[ i*hadronEnergies.size() + j ]);
130  // }
131  //}
132  std::vector< std::vector<double> > ratios;
133  ratios.resize(hadronTypes.size());
134  for ( unsigned i=0; i<hadronTypes.size(); ++i ) {
135  for ( unsigned j=0; j<hadronEnergies.size(); ++j ) {
136  ratios[i].push_back(theRatios[ i*hadronEnergies.size() + j ]);
137  }
138  }
139 
140  // The smallest momentum for elastic interactions
141  double pionEnergy
142  = matEff.getParameter<double>("pionEnergy");
143 
144  // The algorithm to compute the distance between primary and secondaries
145  // when a nuclear interaction occurs
146  unsigned distAlgo
147  = matEff.getParameter<unsigned>("distAlgo");
148  double distCut
149  = matEff.getParameter<double>("distCut");
150 
151  // The file to read the starting interaction in each files
152  // (random reproducibility in case of a crash)
154  = matEff.getUntrackedParameter<std::string>("inputFile");
155 
156  // Build the ID map (i.e., what is to be considered as a proton, etc...)
157  std::map<int,int> idMap;
158  // Protons
159  std::vector<int> idProtons
160  = matEff.getUntrackedParameter<std::vector<int> >("protons");
161  for ( unsigned i=0; i<idProtons.size(); ++i )
162  idMap[idProtons[i]] = 2212;
163  // Anti-Protons
164  std::vector<int> idAntiProtons
165  = matEff.getUntrackedParameter<std::vector<int> >("antiprotons");
166  for ( unsigned i=0; i<idAntiProtons.size(); ++i )
167  idMap[idAntiProtons[i]] = -2212;
168  // Neutrons
169  std::vector<int> idNeutrons
170  = matEff.getUntrackedParameter<std::vector<int> >("neutrons");
171  for ( unsigned i=0; i<idNeutrons.size(); ++i )
172  idMap[idNeutrons[i]] = 2112;
173  // Anti-Neutrons
174  std::vector<int> idAntiNeutrons
175  = matEff.getUntrackedParameter<std::vector<int> >("antineutrons");
176  for ( unsigned i=0; i<idAntiNeutrons.size(); ++i )
177  idMap[idAntiNeutrons[i]] = -2112;
178  // K0L's
179  std::vector<int> idK0Ls
180  = matEff.getUntrackedParameter<std::vector<int> >("K0Ls");
181  for ( unsigned i=0; i<idK0Ls.size(); ++i )
182  idMap[idK0Ls[i]] = 130;
183  // K+'s
184  std::vector<int> idKplusses
185  = matEff.getUntrackedParameter<std::vector<int> >("Kplusses");
186  for ( unsigned i=0; i<idKplusses.size(); ++i )
187  idMap[idKplusses[i]] = 321;
188  // K-'s
189  std::vector<int> idKminusses
190  = matEff.getUntrackedParameter<std::vector<int> >("Kminusses");
191  for ( unsigned i=0; i<idKminusses.size(); ++i )
192  idMap[idKminusses[i]] = -321;
193  // pi+'s
194  std::vector<int> idPiplusses
195  = matEff.getUntrackedParameter<std::vector<int> >("Piplusses");
196  for ( unsigned i=0; i<idPiplusses.size(); ++i )
197  idMap[idPiplusses[i]] = 211;
198  // pi-'s
199  std::vector<int> idPiminusses
200  = matEff.getUntrackedParameter<std::vector<int> >("Piminusses");
201  for ( unsigned i=0; i<idPiminusses.size(); ++i )
202  idMap[idPiminusses[i]] = -211;
203 
204  // Construction
205  if ( doG4NuclInteraction ) {
206  double elimit = matEff.getParameter<double>("EkinBertiniGeV")*CLHEP::GeV;
207  double eth = matEff.getParameter<double>("EkinLimitGeV")*CLHEP::GeV;
209  new NuclearInteractionFTFSimulator(distAlgo, distCut, elimit, eth);
210  } else {
212  new NuclearInteractionSimulator(hadronEnergies, hadronTypes, hadronNames,
213  hadronMasses, hadronPMin, pionEnergy,
214  lengthRatio, ratios, idMap,
215  inputFile, distAlgo, distCut);
216  }
217  }
218 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const double GeV
Definition: MathUtil.h:16
PairProductionSimulator * PairProduction
MaterialEffectsSimulator * NuclearInteraction
MuonBremsstrahlungSimulator * MuonBremsstrahlung
BremsstrahlungSimulator * Bremsstrahlung
EnergyLossSimulator * EnergyLoss
MultipleScatteringSimulator * MultipleScattering
MaterialEffects::~MaterialEffects ( )

Default destructor.

Definition at line 220 of file MaterialEffects.cc.

References Bremsstrahlung, EnergyLoss, MultipleScattering, MuonBremsstrahlung, NuclearInteraction, and PairProduction.

220  {
221 
222  if ( PairProduction ) delete PairProduction;
223  if ( Bremsstrahlung ) delete Bremsstrahlung;
224  if ( EnergyLoss ) delete EnergyLoss;
227  //Muon Brem
229 }
PairProductionSimulator * PairProduction
MaterialEffectsSimulator * NuclearInteraction
MuonBremsstrahlungSimulator * MuonBremsstrahlung
BremsstrahlungSimulator * Bremsstrahlung
EnergyLossSimulator * EnergyLoss
MultipleScatteringSimulator * MultipleScattering

Member Function Documentation

double MaterialEffects::energyLoss ( ) const
inline

Return the energy loss by ionization in the current layer.

Definition at line 77 of file MaterialEffects.h.

Referenced by TrajectoryManager::createPSimHits().

77 { return theEnergyLoss; }
EnergyLossSimulator* MaterialEffects::energyLossSimulator ( ) const
inline

Return the Energy Loss engine.

Definition at line 85 of file MaterialEffects.h.

Referenced by MuonSimHitProducer::applyMaterialEffects(), and CalorimetryManager::MuonMipSimulation().

85  {
86  return EnergyLoss;
87  }
EnergyLossSimulator * EnergyLoss
void MaterialEffects::interact ( FSimEvent simEvent,
const TrackerLayer layer,
ParticlePropagator PP,
unsigned  i,
RandomEngineAndDistribution const *  random 
)

Steer the various interaction processes in the Tracker Material and update the FSimEvent

Energy loss

Multiple scattering

Definition at line 231 of file MaterialEffects.cc.

References funct::abs(), FBaseSimEvent::addSimTrack(), FBaseSimEvent::addSimVertex(), MaterialEffectsSimulator::beginDaughters(), FSimVertexType::BREM_VERTEX, Bremsstrahlung, RawParticle::charge(), MaterialEffectsSimulator::closestDaughterId(), MaterialEffectsSimulator::endDaughters(), EnergyLoss, TrackerLayer::layerNumber(), MultipleScattering, MuonBremsstrahlung, MaterialEffectsSimulator::nDaughters(), normalVector(), FSimVertexType::NUCL_VERTEX, NuclearInteraction, FSimVertexType::PAIR_VERTEX, PairProduction, RawParticle::pid(), FSimVertex::position(), pTmin, radLengths(), FSimTrack::setClosestDaughterId(), MaterialEffectsSimulator::setNormalVector(), theEnergyLoss, theNormalVector, theTECFudgeFactor, FBaseSimEvent::track(), MaterialEffectsSimulator::updateState(), use_hardcoded, FSimTrack::vertex(), and RawParticle::vertex().

Referenced by TrajectoryManager::reconstruct().

235  {
236 
238  double radlen;
239  theEnergyLoss = 0;
240  theNormalVector = normalVector(layer,myTrack);
241  radlen = radLengths(layer,myTrack);
242 
243  //std::cout << "### MaterialEffects: for Track= " << itrack << " in layer #"
244  // << layer.layerNumber() << std::endl;
245  //std::cout << myTrack << std::endl;
246 
247 //-------------------
248 // Photon Conversion
249 //-------------------
250 
251  if ( PairProduction && myTrack.pid()==22 ) {
252 
253  //
254  PairProduction->updateState(myTrack, radlen, random);
255 
256  if ( PairProduction->nDaughters() ) {
257  //add a vertex to the mother particle
258  int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
260 
261  // Check if it is a valid vertex first:
262  if (ivertex>=0) {
263  // This was a photon that converted
264  for ( DaughterIter = PairProduction->beginDaughters();
265  DaughterIter != PairProduction->endDaughters();
266  ++DaughterIter) {
267 
268  mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
269 
270  }
271  // The photon converted. Return.
272  return;
273  }
274  else {
275  edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in photon conv. -> " << ivertex << std::endl;
276  }
277 
278  }
279 
280  }
281 
282  if ( myTrack.pid() == 22 ) return;
283 
284 //------------------------
285 // Nuclear interactions
286 //------------------------
287 
288  if ( NuclearInteraction && abs(myTrack.pid()) > 100
289  && abs(myTrack.pid()) < 1000000) {
290 
291  // Simulate a nuclear interaction
292  double factor = 1.0;
293  if(use_hardcoded){
294  if (layer.layerNumber() >= 19 && layer.layerNumber() <= 27 )
295  factor = theTECFudgeFactor;
296  }
297  NuclearInteraction->updateState(myTrack, radlen*factor, random);
298 
299  //std::cout << "MaterialEffects: nDaughters= "
300  // << NuclearInteraction->nDaughters() << std::endl;
301  if ( NuclearInteraction->nDaughters() ) {
302 
303  //add a end vertex to the mother particle
304  int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
306  //std::cout << "ivertex= " << ivertex << " nDaughters= "
307  // << NuclearInteraction->nDaughters() << std::endl;
308  // Check if it is a valid vertex first:
309  if (ivertex>=0) {
310  // This was a hadron that interacted inelastically
311  int idaugh = 0;
312  for ( DaughterIter = NuclearInteraction->beginDaughters();
313  DaughterIter != NuclearInteraction->endDaughters();
314  ++DaughterIter) {
315 
316  // The daughter in the event
317  int daughId = mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
318 
319  // Store the closest daughter in the mother info (for later tracking purposes)
320  if ( NuclearInteraction->closestDaughterId() == idaugh ) {
321  if ( mySimEvent.track(itrack).vertex().position().Pt() < 4.0 )
322  mySimEvent.track(itrack).setClosestDaughterId(daughId);
323  }
324  ++idaugh;
325  }
326  // The hadron is destroyed. Return.
327  return;
328  }
329  else {
330  edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in nucl. int. -> " << ivertex << std::endl;
331  }
332 
333  }
334 
335  }
336 
337  if ( myTrack.charge() == 0 ) return;
338 
340 
341 //----------------
342 // Bremsstrahlung
343 //----------------
344 
345  if ( Bremsstrahlung && abs(myTrack.pid())==11 ) {
346 
347  Bremsstrahlung->updateState(myTrack,radlen, random);
348 
349  if ( Bremsstrahlung->nDaughters() ) {
350 
351  // Add a vertex, but do not attach it to the electron, because it
352  // continues its way...
353  int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
355 
356  // Check if it is a valid vertex first:
357  if (ivertex>=0) {
358  for ( DaughterIter = Bremsstrahlung->beginDaughters();
359  DaughterIter != Bremsstrahlung->endDaughters();
360  ++DaughterIter) {
361  mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
362  }
363  }
364  else {
365  edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in brem -> " << ivertex << std::endl;
366  }
367 
368  }
369 
370  }
371 
372 //---------------------------
373 // Muon_Bremsstrahlung
374 //--------------------------
375 
376  if ( MuonBremsstrahlung && abs(myTrack.pid())==13 ) {
377 
378  MuonBremsstrahlung->updateState(myTrack, radlen, random);
379 
380  if ( MuonBremsstrahlung->nDaughters() ) {
381 
382  // Add a vertex, but do not attach it to the muon, because it
383  // continues its way...
384  int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
386 
387  // Check if it is a valid vertex first:
388  if (ivertex>=0) {
389  for ( DaughterIter = MuonBremsstrahlung->beginDaughters();
390  DaughterIter != MuonBremsstrahlung->endDaughters();
391  ++DaughterIter) {
392  mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
393  }
394  }
395  else {
396  edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in muon brem -> " << ivertex << std::endl;
397  }
398 
399  }
400 
401  }
402 
406 
407  if ( EnergyLoss )
408  {
409  theEnergyLoss = myTrack.E();
410  EnergyLoss->updateState(myTrack, radlen, random);
411  theEnergyLoss -= myTrack.E();
412  }
413 
414 
418 
419  if ( MultipleScattering && myTrack.Pt() > pTmin ) {
420  // MultipleScattering->setNormalVector(normalVector(layer,myTrack));
422  MultipleScattering->updateState(myTrack,radlen, random);
423  }
424 
425 }
RHEP_const_iter beginDaughters() const
Returns const iterator to the beginning of the daughters list.
unsigned nDaughters() const
Returns the number of daughters.
unsigned int layerNumber() const
Returns the layer number.
Definition: TrackerLayer.h:82
int closestDaughterId()
The id of the closest charged daughter (filled for nuclear interactions only)
void setNormalVector(const GlobalVector &normal)
Sets the vector normal to the surface traversed.
PairProductionSimulator * PairProduction
GlobalVector normalVector(const TrackerLayer &layer, ParticlePropagator &myTrack) const
The vector normal to the surface traversed.
void updateState(ParticlePropagator &myTrack, double radlen, RandomEngineAndDistribution const *)
Compute the material effect (calls the sub class)
TRandom random
Definition: MVATrainer.cc:138
MaterialEffectsSimulator * NuclearInteraction
MuonBremsstrahlungSimulator * MuonBremsstrahlung
std::vector< RawParticle >::const_iterator RHEP_const_iter
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
RHEP_const_iter endDaughters() const
Returns const iterator to the end of the daughters list.
double radLengths(const TrackerLayer &layer, ParticlePropagator &myTrack)
The number of radiation lengths traversed.
BremsstrahlungSimulator * Bremsstrahlung
GlobalVector theNormalVector
EnergyLossSimulator * EnergyLoss
MultipleScatteringSimulator * MultipleScattering
MultipleScatteringSimulator* MaterialEffects::multipleScatteringSimulator ( ) const
inline

Return the Multiple Scattering engine.

Definition at line 80 of file MaterialEffects.h.

Referenced by MuonSimHitProducer::applyMaterialEffects().

80  {
81  return MultipleScattering;
82  }
MultipleScatteringSimulator * MultipleScattering
MuonBremsstrahlungSimulator* MaterialEffects::muonBremsstrahlungSimulator ( ) const
inline

Return the Muon Bremsstrahlung engine.

Definition at line 90 of file MaterialEffects.h.

Referenced by MuonSimHitProducer::applyMaterialEffects().

90  {
91  return MuonBremsstrahlung;
92  }
MuonBremsstrahlungSimulator * MuonBremsstrahlung
GlobalVector MaterialEffects::normalVector ( const TrackerLayer layer,
ParticlePropagator myTrack 
) const
private

The vector normal to the surface traversed.

Definition at line 469 of file MaterialEffects.cc.

References TrackerLayer::disk(), TrackerLayer::forward(), RawParticle::R(), RawParticle::X(), and RawParticle::Y().

Referenced by interact().

470  {
471  return layer.forward() ?
472  layer.disk()->normalVector() :
473  GlobalVector(myTrack.X(),myTrack.Y(),0.)/myTrack.R();
474 }
bool forward() const
Is the layer forward ?
Definition: TrackerLayer.h:70
double R() const
vertex radius
Definition: RawParticle.h:278
double Y() const
y of vertex
Definition: RawParticle.h:275
BoundDisk const * disk() const
Returns the surface.
Definition: TrackerLayer.h:79
double X() const
x of vertex
Definition: RawParticle.h:274
Global3DVector GlobalVector
Definition: GlobalVector.h:10
double MaterialEffects::radLengths ( const TrackerLayer layer,
ParticlePropagator myTrack 
)
private

The number of radiation lengths traversed.

Definition at line 428 of file MaterialEffects.cc.

References TrackerLayer::forward(), TrackerLayer::fudgeFactor(), TrackerLayer::fudgeMax(), TrackerLayer::fudgeMin(), TrackerLayer::fudgeNumber(), Surface::mediumProperties(), RawParticle::R(), MediumProperties::radLen(), TrackerLayer::surface(), theNormalVector, theThickness, and RawParticle::Z().

Referenced by interact().

429  {
430 
431  // Thickness of layer
433 
434  GlobalVector P(myTrack.Px(),myTrack.Py(),myTrack.Pz());
435 
436  // Effective length of track inside layer (considering crossing angle)
437  // double radlen = theThickness / fabs(P.dot(theNormalVector)/(P.mag()*theNormalVector.mag()));
438  double radlen = theThickness / fabs(P.dot(theNormalVector)) * P.mag();
439 
440  // This is a series of fudge factors (from the geometry description),
441  // to describe the layer inhomogeneities (services, cables, supports...)
442  double rad = myTrack.R();
443  double zed = fabs(myTrack.Z());
444 
445  double factor = 1;
446 
447  // Are there fudge factors for this layer
448  if ( layer.fudgeNumber() )
449 
450  // If yes, loop on them
451  for ( unsigned int iLayer=0; iLayer < layer.fudgeNumber(); ++iLayer ) {
452 
453  // Apply to R if forward layer, to Z if barrel layer
454  if ( ( layer.forward() && layer.fudgeMin(iLayer) < rad && rad < layer.fudgeMax(iLayer) ) ||
455  ( !layer.forward() && layer.fudgeMin(iLayer) < zed && zed < layer.fudgeMax(iLayer) ) ) {
456  factor = layer.fudgeFactor(iLayer);
457  break;
458  }
459 
460  }
461 
462  theThickness *= factor;
463 
464  return radlen * factor;
465 
466 }
float radLen() const
bool forward() const
Is the layer forward ?
Definition: TrackerLayer.h:70
double fudgeFactor(unsigned iFudge) const
Definition: TrackerLayer.h:108
double R() const
vertex radius
Definition: RawParticle.h:278
double fudgeMax(unsigned iFudge) const
Definition: TrackerLayer.h:105
double Z() const
z of vertex
Definition: RawParticle.h:276
unsigned int fudgeNumber() const
Set a fudge factor for material inhomogeneities in this layer.
Definition: TrackerLayer.h:101
const BoundSurface & surface() const
Returns the surface.
Definition: TrackerLayer.h:73
std::pair< OmniClusterRef, TrackingParticleRef > P
GlobalVector theNormalVector
double fudgeMin(unsigned iFudge) const
Definition: TrackerLayer.h:102
const MediumProperties & mediumProperties() const
Definition: Surface.h:112
void MaterialEffects::save ( )

Save nuclear interaction information.

Definition at line 477 of file MaterialEffects.cc.

References NuclearInteraction, and MaterialEffectsSimulator::save().

Referenced by Vispa.Main.TabController.TabController::allowClose(), Vispa.Main.TabController.TabController::checkModificationTimestamp(), and TrajectoryManager::reconstruct().

477  {
478 
479  // Save current nuclear interactions in the event libraries.
481 
482 }
MaterialEffectsSimulator * NuclearInteraction
virtual void save()
Used by NuclearInteractionSimulator to save last sampled event.
double MaterialEffects::thickness ( ) const
inline

Return the thickness of the current layer.

Definition at line 74 of file MaterialEffects.h.

Referenced by TrajectoryManager::createPSimHits().

74 { return theThickness; }

Member Data Documentation

BremsstrahlungSimulator* MaterialEffects::Bremsstrahlung
private

Definition at line 107 of file MaterialEffects.h.

Referenced by interact(), MaterialEffects(), and ~MaterialEffects().

EnergyLossSimulator* MaterialEffects::EnergyLoss
private

Definition at line 111 of file MaterialEffects.h.

Referenced by interact(), MaterialEffects(), and ~MaterialEffects().

MultipleScatteringSimulator* MaterialEffects::MultipleScattering
private

Definition at line 110 of file MaterialEffects.h.

Referenced by interact(), MaterialEffects(), and ~MaterialEffects().

MuonBremsstrahlungSimulator* MaterialEffects::MuonBremsstrahlung
private

Definition at line 109 of file MaterialEffects.h.

Referenced by interact(), MaterialEffects(), and ~MaterialEffects().

MaterialEffectsSimulator* MaterialEffects::NuclearInteraction
private

Definition at line 112 of file MaterialEffects.h.

Referenced by interact(), MaterialEffects(), save(), and ~MaterialEffects().

PairProductionSimulator* MaterialEffects::PairProduction
private

Definition at line 106 of file MaterialEffects.h.

Referenced by interact(), MaterialEffects(), and ~MaterialEffects().

double MaterialEffects::pTmin
private

Definition at line 115 of file MaterialEffects.h.

Referenced by interact(), and MaterialEffects().

double MaterialEffects::theEnergyLoss
private

Definition at line 118 of file MaterialEffects.h.

Referenced by interact().

GlobalVector MaterialEffects::theNormalVector
private

Definition at line 116 of file MaterialEffects.h.

Referenced by interact(), and radLengths().

double MaterialEffects::theTECFudgeFactor
private

Definition at line 119 of file MaterialEffects.h.

Referenced by interact(), and MaterialEffects().

double MaterialEffects::theThickness
private

Definition at line 117 of file MaterialEffects.h.

Referenced by radLengths().

bool MaterialEffects::use_hardcoded
private

Definition at line 124 of file MaterialEffects.h.

Referenced by interact(), and MaterialEffects().