CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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
 
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 funct::A, Bremsstrahlung, HLT_FULL_cff::distCut, EnergyLoss, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), GeV, mps_fire::i, writeEcalDQMStatus::inputFile, dqmiolumiharvest::j, MultipleScattering, MuonBremsstrahlung, NuclearInteraction, PairProduction, pTmin, plotFactory::ratios, AlCaHLTBitMon_QueryRunRegistry::string, theTECFudgeFactor, use_hardcoded, and BeamSpotPI::Z.

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

Default destructor.

Definition at line 197 of file MaterialEffects.cc.

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

197  {
198  if (PairProduction)
199  delete PairProduction;
200  if (Bremsstrahlung)
201  delete Bremsstrahlung;
202  if (EnergyLoss)
203  delete EnergyLoss;
204  if (MultipleScattering)
205  delete MultipleScattering;
206  if (NuclearInteraction)
207  delete NuclearInteraction;
208  //Muon Brem
209  if (MuonBremsstrahlung)
210  delete MuonBremsstrahlung;
211 }
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 74 of file MaterialEffects.h.

References theEnergyLoss.

Referenced by TrajectoryManager::createPSimHits().

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

Return the Energy Loss engine.

Definition at line 80 of file MaterialEffects.h.

References EnergyLoss.

Referenced by CalorimetryManager::MuonMipSimulation().

80 { return EnergyLoss; }
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 213 of file MaterialEffects.cc.

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

Referenced by TrajectoryManager::reconstruct().

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

Return the Multiple Scattering engine.

Definition at line 77 of file MaterialEffects.h.

References MultipleScattering.

77 { return MultipleScattering; }
MultipleScatteringSimulator * MultipleScattering
MuonBremsstrahlungSimulator* MaterialEffects::muonBremsstrahlungSimulator ( ) const
inline

Return the Muon Bremsstrahlung engine.

Definition at line 83 of file MaterialEffects.h.

References MuonBremsstrahlung.

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

The vector normal to the surface traversed.

Definition at line 416 of file MaterialEffects.cc.

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

Referenced by interact().

416  {
417  return layer.forward() ? layer.disk()->normalVector()
418  : GlobalVector(myTrack.particle().X(), myTrack.particle().Y(), 0.) / myTrack.particle().R();
419 }
bool forward() const
Is the layer forward ?
Definition: TrackerLayer.h:66
RawParticle const & particle() const
The particle being propagated.
double R() const
vertex radius
Definition: RawParticle.h:290
double Y() const
y of vertex
Definition: RawParticle.h:287
BoundDisk const * disk() const
Returns the surface.
Definition: TrackerLayer.h:75
double X() const
x of vertex
Definition: RawParticle.h:286
Global3DVector GlobalVector
Definition: GlobalVector.h:10
double MaterialEffects::radLengths ( const TrackerLayer layer,
ParticlePropagator myTrack 
)
private

The number of radiation lengths traversed.

Definition at line 381 of file MaterialEffects.cc.

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

Referenced by interact().

381  {
382  // Thickness of layer
384 
385  GlobalVector P(myTrack.particle().Px(), myTrack.particle().Py(), myTrack.particle().Pz());
386 
387  // Effective length of track inside layer (considering crossing angle)
388  // double radlen = theThickness / fabs(P.dot(theNormalVector)/(P.mag()*theNormalVector.mag()));
389  double radlen = theThickness / fabs(P.dot(theNormalVector)) * P.mag();
390 
391  // This is a series of fudge factors (from the geometry description),
392  // to describe the layer inhomogeneities (services, cables, supports...)
393  double rad = myTrack.particle().R();
394  double zed = fabs(myTrack.particle().Z());
395 
396  double factor = 1;
397 
398  // Are there fudge factors for this layer
399  if (layer.fudgeNumber())
400 
401  // If yes, loop on them
402  for (unsigned int iLayer = 0; iLayer < layer.fudgeNumber(); ++iLayer) {
403  // Apply to R if forward layer, to Z if barrel layer
404  if ((layer.forward() && layer.fudgeMin(iLayer) < rad && rad < layer.fudgeMax(iLayer)) ||
405  (!layer.forward() && layer.fudgeMin(iLayer) < zed && zed < layer.fudgeMax(iLayer))) {
406  factor = layer.fudgeFactor(iLayer);
407  break;
408  }
409  }
410 
411  theThickness *= factor;
412 
413  return radlen * factor;
414 }
float radLen() const
bool forward() const
Is the layer forward ?
Definition: TrackerLayer.h:66
double fudgeFactor(unsigned iFudge) const
Definition: TrackerLayer.h:104
RawParticle const & particle() const
The particle being propagated.
double R() const
vertex radius
Definition: RawParticle.h:290
double fudgeMax(unsigned iFudge) const
Definition: TrackerLayer.h:101
double Py() const
y of the momentum
Definition: RawParticle.h:300
double Z() const
z of vertex
Definition: RawParticle.h:288
double Pz() const
z of the momentum
Definition: RawParticle.h:303
unsigned int fudgeNumber() const
Set a fudge factor for material inhomogeneities in this layer.
Definition: TrackerLayer.h:97
const BoundSurface & surface() const
Returns the surface.
Definition: TrackerLayer.h:69
std::pair< OmniClusterRef, TrackingParticleRef > P
double Px() const
x of the momentum
Definition: RawParticle.h:297
GlobalVector theNormalVector
double fudgeMin(unsigned iFudge) const
Definition: TrackerLayer.h:98
const MediumProperties & mediumProperties() const
Definition: Surface.h:83
void MaterialEffects::save ( )

Save nuclear interaction information.

Definition at line 421 of file MaterialEffects.cc.

References NuclearInteraction, and MaterialEffectsSimulator::save().

Referenced by TrajectoryManager::reconstruct().

421  {
422  // Save current nuclear interactions in the event libraries.
423  if (NuclearInteraction)
425 }
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 71 of file MaterialEffects.h.

References theThickness.

Referenced by TrajectoryManager::createPSimHits().

71 { return theThickness; }

Member Data Documentation

BremsstrahlungSimulator* MaterialEffects::Bremsstrahlung
private

Definition at line 94 of file MaterialEffects.h.

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

EnergyLossSimulator* MaterialEffects::EnergyLoss
private

Definition at line 98 of file MaterialEffects.h.

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

MultipleScatteringSimulator* MaterialEffects::MultipleScattering
private
MuonBremsstrahlungSimulator* MaterialEffects::MuonBremsstrahlung
private
MaterialEffectsSimulator* MaterialEffects::NuclearInteraction
private

Definition at line 99 of file MaterialEffects.h.

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

PairProductionSimulator* MaterialEffects::PairProduction
private

Definition at line 93 of file MaterialEffects.h.

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

double MaterialEffects::pTmin
private

Definition at line 102 of file MaterialEffects.h.

Referenced by interact(), and MaterialEffects().

double MaterialEffects::theEnergyLoss
private

Definition at line 105 of file MaterialEffects.h.

Referenced by energyLoss(), and interact().

GlobalVector MaterialEffects::theNormalVector
private

Definition at line 103 of file MaterialEffects.h.

Referenced by interact(), and radLengths().

double MaterialEffects::theTECFudgeFactor
private

Definition at line 106 of file MaterialEffects.h.

Referenced by interact(), and MaterialEffects().

double MaterialEffects::theThickness
private

Definition at line 104 of file MaterialEffects.h.

Referenced by radLengths(), and thickness().

bool MaterialEffects::use_hardcoded
private

Definition at line 111 of file MaterialEffects.h.

Referenced by interact(), and MaterialEffects().