CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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
 
NuclearInteractionSimulatorNuclearInteraction
 
PairProductionSimulatorPairProduction
 
double pTmin
 
double theEnergyLoss
 
GlobalVector theNormalVector
 
double theTECFudgeFactor
 
double theThickness
 
bool use_hardcoded
 

Detailed Description

Definition at line 50 of file MaterialEffects.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 25 of file MaterialEffects.cc.

References funct::A, Bremsstrahlung, EnergyLoss, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), i, analyzePatCleaning_cfg::inputFile, j, MultipleScattering, MuonBremsstrahlung, NuclearInteraction, PairProduction, pTmin, plotFactory::ratios, AlCaHLTBitMon_QueryRunRegistry::string, theTECFudgeFactor, use_hardcoded, and Gflash::Z.

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

Default destructor.

Definition at line 211 of file MaterialEffects.cc.

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

211  {
212 
213  if ( PairProduction ) delete PairProduction;
214  if ( Bremsstrahlung ) delete Bremsstrahlung;
215  if ( EnergyLoss ) delete EnergyLoss;
218 //Muon Brem
220 }
PairProductionSimulator * PairProduction
MuonBremsstrahlungSimulator * MuonBremsstrahlung
BremsstrahlungSimulator * Bremsstrahlung
EnergyLossSimulator * EnergyLoss
NuclearInteractionSimulator * NuclearInteraction
MultipleScatteringSimulator * MultipleScattering

Member Function Documentation

double MaterialEffects::energyLoss ( ) const
inline

Return the energy loss by ionization in the current layer.

Definition at line 76 of file MaterialEffects.h.

References theEnergyLoss.

Referenced by TrajectoryManager::createPSimHits().

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

Return the Energy Loss engine.

Definition at line 84 of file MaterialEffects.h.

References EnergyLoss.

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

84  {
85  return EnergyLoss;
86  }
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 222 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().

226  {
227 
229  double radlen;
230  theEnergyLoss = 0;
231  theNormalVector = normalVector(layer,myTrack);
232  radlen = radLengths(layer,myTrack);
233 
234 //-------------------
235 // Photon Conversion
236 //-------------------
237 
238  if ( PairProduction && myTrack.pid()==22 ) {
239 
240  //
241  PairProduction->updateState(myTrack, radlen, random);
242 
243  if ( PairProduction->nDaughters() ) {
244  //add a vertex to the mother particle
245  int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
247 
248  // Check if it is a valid vertex first:
249  if (ivertex>=0) {
250  // This was a photon that converted
251  for ( DaughterIter = PairProduction->beginDaughters();
252  DaughterIter != PairProduction->endDaughters();
253  ++DaughterIter) {
254 
255  mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
256 
257  }
258  // The photon converted. Return.
259  return;
260  }
261  else {
262  edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in photon conv. -> " << ivertex << std::endl;
263  }
264 
265  }
266 
267  }
268 
269  if ( myTrack.pid() == 22 ) return;
270 
271 //------------------------
272 // Nuclear interactions
273 //------------------------
274 
275  if ( NuclearInteraction && abs(myTrack.pid()) > 100
276  && abs(myTrack.pid()) < 1000000) {
277 
278  // Simulate a nuclear interaction
279  double factor = 1.0;
280  if(use_hardcoded){
281  if (layer.layerNumber() >= 19 && layer.layerNumber() <= 27 )
282  factor = theTECFudgeFactor;
283  }
284  NuclearInteraction->updateState(myTrack, radlen*factor, random);
285 
286  if ( NuclearInteraction->nDaughters() ) {
287 
288  //add a end vertex to the mother particle
289  int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
291 
292  // Check if it is a valid vertex first:
293  if (ivertex>=0) {
294  // This was a hadron that interacted inelastically
295  int idaugh = 0;
296  for ( DaughterIter = NuclearInteraction->beginDaughters();
297  DaughterIter != NuclearInteraction->endDaughters();
298  ++DaughterIter) {
299 
300  // The daughter in the event
301  int daughId = mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
302 
303  // Store the closest daughter in the mother info (for later tracking purposes)
304  if ( NuclearInteraction->closestDaughterId() == idaugh++ ) {
305  if ( mySimEvent.track(itrack).vertex().position().Pt() < 4.0 )
306  mySimEvent.track(itrack).setClosestDaughterId(daughId);
307  }
308 
309  }
310  // The hadron is destroyed. Return.
311  return;
312  }
313  else {
314  edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in nucl. int. -> " << ivertex << std::endl;
315  }
316 
317  }
318 
319  }
320 
321  if ( myTrack.charge() == 0 ) return;
322 
324 
325 //----------------
326 // Bremsstrahlung
327 //----------------
328 
329  if ( Bremsstrahlung && abs(myTrack.pid())==11 ) {
330 
331  Bremsstrahlung->updateState(myTrack,radlen, random);
332 
333  if ( Bremsstrahlung->nDaughters() ) {
334 
335  // Add a vertex, but do not attach it to the electron, because it
336  // continues its way...
337  int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
339 
340  // Check if it is a valid vertex first:
341  if (ivertex>=0) {
342  for ( DaughterIter = Bremsstrahlung->beginDaughters();
343  DaughterIter != Bremsstrahlung->endDaughters();
344  ++DaughterIter) {
345  mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
346  }
347  }
348  else {
349  edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in brem -> " << ivertex << std::endl;
350  }
351 
352  }
353 
354  }
355 
356 //---------------------------
357 // Muon_Bremsstrahlung
358 //--------------------------
359 
360  if ( MuonBremsstrahlung && abs(myTrack.pid())==13 ) {
361 
362  MuonBremsstrahlung->updateState(myTrack, radlen, random);
363 
364  if ( MuonBremsstrahlung->nDaughters() ) {
365 
366  // Add a vertex, but do not attach it to the muon, because it
367  // continues its way...
368  int ivertex = mySimEvent.addSimVertex(myTrack.vertex(),itrack,
370 
371  // Check if it is a valid vertex first:
372  if (ivertex>=0) {
373  for ( DaughterIter = MuonBremsstrahlung->beginDaughters();
374  DaughterIter != MuonBremsstrahlung->endDaughters();
375  ++DaughterIter) {
376  mySimEvent.addSimTrack(&(*DaughterIter), ivertex);
377  }
378  }
379  else {
380  edm::LogWarning("MaterialEffects") << " WARNING: A non valid vertex was found in muon brem -> " << ivertex << std::endl;
381  }
382 
383  }
384 
385  }
386 
390 
391  if ( EnergyLoss )
392  {
393  theEnergyLoss = myTrack.E();
394  EnergyLoss->updateState(myTrack, radlen, random);
395  theEnergyLoss -= myTrack.E();
396  }
397 
398 
402 
403  if ( MultipleScattering && myTrack.Pt() > pTmin ) {
404  // MultipleScattering->setNormalVector(normalVector(layer,myTrack));
406  MultipleScattering->updateState(myTrack,radlen, random);
407  }
408 
409 }
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
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
NuclearInteractionSimulator * NuclearInteraction
MultipleScatteringSimulator * MultipleScattering
MultipleScatteringSimulator* MaterialEffects::multipleScatteringSimulator ( ) const
inline

Return the Multiple Scattering engine.

Definition at line 79 of file MaterialEffects.h.

References MultipleScattering.

Referenced by MuonSimHitProducer::applyMaterialEffects().

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

Return the Muon Bremsstrahlung engine.

Definition at line 89 of file MaterialEffects.h.

References MuonBremsstrahlung.

Referenced by MuonSimHitProducer::applyMaterialEffects().

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

The vector normal to the surface traversed.

Definition at line 453 of file MaterialEffects.cc.

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

Referenced by interact().

454  {
455  return layer.forward() ?
456  layer.disk()->normalVector() :
457  GlobalVector(myTrack.X(),myTrack.Y(),0.)/myTrack.R();
458 }
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 412 of file MaterialEffects.cc.

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

Referenced by interact().

413  {
414 
415  // Thickness of layer
417 
418  GlobalVector P(myTrack.Px(),myTrack.Py(),myTrack.Pz());
419 
420  // Effective length of track inside layer (considering crossing angle)
421  // double radlen = theThickness / fabs(P.dot(theNormalVector)/(P.mag()*theNormalVector.mag()));
422  double radlen = theThickness / fabs(P.dot(theNormalVector)) * P.mag();
423 
424  // This is a series of fudge factors (from the geometry description),
425  // to describe the layer inhomogeneities (services, cables, supports...)
426  double rad = myTrack.R();
427  double zed = fabs(myTrack.Z());
428 
429  double factor = 1;
430 
431  // Are there fudge factors for this layer
432  if ( layer.fudgeNumber() )
433 
434  // If yes, loop on them
435  for ( unsigned int iLayer=0; iLayer < layer.fudgeNumber(); ++iLayer ) {
436 
437  // Apply to R if forward layer, to Z if barrel layer
438  if ( ( layer.forward() && layer.fudgeMin(iLayer) < rad && rad < layer.fudgeMax(iLayer) ) ||
439  ( !layer.forward() && layer.fudgeMin(iLayer) < zed && zed < layer.fudgeMax(iLayer) ) ) {
440  factor = layer.fudgeFactor(iLayer);
441  break;
442  }
443 
444  }
445 
446  theThickness *= factor;
447 
448  return radlen * factor;
449 
450 }
float radLen() const
bool forward() const
Is the layer forward ?
Definition: TrackerLayer.h:70
double fudgeFactor(unsigned iFudge) const
Definition: TrackerLayer.h:108
#define P
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
GlobalVector theNormalVector
double fudgeMin(unsigned iFudge) const
Definition: TrackerLayer.h:102
const MediumProperties & mediumProperties() const
Definition: Surface.h:120
void MaterialEffects::save ( )

Save nuclear interaction information.

Definition at line 461 of file MaterialEffects.cc.

References NuclearInteraction, and NuclearInteractionSimulator::save().

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

461  {
462 
463  // Save current nuclear interactions in the event libraries.
465 
466 }
NuclearInteractionSimulator * NuclearInteraction
void save()
Save current nuclear interaction (for later use)
double MaterialEffects::thickness ( ) const
inline

Return the thickness of the current layer.

Definition at line 73 of file MaterialEffects.h.

References theThickness.

Referenced by TrajectoryManager::createPSimHits().

73 { return theThickness; }

Member Data Documentation

BremsstrahlungSimulator* MaterialEffects::Bremsstrahlung
private

Definition at line 106 of file MaterialEffects.h.

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

EnergyLossSimulator* MaterialEffects::EnergyLoss
private
MultipleScatteringSimulator* MaterialEffects::MultipleScattering
private
MuonBremsstrahlungSimulator* MaterialEffects::MuonBremsstrahlung
private
NuclearInteractionSimulator* MaterialEffects::NuclearInteraction
private

Definition at line 111 of file MaterialEffects.h.

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

PairProductionSimulator* MaterialEffects::PairProduction
private

Definition at line 105 of file MaterialEffects.h.

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

double MaterialEffects::pTmin
private

Definition at line 114 of file MaterialEffects.h.

Referenced by interact(), and MaterialEffects().

double MaterialEffects::theEnergyLoss
private

Definition at line 117 of file MaterialEffects.h.

Referenced by energyLoss(), and interact().

GlobalVector MaterialEffects::theNormalVector
private

Definition at line 115 of file MaterialEffects.h.

Referenced by interact(), and radLengths().

double MaterialEffects::theTECFudgeFactor
private

Definition at line 118 of file MaterialEffects.h.

Referenced by interact(), and MaterialEffects().

double MaterialEffects::theThickness
private

Definition at line 116 of file MaterialEffects.h.

Referenced by radLengths(), and thickness().

bool MaterialEffects::use_hardcoded
private

Definition at line 123 of file MaterialEffects.h.

Referenced by interact(), and MaterialEffects().