CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
G4MonopoleTransportation.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4MonopoleTransportation.cc,v 1.1 2010/07/29 23:05:19 sunanda Exp $
27 // GEANT4 tag $Name: CMSSW_5_3_11_patch3 $
28 //
29 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
31 //
32 // This class is a process responsible for the transportation of
33 // magnetic monopoles, ie the geometrical propagation that encounters the
34 // geometrical sub-volumes of the detectors.
35 //
36 // For monopoles, uses a different equation of motion and ignores energy
37 // conservation.
38 //
39 
40 // =======================================================================
41 // Created: 3 May 2010, J. Apostolakis, B. Bozsogi
42 // =======================================================================
43 
44 #include "G4MonopoleTransportation.hh"
45 #include "G4ProductionCutsTable.hh"
46 #include "G4ParticleTable.hh"
47 #include "G4ChordFinder.hh"
48 #include "G4SafetyHelper.hh"
49 #include "G4FieldManagerStore.hh"
50 #include "SimG4Core/Physics/interface/G4Monopole.hh"
51 
52 class G4VSensitiveDetector;
53 
55 //
56 // Constructor
57 
58 G4MonopoleTransportation::G4MonopoleTransportation(const G4Monopole* mpl,
59  sim::FieldBuilder* fieldBuilder,
60  G4int verb)
61  : G4VProcess( G4String("MonopoleTransportation"), fTransportation ),
62  fFieldBuilder(fieldBuilder),
63  fParticleIsLooping( false ),
64  fPreviousSftOrigin (0.,0.,0.),
65  fPreviousSafety ( 0.0 ),
66  fThreshold_Warning_Energy( 100 * MeV ),
67  fThreshold_Important_Energy( 250 * MeV ),
68  fThresholdTrials( 10 ),
69  fUnimportant_Energy( 1 * MeV ),
70  fNoLooperTrials(0),
71  fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
72  fShortStepOptimisation(false), // Old default: true (=fast short steps)
73  fVerboseLevel( verb )
74 {
75  pParticleDef = mpl;
76 
77  //magSetup = G4MonopoleFieldSetup::GetMonopoleFieldSetup();
78 
79  G4TransportationManager* transportMgr ;
80 
81  transportMgr = G4TransportationManager::GetTransportationManager() ;
82 
83  fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
84 
85  // fGlobalFieldMgr = transportMgr->GetFieldManager() ;
86 
87  fFieldPropagator = transportMgr->GetPropagatorInField() ;
88 
89  fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
90 
91  // Cannot determine whether a field exists here,
92  // because it would only work if the field manager has informed
93  // about the detector's field before this transportation process
94  // is constructed.
95  // Instead later the method DoesGlobalFieldExist() is called
96 
97  static G4TouchableHandle nullTouchableHandle; // Points to (G4VTouchable*) 0
98  fCurrentTouchableHandle = nullTouchableHandle;
99 
100  fEndGlobalTimeComputed = false;
101  fCandidateEndGlobalTime = 0;
102 }
103 
105 
106 G4MonopoleTransportation::~G4MonopoleTransportation()
107 {
108  if( (fVerboseLevel > 0) && (fSumEnergyKilled > 0.0 ) ){
109  G4cout << " G4MonopoleTransportation: Statistics for looping particles " << G4endl;
110  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
111  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
112  }
113 }
114 
116 //
117 // Responsibilities:
118 // Find whether the geometry limits the Step, and to what length
119 // Calculate the new value of the safety and return it.
120 // Store the final time, position and momentum.
121 
122 G4double G4MonopoleTransportation::
123 AlongStepGetPhysicalInteractionLength( const G4Track& track,
124  G4double, // previousStepSize
125  G4double currentMinimumStep,
126  G4double& currentSafety,
127  G4GPILSelection* selection )
128 {
129 
130  //magSetup->SetStepperAndChordFinder(1);
131  // change to monopole equation
132  G4FieldManager* fieldMgrX=fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
133  fFieldBuilder->setStepperAndChordFinder (fieldMgrX, 1);
134 
135  G4double geometryStepLength, newSafety ;
136  fParticleIsLooping = false ;
137 
138  // Initial actions moved to StartTrack()
139  // --------------------------------------
140  // Note: in case another process changes touchable handle
141  // it will be necessary to add here (for all steps)
142  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
143 
144  // GPILSelection is set to defaule value of CandidateForSelection
145  // It is a return value
146  //
147  *selection = CandidateForSelection ;
148 
149  // Get initial Energy/Momentum of the track
150  //
151  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
152  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
153  G4ThreeVector startPosition = track.GetPosition() ;
154 
155  // G4double theTime = track.GetGlobalTime() ;
156 
157  // The Step Point safety can be limited by other geometries and/or the
158  // assumptions of any process - it's not always the geometrical safety.
159  // We calculate the starting point's isotropic safety here.
160  //
161  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
162  G4double MagSqShift = OriginShift.mag2() ;
163  if( MagSqShift >= sqr(fPreviousSafety) )
164  {
165  currentSafety = 0.0 ;
166  }
167  else
168  {
169  currentSafety = fPreviousSafety - std::sqrt(MagSqShift) ;
170  }
171 
172  // Is the monopole charged ?
173  //
174  G4double particleMagneticCharge = pParticleDef->MagneticCharge() ;
175  G4double particleElectricCharge = pParticle->GetCharge();
176 
177  fGeometryLimitedStep = false ;
178  // fEndGlobalTimeComputed = false ;
179 
180  // There is no need to locate the current volume. It is Done elsewhere:
181  // On track construction
182  // By the tracking, after all AlongStepDoIts, in "Relocation"
183 
184  // Check whether the particle have an (EM) field force exerting upon it
185  //
186  G4FieldManager* fieldMgr=0;
187  G4bool fieldExertsForce = false ;
188 
189  if( (particleMagneticCharge != 0.0) )
190  {
191  fieldMgr= fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
192  if (fieldMgr != 0) {
193  // Message the field Manager, to configure it for this track
194  fieldMgr->ConfigureForTrack( &track );
195  // Moved here, in order to allow a transition
196  // from a zero-field status (with fieldMgr->(field)0
197  // to a finite field status
198 
199  // If the field manager has no field, there is no field !
200  fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
201  }
202  }
203 
204  // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
205  // << " fieldMgr= " << fieldMgr << G4endl;
206 
207  // Choose the calculation of the transportation: Field or not
208  //
209  if( !fieldExertsForce )
210  {
211  G4double linearStepLength ;
212  if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
213  {
214  // The Step is guaranteed to be taken
215  //
216  geometryStepLength = currentMinimumStep ;
217  fGeometryLimitedStep = false ;
218  }
219  else
220  {
221  // Find whether the straight path intersects a volume
222  //
223  linearStepLength = fLinearNavigator->ComputeStep( startPosition,
224  startMomentumDir,
225  currentMinimumStep,
226  newSafety) ;
227  // Remember last safety origin & value.
228  //
229  fPreviousSftOrigin = startPosition ;
230  fPreviousSafety = newSafety ;
231  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
232 
233  // The safety at the initial point has been re-calculated:
234  //
235  currentSafety = newSafety ;
236 
237  fGeometryLimitedStep= (linearStepLength <= currentMinimumStep);
238  if( fGeometryLimitedStep )
239  {
240  // The geometry limits the Step size (an intersection was found.)
241  geometryStepLength = linearStepLength ;
242  }
243  else
244  {
245  // The full Step is taken.
246  geometryStepLength = currentMinimumStep ;
247  }
248  }
249  endpointDistance = geometryStepLength ;
250 
251  // Calculate final position
252  //
253  fTransportEndPosition = startPosition+geometryStepLength*startMomentumDir ;
254 
255  // Momentum direction, energy and polarisation are unchanged by transport
256  //
257  fTransportEndMomentumDir = startMomentumDir ;
258  fTransportEndKineticEnergy = track.GetKineticEnergy() ;
259  fTransportEndSpin = track.GetPolarization();
260  fParticleIsLooping = false ;
261  fMomentumChanged = false ;
262  fEndGlobalTimeComputed = false ;
263  }
264  else // A field exerts force
265  {
266  // G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
267  G4ThreeVector EndUnitMomentum ;
268  G4double lengthAlongCurve ;
269  G4double restMass = pParticleDef->GetPDGMass() ;
270 
271  fFieldPropagator->SetChargeMomentumMass( particleMagneticCharge, // in Mev/c
272  particleElectricCharge, // in e+ units
273  restMass ) ;
274 
275  // SetChargeMomentumMass is _not_ used here as it would in everywhere else,
276  // it's just a workaround to pass the electric charge as well.
277 
278 
279  G4ThreeVector spin = track.GetPolarization() ;
280  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
281  track.GetMomentumDirection(),
282  0.0,
283  track.GetKineticEnergy(),
284  restMass,
285  track.GetVelocity(),
286  track.GetGlobalTime(), // Lab.
287  track.GetProperTime(), // Part.
288  &spin ) ;
289  if( currentMinimumStep > 0 )
290  {
291  // Do the Transport in the field (non recti-linear)
292  //
293  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
294  currentMinimumStep,
295  currentSafety,
296  track.GetVolume() ) ;
297  fGeometryLimitedStep= lengthAlongCurve < currentMinimumStep;
298  if( fGeometryLimitedStep ) {
299  geometryStepLength = lengthAlongCurve ;
300  } else {
301  geometryStepLength = currentMinimumStep ;
302  }
303  }
304  else
305  {
306  geometryStepLength = lengthAlongCurve= 0.0 ;
307  fGeometryLimitedStep = false ;
308  }
309 
310  // Remember last safety origin & value.
311  //
312  fPreviousSftOrigin = startPosition ;
313  fPreviousSafety = currentSafety ;
314  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
315 
316  // Get the End-Position and End-Momentum (Dir-ection)
317  //
318  fTransportEndPosition = aFieldTrack.GetPosition() ;
319 
320  // Momentum: Magnitude and direction can be changed too now ...
321  //
322  fMomentumChanged = true ;
323  fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
324 
325  fTransportEndKineticEnergy = aFieldTrack.GetKineticEnergy() ;
326 
327 // if( fFieldPropagator->GetCurrentFieldManager()->DoesFieldChangeEnergy() )
328 // {
329 // // If the field can change energy, then the time must be integrated
330 // // - so this should have been updated
331 
332  fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
333  fEndGlobalTimeComputed = true;
334 
335 // // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
336 // // a cleaner way is to have FieldTrack knowing whether time is updated.
337 // }
338 // else
339 // {
340 // // The energy should be unchanged by field transport,
341 // // - so the time changed will be calculated elsewhere
342 // //
343 // fEndGlobalTimeComputed = false;
344 //
345 // // Check that the integration preserved the energy
346 // // - and if not correct this!
347 // G4double startEnergy= track.GetKineticEnergy();
348 // G4double endEnergy= fTransportEndKineticEnergy;
349 //
350 // static G4int no_inexact_steps=0, no_large_ediff;
351 // G4double absEdiff = std::fabs(startEnergy- endEnergy);
352 // if( absEdiff > perMillion * endEnergy )
353 // {
354 // no_inexact_steps++;
355 // // Possible statistics keeping here ...
356 // }
357 // if( fVerboseLevel > 1 )
358 // {
359 // if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
360 // {
361 // static G4int no_warnings= 0, warnModulo=1, moduloFactor= 10;
362 // no_large_ediff ++;
363 // if( (no_large_ediff% warnModulo) == 0 )
364 // {
365 // no_warnings++;
366 // G4cout << "WARNING - G4MonopoleTransportation::AlongStepGetPIL() "
367 // << " Energy change in Step is above 1^-3 relative value. " << G4endl
368 // << " Relative change in 'tracking' step = "
369 // << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
370 // << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV " << G4endl
371 // << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV " << G4endl;
372 // G4cout << " Energy has been corrected -- however, review"
373 // << " field propagation parameters for accuracy." << G4endl;
374 // if( (fVerboseLevel > 2 ) || (no_warnings<4) || (no_large_ediff == warnModulo * moduloFactor) ){
375 // G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
376 // << " which determine fractional error per step for integrated quantities. " << G4endl
377 // << " Note also the influence of the permitted number of integration steps."
378 // << G4endl;
379 // }
380 // G4cerr << "ERROR - G4MonopoleTransportation::AlongStepGetPIL()" << G4endl
381 // << " Bad 'endpoint'. Energy change detected"
382 // << " and corrected. "
383 // << " Has occurred already "
384 // << no_large_ediff << " times." << G4endl;
385 // if( no_large_ediff == warnModulo * moduloFactor )
386 // {
387 // warnModulo *= moduloFactor;
388 // }
389 // }
390 // }
391 // } // end of if (fVerboseLevel)
392 //
393 // // Correct the energy for fields that conserve it
394 // // This - hides the integration error
395 // // - but gives a better physical answer
396 // fTransportEndKineticEnergy= track.GetKineticEnergy();
397 // }
398 
399 
400  fTransportEndSpin = aFieldTrack.GetSpin();
401  fParticleIsLooping = fFieldPropagator->IsParticleLooping() ;
402  endpointDistance = (fTransportEndPosition - startPosition).mag() ;
403  }
404 
405  // If we are asked to go a step length of 0, and we are on a boundary
406  // then a boundary will also limit the step -> we must flag this.
407  //
408  if( currentMinimumStep == 0.0 )
409  {
410  if( currentSafety == 0.0 ) fGeometryLimitedStep = true ;
411  }
412 
413  // Update the safety starting from the end-point,
414  // if it will become negative at the end-point.
415  //
416  if( currentSafety < endpointDistance )
417  {
418  // if( particleMagneticCharge == 0.0 )
419  // G4cout << " Avoiding call to ComputeSafety : charge = 0.0 " << G4endl;
420 
421  if( particleMagneticCharge != 0.0 ) {
422 
423  G4double endSafety =
424  fLinearNavigator->ComputeSafety( fTransportEndPosition) ;
425  currentSafety = endSafety ;
426  fPreviousSftOrigin = fTransportEndPosition ;
427  fPreviousSafety = currentSafety ;
428  fpSafetyHelper->SetCurrentSafety( currentSafety, fTransportEndPosition);
429 
430  // Because the Stepping Manager assumes it is from the start point,
431  // add the StepLength
432  //
433  currentSafety += endpointDistance ;
434 
435 #ifdef G4DEBUG_TRANSPORT
436  G4cout.precision(12) ;
437  G4cout << "***G4MonopoleTransportation::AlongStepGPIL ** " << G4endl ;
438  G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
439  << " and it returned safety= " << endSafety << G4endl ;
440  G4cout << " Adding endpoint distance " << endpointDistance
441  << " to obtain pseudo-safety= " << currentSafety << G4endl ;
442 #endif
443  }
444  }
445 
446  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
447 
448  return geometryStepLength ;
449 }
450 
452 //
453 // Initialize ParticleChange (by setting all its members equal
454 // to corresponding members in G4Track)
455 
456 G4VParticleChange* G4MonopoleTransportation::AlongStepDoIt( const G4Track& track,
457  const G4Step& stepData )
458 {
459  static G4int noCalls=0;
460  static const G4ParticleDefinition* fOpticalPhoton =
461  G4ParticleTable::GetParticleTable()->FindParticle("opticalphoton");
462 
463  noCalls++;
464 
465  fParticleChange.Initialize(track) ;
466 
467  // Code for specific process
468  //
469  fParticleChange.ProposePosition(fTransportEndPosition) ;
470  fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
471  fParticleChange.ProposeEnergy(fTransportEndKineticEnergy) ;
472  fParticleChange.SetMomentumChanged(fMomentumChanged) ;
473 
474  fParticleChange.ProposePolarization(fTransportEndSpin);
475 
476  G4double deltaTime = 0.0 ;
477 
478  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
479  // G4double endTime = fCandidateEndGlobalTime;
480  // G4double delta_time = endTime - startTime;
481 
482  G4double startTime = track.GetGlobalTime() ;
483 
484  if (!fEndGlobalTimeComputed)
485  {
486  // The time was not integrated .. make the best estimate possible
487  //
488  G4double finalVelocity = track.GetVelocity() ;
489  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
490  G4double stepLength = track.GetStepLength() ;
491 
492  deltaTime= 0.0; // in case initialVelocity = 0
493  const G4DynamicParticle* fpDynamicParticle = track.GetDynamicParticle();
494  if (fpDynamicParticle->GetDefinition()== fOpticalPhoton)
495  {
496  // A photon is in the medium of the final point
497  // during the step, so it has the final velocity.
498  deltaTime = stepLength/finalVelocity ;
499  }
500  else if (finalVelocity > 0.0)
501  {
502  G4double meanInverseVelocity ;
503  // deltaTime = stepLength/finalVelocity ;
504  meanInverseVelocity = 0.5
505  * ( 1.0 / initialVelocity + 1.0 / finalVelocity ) ;
506  deltaTime = stepLength * meanInverseVelocity ;
507  }
508  else if( initialVelocity > 0.0 )
509  {
510  deltaTime = stepLength/initialVelocity ;
511  }
512  fCandidateEndGlobalTime = startTime + deltaTime ;
513  }
514  else
515  {
516  deltaTime = fCandidateEndGlobalTime - startTime ;
517  }
518 
519  fParticleChange.ProposeGlobalTime( fCandidateEndGlobalTime ) ;
520 
521  // Now Correct by Lorentz factor to get "proper" deltaTime
522 
523  G4double restMass = track.GetDynamicParticle()->GetMass() ;
524  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
525 
526  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
527  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
528 
529  // If the particle is caught looping or is stuck (in very difficult
530  // boundaries) in a magnetic field (doing many steps)
531  // THEN this kills it ...
532  //
533  if ( fParticleIsLooping )
534  {
535  G4double endEnergy= fTransportEndKineticEnergy;
536 
537  if( (endEnergy < fThreshold_Important_Energy)
538  || (fNoLooperTrials >= fThresholdTrials ) ){
539  // Kill the looping particle
540  //
541  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
542 
543  // 'Bare' statistics
544  fSumEnergyKilled += endEnergy;
545  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
546 
547 #ifdef G4VERBOSE
548  if( (fVerboseLevel > 1) ||
549  ( endEnergy > fThreshold_Warning_Energy ) ) {
550  G4cout << " G4MonopoleTransportation is killing track that is looping or stuck "
551  << G4endl
552  << " This track has " << track.GetKineticEnergy() / MeV
553  << " MeV energy." << G4endl;
554  G4cout << " Number of trials = " << fNoLooperTrials
555  << " No of calls to AlongStepDoIt = " << noCalls
556  << G4endl;
557  }
558 #endif
559  fNoLooperTrials=0;
560  }
561  else{
562  fNoLooperTrials ++;
563 #ifdef G4VERBOSE
564  if( (fVerboseLevel > 2) ){
565  G4cout << " G4MonopoleTransportation::AlongStepDoIt(): Particle looping - "
566  << " Number of trials = " << fNoLooperTrials
567  << " No of calls to = " << noCalls
568  << G4endl;
569  }
570 #endif
571  }
572  }else{
573  fNoLooperTrials=0;
574  }
575 
576  // Another (sometimes better way) is to use a user-limit maximum Step size
577  // to alleviate this problem ..
578 
579  // Introduce smooth curved trajectories to particle-change
580  //
581  fParticleChange.SetPointerToVectorOfAuxiliaryPoints
582  (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
583 
584  return &fParticleChange ;
585 }
586 
588 //
589 // This ensures that the PostStep action is always called,
590 // so that it can do the relocation if it is needed.
591 //
592 
593 G4double G4MonopoleTransportation::
594 PostStepGetPhysicalInteractionLength( const G4Track& track,
595  G4double, // previousStepSize
596  G4ForceCondition* pForceCond )
597 {
598  //magSetup->SetStepperAndChordFinder(0);
599  // change back to usual equation
600  G4FieldManager* fieldMgr=fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
601  fFieldBuilder->setStepperAndChordFinder (fieldMgr, 0);
602 
603  *pForceCond = Forced ;
604  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
605 }
606 
608 //
609 
610 G4VParticleChange* G4MonopoleTransportation::PostStepDoIt( const G4Track& track,
611  const G4Step& )
612 {
613  G4TouchableHandle retCurrentTouchable ; // The one to return
614 
615  // Initialize ParticleChange (by setting all its members equal
616  // to corresponding members in G4Track)
617  // fParticleChange.Initialize(track) ; // To initialise TouchableChange
618 
619  fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
620 
621  // If the Step was determined by the volume boundary,
622  // logically relocate the particle
623 
624  if(fGeometryLimitedStep)
625  {
626  // fCurrentTouchable will now become the previous touchable,
627  // and what was the previous will be freed.
628  // (Needed because the preStepPoint can point to the previous touchable)
629 
630  fLinearNavigator->SetGeometricallyLimitedStep() ;
631  fLinearNavigator->
632  LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
633  track.GetMomentumDirection(),
634  fCurrentTouchableHandle,
635  true ) ;
636  // Check whether the particle is out of the world volume
637  // If so it has exited and must be killed.
638  //
639  if( fCurrentTouchableHandle->GetVolume() == 0 )
640  {
641  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
642  }
643  retCurrentTouchable = fCurrentTouchableHandle ;
644  fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
645  }
646  else // fGeometryLimitedStep is false
647  {
648  // This serves only to move the Navigator's location
649  //
650  fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
651 
652  // The value of the track's current Touchable is retained.
653  // (and it must be correct because we must use it below to
654  // overwrite the (unset) one in particle change)
655  // It must be fCurrentTouchable too ??
656  //
657  fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
658  retCurrentTouchable = track.GetTouchableHandle() ;
659  } // endif ( fGeometryLimitedStep )
660 
661  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
662  const G4Material* pNewMaterial = 0 ;
663  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
664 
665  if( pNewVol != 0 )
666  {
667  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
668  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
669  }
670 
671  // ( <const_cast> pNewMaterial ) ;
672  // ( <const_cast> pNewSensitiveDetector) ;
673 
674  fParticleChange.SetMaterialInTouchable( (G4Material *) pNewMaterial ) ;
675  fParticleChange.SetSensitiveDetectorInTouchable( (G4VSensitiveDetector *) pNewSensitiveDetector ) ;
676 
677  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
678  if( pNewVol != 0 )
679  {
680  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
681  }
682 
683  if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
684  {
685  // for parametrized volume
686  //
687  pNewMaterialCutsCouple =
688  G4ProductionCutsTable::GetProductionCutsTable()
689  ->GetMaterialCutsCouple(pNewMaterial,
690  pNewMaterialCutsCouple->GetProductionCuts());
691  }
692  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
693 
694  // temporarily until Get/Set Material of ParticleChange,
695  // and StepPoint can be made const.
696  // Set the touchable in ParticleChange
697  // this must always be done because the particle change always
698  // uses this value to overwrite the current touchable pointer.
699  //
700  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
701 
702  return &fParticleChange ;
703 }
704 
705 // New method takes over the responsibility to reset the state of G4MonopoleTransportation
706 // object at the start of a new track or the resumption of a suspended track.
707 
708 void
709 G4MonopoleTransportation::StartTracking(G4Track* aTrack)
710 {
711  G4VProcess::StartTracking(aTrack);
712 
713 // The actions here are those that were taken in AlongStepGPIL
714 // when track.GetCurrentStepNumber()==1
715 
716  // reset safety value and center
717  //
718  fPreviousSafety = 0.0 ;
719  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
720 
721  // reset looping counter -- for motion in field
722  fNoLooperTrials= 0;
723  // Must clear this state .. else it depends on last track's value
724  // --> a better solution would set this from state of suspended track TODO ?
725  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
726 
727  // ChordFinder reset internal state
728  //
729  if( DoesGlobalFieldExist() ) {
730  fFieldPropagator->ClearPropagatorState();
731  // Resets all state of field propagator class (ONLY)
732  // including safety values (in case of overlaps and to wipe for first track).
733 
734  // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
735  // if( chordF ) chordF->ResetStepEstimate();
736  }
737 
738  // Make sure to clear the chord finders of all fields (ie managers)
739  static G4FieldManagerStore* fieldMgrStore= G4FieldManagerStore::GetInstance();
740  fieldMgrStore->ClearAllChordFindersState();
741 
742  // Update the current touchable handle (from the track's)
743  //
744  fCurrentTouchableHandle = aTrack->GetTouchableHandle();
745 }
746 
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
selection
main part
Definition: corrVsCorr.py:98
T sqrt(T t)
Definition: SSEVec.h:46
Square< F >::type sqr(const F &f)
Definition: Square.h:13
float spin(float ph)