CMS 3D CMS Logo

MonopoleTransportation.cc
Go to the documentation of this file.
1 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
2 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
3 //
4 // This class is a process responsible for the transportation of
5 // magnetic monopoles, ie the geometrical propagation that encounters the
6 // geometrical sub-volumes of the detectors.
7 //
8 // =======================================================================
9 // Created: 3 May 2010, J. Apostolakis, B. Bozsogi
10 // G4MonopoleTransportation class for
11 // Geant4 extended example "monopole"
12 //
13 // Adopted for CMSSW by V.Ivanchenko 30 April 2018
14 //
15 // =======================================================================
16 
20 
21 #include "G4ProductionCutsTable.hh"
22 #include "G4ParticleTable.hh"
23 #include "G4ChordFinder.hh"
24 #include "G4SafetyHelper.hh"
25 #include "G4FieldManagerStore.hh"
26 #include "G4TransportationProcessType.hh"
27 #include "G4SystemOfUnits.hh"
28 
29 class G4VSensitiveDetector;
30 
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 
34  : G4VProcess(G4String("MonopoleTransportation"), fTransportation),
35  fParticleDef(mpl),
36  fieldMgrCMS(nullptr),
37  fLinearNavigator(nullptr),
38  fFieldPropagator(nullptr),
39  fParticleIsLooping(false),
40  fPreviousSftOrigin(0., 0., 0.),
41  fPreviousSafety(0.0),
42  fThreshold_Warning_Energy(100 * MeV),
43  fThreshold_Important_Energy(250 * MeV),
44  fThresholdTrials(10),
45  fNoLooperTrials(0),
46  fSumEnergyKilled(0.0),
47  fMaxEnergyKilled(0.0),
48  fShortStepOptimisation(false), // Old default: true (=fast short steps)
49  fpSafetyHelper(nullptr) {
50  verboseLevel = verb;
51 
52  // set Process Sub Type
53  SetProcessSubType(TRANSPORTATION);
54 
55 #ifdef G4MULTITHREADED
56  // Do not finalize the MonopoleTransportation class
57  if (G4Threading::IsMasterThread()) {
58  return;
59  }
60 #endif
61 
62  G4TransportationManager* transportMgr = G4TransportationManager::GetTransportationManager();
63 
64  fLinearNavigator = transportMgr->GetNavigatorForTracking();
65 
66  fFieldPropagator = transportMgr->GetPropagatorInField();
67  fpSafetyHelper = transportMgr->GetSafetyHelper();
68 
69  // Cannot determine whether a field exists here,
70  // because it would only work if the field manager has informed
71  // about the detector's field before this transportation process
72  // is constructed.
73  // Instead later the method DoesGlobalFieldExist() is called
74 
75  fCurrentTouchableHandle = nullptr;
76 
77  fEndGlobalTimeComputed = false;
79 }
80 
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82 
84  if ((verboseLevel > 0) && (fSumEnergyKilled > 0.0)) {
85  /*
86  G4cout << " MonopoleTransportation: Statistics for looping particles "
87  << G4endl;
88  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
89  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
90  */
91  }
92 }
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95 //
96 // Responsibilities:
97 // Find whether the geometry limits the Step, and to what length
98 // Calculate the new value of the safety and return it.
99 // Store the final time, position and momentum.
100 
102  G4double, // previousStepSize
103  G4double currentMinimumStep,
104  G4double& currentSafety,
105  G4GPILSelection* selection) {
106  // change to monopole equation
108 
109  G4double geometryStepLength, newSafety;
110  fParticleIsLooping = false;
111 
112  // Initial actions moved to StartTrack()
113  // --------------------------------------
114  // Note: in case another process changes touchable handle
115  // it will be necessary to add here (for all steps)
116  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
117 
118  // GPILSelection is set to defaule value of CandidateForSelection
119  // It is a return value
120  //
121  *selection = CandidateForSelection;
122 
123  // Get initial Energy/Momentum of the track
124  //
125  const G4DynamicParticle* pParticle = track.GetDynamicParticle();
126  const G4ThreeVector& startMomentumDir = pParticle->GetMomentumDirection();
127  G4ThreeVector startPosition = track.GetPosition();
128 
129  // The Step Point safety can be limited by other geometries and/or the
130  // assumptions of any process - it's not always the geometrical safety.
131  // We calculate the starting point's isotropic safety here.
132  //
133  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin;
134  G4double MagSqShift = OriginShift.mag2();
135  if (MagSqShift >= sqr(fPreviousSafety)) {
136  currentSafety = 0.0;
137  } else {
138  currentSafety = fPreviousSafety - std::sqrt(MagSqShift);
139  }
140 
141  // Is the monopole charged ?
142  //
143  G4double particleMagneticCharge = fParticleDef->MagneticCharge();
144  G4double particleElectricCharge = pParticle->GetCharge();
145 
146  fGeometryLimitedStep = false;
147 
148  // There is no need to locate the current volume. It is Done elsewhere:
149  // On track construction
150  // By the tracking, after all AlongStepDoIts, in "Relocation"
151 
152  // Check whether the particle have an (EM) field force exerting upon it
153  G4bool fieldExertsForce = false;
154 
155  if (particleMagneticCharge != 0.0 && fieldMgrCMS) {
156  // Message the field Manager, to configure it for this track
158  // Moved here, in order to allow a transition
159  // from a zero-field status (with fieldMgr->(field)0
160  // to a finite field status
161 
162  // If the field manager has no field, there is no field !
163  fieldExertsForce = (fieldMgrCMS->GetDetectorField() != nullptr);
164  }
165 
166  // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
167  // << " fieldMgr= " << fieldMgr << G4endl;
168 
169  // Choose the calculation of the transportation: Field or not
170  //
171  if (!fieldExertsForce) {
172  G4double linearStepLength;
173  if (fShortStepOptimisation && (currentMinimumStep <= currentSafety)) {
174  // The Step is guaranteed to be taken
175  //
176  geometryStepLength = currentMinimumStep;
177  fGeometryLimitedStep = false;
178  } else {
179  // Find whether the straight path intersects a volume
180  //
181  linearStepLength = fLinearNavigator->ComputeStep(startPosition, startMomentumDir, currentMinimumStep, newSafety);
182  // Remember last safety origin & value.
183  //
184  fPreviousSftOrigin = startPosition;
185  fPreviousSafety = newSafety;
186  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
187 
188  // The safety at the initial point has been re-calculated:
189  //
190  currentSafety = newSafety;
191 
192  fGeometryLimitedStep = (linearStepLength <= currentMinimumStep);
193  if (fGeometryLimitedStep) {
194  // The geometry limits the Step size (an intersection was found.)
195  geometryStepLength = linearStepLength;
196  } else {
197  // The full Step is taken.
198  geometryStepLength = currentMinimumStep;
199  }
200  }
201  endpointDistance = geometryStepLength;
202 
203  // Calculate final position
204  //
205  fTransportEndPosition = startPosition + geometryStepLength * startMomentumDir;
206 
207  // Momentum direction, energy and polarisation are unchanged by transport
208  //
209  fTransportEndMomentumDir = startMomentumDir;
210  fTransportEndKineticEnergy = track.GetKineticEnergy();
211  fTransportEndSpin = track.GetPolarization();
212  fParticleIsLooping = false;
213  fMomentumChanged = false;
214  fEndGlobalTimeComputed = false;
215  } else // A field exerts force
216  {
217  G4double momentumMagnitude = pParticle->GetTotalMomentum();
218  G4ThreeVector EndUnitMomentum;
219  G4double lengthAlongCurve;
220  G4double restMass = fParticleDef->GetPDGMass();
221 
222  G4ChargeState chargeState(particleElectricCharge, // The charge can change (dynamic)
223  fParticleDef->GetPDGSpin(),
224  0, // Magnetic moment: pParticleDef->GetMagneticMoment(),
225  0, // Electric Dipole moment - not in Particle Definition
226  particleMagneticCharge); // in Mev/c
227 
228  G4EquationOfMotion* equationOfMotion =
229  fFieldPropagator->GetChordFinder()->GetIntegrationDriver()->GetEquationOfMotion();
230 
231  equationOfMotion->SetChargeMomentumMass(chargeState, // Was particleMagneticCharge - in Mev/c
232  momentumMagnitude, // Was particleElectricCharge
233  restMass);
234  // SetChargeMomentumMass now passes both the electric and magnetic charge - in chargeState
235 
236  G4ThreeVector spin = track.GetPolarization();
237  G4FieldTrack aFieldTrack = G4FieldTrack(startPosition,
238  track.GetMomentumDirection(),
239  0.0,
240  track.GetKineticEnergy(),
241  restMass,
242  track.GetVelocity(),
243  track.GetGlobalTime(), // Lab.
244  track.GetProperTime(), // Part.
245  &spin);
246  if (currentMinimumStep > 0) {
247  // Do the Transport in the field (non recti-linear)
248  //
249  lengthAlongCurve =
250  fFieldPropagator->ComputeStep(aFieldTrack, currentMinimumStep, currentSafety, track.GetVolume());
251  fGeometryLimitedStep = lengthAlongCurve < currentMinimumStep;
252  if (fGeometryLimitedStep) {
253  geometryStepLength = lengthAlongCurve;
254  } else {
255  geometryStepLength = currentMinimumStep;
256  }
257  } else {
258  geometryStepLength = lengthAlongCurve = 0.0;
259  fGeometryLimitedStep = false;
260  }
261 
262  // Remember last safety origin & value.
263  //
264  fPreviousSftOrigin = startPosition;
265  fPreviousSafety = currentSafety;
266  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
267 
268  // Get the End-Position and End-Momentum (Dir-ection)
269  //
270  fTransportEndPosition = aFieldTrack.GetPosition();
271 
272  // Momentum: Magnitude and direction can be changed too now ...
273  //
274  fMomentumChanged = true;
275  fTransportEndMomentumDir = aFieldTrack.GetMomentumDir();
276 
277  fTransportEndKineticEnergy = aFieldTrack.GetKineticEnergy();
278 
279  fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
280  fEndGlobalTimeComputed = true;
281 
282  fTransportEndSpin = aFieldTrack.GetSpin();
283  fParticleIsLooping = fFieldPropagator->IsParticleLooping();
284  endpointDistance = (fTransportEndPosition - startPosition).mag();
285  }
286 
287  // If we are asked to go a step length of 0, and we are on a boundary
288  // then a boundary will also limit the step -> we must flag this.
289  //
290  if (currentMinimumStep == 0.0) {
291  if (currentSafety == 0.0)
292  fGeometryLimitedStep = true;
293  }
294 
295  // Update the safety starting from the end-point,
296  // if it will become negative at the end-point.
297  //
298  if (currentSafety < endpointDistance) {
299  if (particleMagneticCharge != 0.0) {
300  G4double endSafety = fLinearNavigator->ComputeSafety(fTransportEndPosition);
301  currentSafety = endSafety;
303  fPreviousSafety = currentSafety;
304  fpSafetyHelper->SetCurrentSafety(currentSafety, fTransportEndPosition);
305 
306  // Because the Stepping Manager assumes it is from the start point,
307  // add the StepLength
308  //
309  currentSafety += endpointDistance;
310 
311 #ifdef G4DEBUG_TRANSPORT
312  G4cout.precision(12);
313  G4cout << "***MonopoleTransportation::AlongStepGPIL ** " << G4endl;
314  G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
315  << " and it returned safety= " << endSafety << G4endl;
316  G4cout << " Adding endpoint distance " << endpointDistance << " to obtain pseudo-safety= " << currentSafety
317  << G4endl;
318 #endif
319  }
320  }
321 
322  fParticleChange.ProposeTrueStepLength(geometryStepLength);
323 
324  return geometryStepLength;
325 }
326 
327 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
328 //
329 // Initialize ParticleChange (by setting all its members equal
330 // to corresponding members in G4Track)
331 
332 G4VParticleChange* MonopoleTransportation::AlongStepDoIt(const G4Track& track, const G4Step& stepData) {
333  fParticleChange.Initialize(track);
334 
335  // Code for specific process
336  //
337  fParticleChange.ProposePosition(fTransportEndPosition);
338  fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir);
340  fParticleChange.SetMomentumChanged(fMomentumChanged);
341 
342  fParticleChange.ProposePolarization(fTransportEndSpin);
343 
344  G4double deltaTime = 0.0;
345 
346  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
347  G4double startTime = track.GetGlobalTime();
348 
349  if (!fEndGlobalTimeComputed) {
350  // The time was not integrated .. make the best estimate possible
351  //
352  G4double finalVelocity = track.GetVelocity();
353  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity();
354  G4double stepLength = track.GetStepLength();
355 
356  deltaTime = 0.0; // in case initialVelocity = 0
357  if (finalVelocity > 0.0) {
358  G4double meanInverseVelocity;
359  meanInverseVelocity = 0.5 * (1.0 / initialVelocity + 1.0 / finalVelocity);
360  deltaTime = stepLength * meanInverseVelocity;
361  } else if (initialVelocity > 0.0) {
362  deltaTime = stepLength / initialVelocity;
363  }
364  fCandidateEndGlobalTime = startTime + deltaTime;
365  } else {
366  deltaTime = fCandidateEndGlobalTime - startTime;
367  }
368 
369  fParticleChange.ProposeGlobalTime(fCandidateEndGlobalTime);
370 
371  // Now Correct by Lorentz factor to get "proper" deltaTime
372 
373  G4double restMass = track.GetDynamicParticle()->GetMass();
374  G4double deltaProperTime = deltaTime * (restMass / track.GetTotalEnergy());
375 
376  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime);
377  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
378 
379  // If the particle is caught looping or is stuck (in very difficult
380  // boundaries) in a magnetic field (doing many steps)
381  // THEN this kills it ...
382  //
383  if (fParticleIsLooping) {
384  G4double endEnergy = fTransportEndKineticEnergy;
385 
387  // Kill the looping particle
388  //
389  fParticleChange.ProposeTrackStatus(fStopAndKill);
390 
391  // 'Bare' statistics
392  fSumEnergyKilled += endEnergy;
393  if (endEnergy > fMaxEnergyKilled) {
394  fMaxEnergyKilled = endEnergy;
395  }
396 
397 #ifdef G4VERBOSE
398  if ((verboseLevel > 1) || (endEnergy > fThreshold_Warning_Energy)) {
399  G4cout << " MonopoleTransportation is killing track that is looping or stuck " << G4endl << " This track has "
400  << track.GetKineticEnergy() / MeV << " MeV energy." << G4endl;
401  G4cout << " Number of trials = " << fNoLooperTrials << " No of calls to AlongStepDoIt = " << noCalls
402  << G4endl;
403  }
404 #endif
405  fNoLooperTrials = 0;
406  } else {
407  fNoLooperTrials++;
408 #ifdef G4VERBOSE
409  if ((verboseLevel > 2)) {
410  G4cout << " MonopoleTransportation::AlongStepDoIt(): Particle looping - "
411  << " Number of trials = " << fNoLooperTrials << " No of calls to = " << noCalls << G4endl;
412  }
413 #endif
414  }
415  } else {
416  fNoLooperTrials = 0;
417  }
418 
419  // Another (sometimes better way) is to use a user-limit maximum Step size
420  // to alleviate this problem ..
421 
422  // Introduce smooth curved trajectories to particle-change
423  //
424  fParticleChange.SetPointerToVectorOfAuxiliaryPoints(fFieldPropagator->GimmeTrajectoryVectorAndForgetIt());
425 
426  return &fParticleChange;
427 }
428 
429 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
430 //
431 // This ensures that the PostStep action is always called,
432 // so that it can do the relocation if it is needed.
433 //
434 
436  G4double, // previousStepSize
437  G4ForceCondition* pForceCond) {
438  *pForceCond = Forced;
439  return DBL_MAX; // was kInfinity ; but convention now is DBL_MAX
440 }
441 
442 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
443 
444 G4VParticleChange* MonopoleTransportation::PostStepDoIt(const G4Track& track, const G4Step&) {
445  G4TouchableHandle retCurrentTouchable; // The one to return
446 
447  fParticleChange.ProposeTrackStatus(track.GetTrackStatus());
448 
449  // If the Step was determined by the volume boundary,
450  // logically relocate the particle
451 
452  if (fGeometryLimitedStep) {
453  // fCurrentTouchable will now become the previous touchable,
454  // and what was the previous will be freed.
455  // (Needed because the preStepPoint can point to the previous touchable)
456 
457  fLinearNavigator->SetGeometricallyLimitedStep();
458  fLinearNavigator->LocateGlobalPointAndUpdateTouchableHandle(
459  track.GetPosition(), track.GetMomentumDirection(), fCurrentTouchableHandle, true);
460  // Check whether the particle is out of the world volume
461  // If so it has exited and must be killed.
462  //
463  if (fCurrentTouchableHandle->GetVolume() == nullptr) {
464  fParticleChange.ProposeTrackStatus(fStopAndKill);
465  }
466  retCurrentTouchable = fCurrentTouchableHandle;
467  fParticleChange.SetTouchableHandle(fCurrentTouchableHandle);
468  } else // fGeometryLimitedStep is false
469  {
470  // This serves only to move the Navigator's location
471  //
472  fLinearNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
473 
474  // The value of the track's current Touchable is retained.
475  // (and it must be correct because we must use it below to
476  // overwrite the (unset) one in particle change)
477  // It must be fCurrentTouchable too ??
478  //
479  fParticleChange.SetTouchableHandle(track.GetTouchableHandle());
480  retCurrentTouchable = track.GetTouchableHandle();
481  } // endif ( fGeometryLimitedStep )
482 
483  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume();
484  const G4Material* pNewMaterial = nullptr;
485  const G4VSensitiveDetector* pNewSensitiveDetector = nullptr;
486 
487  if (pNewVol != nullptr) {
488  pNewMaterial = pNewVol->GetLogicalVolume()->GetMaterial();
489  pNewSensitiveDetector = pNewVol->GetLogicalVolume()->GetSensitiveDetector();
490  }
491 
492  fParticleChange.SetMaterialInTouchable((G4Material*)pNewMaterial);
493  fParticleChange.SetSensitiveDetectorInTouchable((G4VSensitiveDetector*)pNewSensitiveDetector);
494 
495  const G4MaterialCutsCouple* pNewMaterialCutsCouple = nullptr;
496  if (pNewVol != nullptr) {
497  pNewMaterialCutsCouple = pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
498  }
499 
500  if (pNewVol != nullptr && pNewMaterialCutsCouple != nullptr &&
501  pNewMaterialCutsCouple->GetMaterial() != pNewMaterial) {
502  // for parametrized volume
503  //
504  pNewMaterialCutsCouple = G4ProductionCutsTable::GetProductionCutsTable()->GetMaterialCutsCouple(
505  pNewMaterial, pNewMaterialCutsCouple->GetProductionCuts());
506  }
507  fParticleChange.SetMaterialCutsCoupleInTouchable(pNewMaterialCutsCouple);
508 
509  // temporarily until Get/Set Material of ParticleChange,
510  // and StepPoint can be made const.
511  // Set the touchable in ParticleChange
512  // this must always be done because the particle change always
513  // uses this value to overwrite the current touchable pointer.
514  //
515  fParticleChange.SetTouchableHandle(retCurrentTouchable);
516 
517  // change to normal equation
519 
520  return &fParticleChange;
521 }
522 
523 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
524 
525 // New method takes over the responsibility to reset the state
526 // of MonopoleTransportation object at the start of a new track
527 // or the resumption of a suspended track.
528 
530  // initialise pointer
531  if (!fieldMgrCMS) {
532  G4FieldManager* fieldMgr = fFieldPropagator->FindAndSetFieldManager(aTrack->GetVolume());
533  fieldMgrCMS = static_cast<CMSFieldManager*>(fieldMgr);
534  }
535 
536  G4VProcess::StartTracking(aTrack);
537 
538  // The actions here are those that were taken in AlongStepGPIL
539  // when track.GetCurrentStepNumber()==1
540 
541  // reset safety value and center
542  //
543  fPreviousSafety = 0.0;
544  fPreviousSftOrigin = G4ThreeVector(0., 0., 0.);
545 
546  // reset looping counter -- for motion in field
547  fNoLooperTrials = 0;
548  // Must clear this state .. else it depends on last track's value
549 
550  // ChordFinder reset internal state
551  //
552  if (DoesGlobalFieldExist()) {
553  fFieldPropagator->ClearPropagatorState();
554  // Resets all state of field propagator class (ONLY)
555  // including safety values (in case of overlaps and to wipe for first track).
556 
557  G4ChordFinder* chordF = fFieldPropagator->GetChordFinder();
558  if (chordF)
559  chordF->ResetStepEstimate();
560  }
561 
562  // Make sure to clear the chord finders of all fields (ie managers)
563  G4FieldManagerStore* fieldMgrStore = G4FieldManagerStore::GetInstance();
564  fieldMgrStore->ClearAllChordFindersState();
565 
566  // Update the current touchable handle (from the track's)
567  //
568  fCurrentTouchableHandle = aTrack->GetTouchableHandle();
569 }
570 
571 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
572 
573 G4double MonopoleTransportation::AtRestGetPhysicalInteractionLength(const G4Track&, G4ForceCondition*) { return -1.0; }
574 
575 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
576 
577 G4VParticleChange* MonopoleTransportation::AtRestDoIt(const G4Track&, const G4Step&) { return nullptr; }
578 
579 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
580 
582 // Statistics for tracks killed (currently due to looping in field)
583 {
584  if (report) {
585  G4cout << " MonopoleTransportation: Statistics for looping particles " << G4endl;
586  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
587  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
588  }
589 
590  fSumEnergyKilled = 0;
591  fMaxEnergyKilled = -1.0 * CLHEP::GeV;
592 }
593 
594 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
CMSFieldManager::ConfigureForTrack
void ConfigureForTrack(const G4Track *) override
Definition: CMSFieldManager.cc:140
Monopole.h
MonopoleTransportation::AlongStepDoIt
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData) override
Definition: MonopoleTransportation.cc:332
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11713
MonopoleTransportation::fpSafetyHelper
G4SafetyHelper * fpSafetyHelper
Definition: MonopoleTransportation.h:158
funct::false
false
Definition: Factorize.h:29
MonopoleTransportation.h
HLT_FULL_cff.verboseLevel
verboseLevel
Definition: HLT_FULL_cff.py:8499
MonopoleTransportation::fTransportEndSpin
G4ThreeVector fTransportEndSpin
Definition: MonopoleTransportation.h:116
MonopoleTransportation::fPreviousSftOrigin
G4ThreeVector fPreviousSftOrigin
Definition: MonopoleTransportation.h:130
CMSFieldManager.h
MonopoleTransportation::fGeometryLimitedStep
G4bool fGeometryLimitedStep
Definition: MonopoleTransportation.h:127
MonopoleTransportation::fMaxEnergyKilled
G4double fMaxEnergyKilled
Definition: MonopoleTransportation.h:152
MonopoleTransportation::fFieldPropagator
G4PropagatorInField * fFieldPropagator
Definition: MonopoleTransportation.h:110
MonopoleTransportation::fieldMgrCMS
CMSFieldManager * fieldMgrCMS
Definition: MonopoleTransportation.h:107
MonopoleTransportation::fCandidateEndGlobalTime
G4double fCandidateEndGlobalTime
Definition: MonopoleTransportation.h:120
sqr
int sqr(const T &t)
Definition: pfalgo_common_ref.h:9
MeV
const double MeV
MonopoleTransportation::fTransportEndPosition
G4ThreeVector fTransportEndPosition
Definition: MonopoleTransportation.h:113
MonopoleTransportation::fTransportEndMomentumDir
G4ThreeVector fTransportEndMomentumDir
Definition: MonopoleTransportation.h:114
MonopoleTransportation::fSumEnergyKilled
G4double fSumEnergyKilled
Definition: MonopoleTransportation.h:151
MonopoleTransportation::fShortStepOptimisation
G4bool fShortStepOptimisation
Definition: MonopoleTransportation.h:156
MonopoleTransportation::fLinearNavigator
G4Navigator * fLinearNavigator
Definition: MonopoleTransportation.h:109
CMSFieldManager::setMonopoleTracking
void setMonopoleTracking(G4bool)
Definition: CMSFieldManager.cc:158
MonopoleTransportation::ResetKilledStatistics
void ResetKilledStatistics(G4int report=1)
Definition: MonopoleTransportation.cc:581
MonopoleTransportation::fParticleChange
G4ParticleChangeForTransport fParticleChange
Definition: MonopoleTransportation.h:134
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
MonopoleTransportation::PostStepDoIt
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData) override
Definition: MonopoleTransportation.cc:444
corrVsCorr.selection
selection
main part
Definition: corrVsCorr.py:100
MonopoleTransportation::DoesGlobalFieldExist
G4bool DoesGlobalFieldExist()
Definition: MonopoleTransportation.h:167
MonopoleTransportation::AtRestDoIt
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &) override
Definition: MonopoleTransportation.cc:577
MonopoleTransportation::AlongStepGetPhysicalInteractionLength
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
Definition: MonopoleTransportation.cc:101
MonopoleTransportation::fCurrentTouchableHandle
G4TouchableHandle fCurrentTouchableHandle
Definition: MonopoleTransportation.h:125
MonopoleTransportation::PostStepGetPhysicalInteractionLength
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond) override
Definition: MonopoleTransportation.cc:435
MonopoleTransportation::fPreviousSafety
G4double fPreviousSafety
Definition: MonopoleTransportation.h:131
MonopoleTransportation::fParticleDef
const Monopole * fParticleDef
Definition: MonopoleTransportation.h:105
Monopole
Definition: Monopole.h:12
MonopoleTransportation::StartTracking
void StartTracking(G4Track *aTrack) override
Definition: MonopoleTransportation.cc:529
Monopole::MagneticCharge
G4double MagneticCharge() const
Definition: Monopole.h:20
MonopoleTransportation::fThreshold_Important_Energy
G4double fThreshold_Important_Energy
Definition: MonopoleTransportation.h:142
GeV
const double GeV
Definition: MathUtil.h:16
MonopoleTransportation::AtRestGetPhysicalInteractionLength
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override
Definition: MonopoleTransportation.cc:573
MonopoleTransportation::fParticleIsLooping
G4bool fParticleIsLooping
Definition: MonopoleTransportation.h:123
MonopoleTransportation::fThresholdTrials
G4int fThresholdTrials
Definition: MonopoleTransportation.h:143
edmIntegrityCheck.report
report
Definition: edmIntegrityCheck.py:349
MonopoleTransportation::fNoLooperTrials
G4int fNoLooperTrials
Definition: MonopoleTransportation.h:149
mag
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Definition: Basic3DVectorLD.h:127
MonopoleTransportation::fThreshold_Warning_Energy
G4double fThreshold_Warning_Energy
Definition: MonopoleTransportation.h:141
MonopoleTransportation::MonopoleTransportation
MonopoleTransportation(const Monopole *p, G4int verbosityLevel=1)
Definition: MonopoleTransportation.cc:33
ecalTB2006H4_GenSimDigiReco_cfg.G4cout
G4cout
Definition: ecalTB2006H4_GenSimDigiReco_cfg.py:303
MonopoleTransportation::endpointDistance
G4double endpointDistance
Definition: MonopoleTransportation.h:137
MonopoleTransportation::fTransportEndKineticEnergy
G4double fTransportEndKineticEnergy
Definition: MonopoleTransportation.h:115
MonopoleTransportation::fMomentumChanged
G4bool fMomentumChanged
Definition: MonopoleTransportation.h:117
MonopoleTransportation::~MonopoleTransportation
~MonopoleTransportation() override
Definition: MonopoleTransportation.cc:83
MonopoleTransportation::fEndGlobalTimeComputed
G4bool fEndGlobalTimeComputed
Definition: MonopoleTransportation.h:119