Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ErrorPropagator Class Reference

#include <G4ErrorPropagator.hh>

Public Member Functions

 G4ErrorPropagator ()
 
 ~G4ErrorPropagator ()
 
G4TrackInitG4Track (G4ErrorTrajState &initialTS)
 
G4int Propagate (G4ErrorTrajState *currentTS, const G4ErrorTarget *target, G4ErrorMode mode=G4ErrorMode_PropForwards)
 
G4int PropagateOneStep (G4ErrorTrajState *currentTS)
 
G4int MakeOneStep (G4ErrorFreeTrajState *currentTS_FREE)
 
G4ErrorFreeTrajStateInitFreeTrajState (G4ErrorTrajState *currentTS)
 
void GetFinalTrajState (G4ErrorTrajState *currentTS, G4ErrorFreeTrajState *currentTS_FREE, const G4ErrorTarget *target)
 
void InvokePreUserTrackingAction (G4Track *fpTrack)
 
void InvokePostUserTrackingAction (G4Track *fpTrack)
 
G4bool CheckIfLastStep (G4Track *aTrack)
 
const G4ErrorTrajStateGetInitialTrajState () const
 
G4double GetStepLength () const
 
void SetStepLength (const G4double sl)
 
void SetStepN (const G4int sn)
 

Detailed Description

Definition at line 53 of file G4ErrorPropagator.hh.

Constructor & Destructor Documentation

◆ G4ErrorPropagator()

G4ErrorPropagator::G4ErrorPropagator ( )

Definition at line 56 of file G4ErrorPropagator.cc.

57 : theStepLength(0.), theInitialTrajState(0), theStepN(0), theG4Track(0)
58{
60#ifdef G4EVERBOSE
61 if(verbose >= 5) { G4cout << "G4ErrorPropagator " << this << G4endl; }
62#endif
63
64 fpSteppingManager = G4EventManager::GetEventManager()
66 thePropIsInitialized = false;
67}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4EventManager * GetEventManager()
G4TrackingManager * GetTrackingManager() const
G4SteppingManager * GetSteppingManager() const

◆ ~G4ErrorPropagator()

G4ErrorPropagator::~G4ErrorPropagator ( )
inline

Definition at line 58 of file G4ErrorPropagator.hh.

58{}

Member Function Documentation

◆ CheckIfLastStep()

G4bool G4ErrorPropagator::CheckIfLastStep ( G4Track aTrack)

Definition at line 526 of file G4ErrorPropagator.cc.

527{
528 G4bool lastG4eStep = false;
529 G4ErrorPropagatorData* g4edata =
531
532#ifdef G4EVERBOSE
533 if( verbose >= 4 )
534 {
535 G4cout << " G4ErrorPropagator::CheckIfLastStep G4ErrorState= "
536 << G4int(g4edata->GetState()) << G4endl;
537 }
538#endif
539
540 //----- Check if this is the last step: track has reached the target
541 // or the end of world
542 //
544 {
545 lastG4eStep = true;
546#ifdef G4EVERBOSE
547 if(verbose >= 5 )
548 {
549 G4cout << " G4ErrorPropagator::CheckIfLastStep " << lastG4eStep
550 << " " << G4int(g4edata->GetState()) << G4endl;
551 }
552#endif
553 }
554 else if( aTrack->GetNextVolume() == 0 )
555 {
556 //----- If particle is out of world, without finding the G4ErrorTarget,
557 // give a n error/warning
558 //
559 lastG4eStep = true;
560 if( verbose >= 1 )
561 {
562 std::ostringstream message;
563 message << "Track extrapolated until end of World" << G4endl
564 << "without finding the defined target.";
565 G4Exception("G4ErrorPropagator::CheckIfLastStep()",
566 "GEANT4e-Notification", JustWarning, message);
567 }
568 } //----- not last step from G4e, but track is stopped (energy exhausted)
569 else if( aTrack->GetTrackStatus() == fStopAndKill )
570 {
571 if( verbose >= 1 )
572 {
573 std::ostringstream message;
574 message << "Track extrapolated until energy is exhausted" << G4endl
575 << "without finding the defined target.";
576 G4Exception("G4ErrorPropagator::CheckIfLastStep()",
577 "GEANT4e-Notification", JustWarning, message);
578 }
579 lastG4eStep = 1;
580 }
581
582#ifdef G4EVERBOSE
583 if( verbose >= 5 )
584 {
585 G4cout << " return CheckIfLastStep " << lastG4eStep << G4endl;
586 }
587#endif
588
589 return lastG4eStep;
590}
@ G4ErrorState_StoppedAtTarget
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ fStopAndKill
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4ErrorState GetState() const
static G4ErrorPropagatorData * GetErrorPropagatorData()
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetNextVolume() const

◆ GetFinalTrajState()

void G4ErrorPropagator::GetFinalTrajState ( G4ErrorTrajState currentTS,
G4ErrorFreeTrajState currentTS_FREE,
const G4ErrorTarget target 
)

Definition at line 480 of file G4ErrorPropagator.cc.

483{
484 G4ErrorPropagatorData* g4edata =
486
487#ifdef G4EVERBOSE
488 if(verbose >= 1 )
489 {
490 G4cout << " G4ErrorPropagator::Propagate: final state "
491 << G4int(g4edata->GetState()) << " TSType "
492 << currentTS->GetTSType() << G4endl;
493 }
494#endif
495
496 if( (currentTS->GetTSType() == G4eTS_FREE) ||
497 (g4edata->GetState() != G4ErrorState_StoppedAtTarget) )
498 {
499 currentTS = currentTS_FREE;
500 }
501 else if( currentTS->GetTSType() == G4eTS_OS )
502 {
503 if( target->GetType() == G4ErrorTarget_TrkL )
504 {
505 G4Exception("G4ErrorPropagator:GetFinalTrajState()",
506 "InvalidSetup", FatalException,
507 "Using a G4ErrorSurfaceTrajState with wrong target");
508 }
509 const G4ErrorTanPlaneTarget* targetWTP =
510 static_cast<const G4ErrorTanPlaneTarget*>(target);
511 *currentTS = G4ErrorSurfaceTrajState(
512 *(static_cast<G4ErrorFreeTrajState*>(currentTS_FREE)),
513 targetWTP->GetTangentPlane( currentTS_FREE->GetPosition() ) );
514#ifdef G4EVERBOSE
515 if(verbose >= 1 )
516 {
517 G4cout << currentTS << " returning tssd " << *currentTS << G4endl;
518 }
519#endif
520 delete currentTS_FREE;
521 }
522}
@ G4ErrorTarget_TrkL
@ G4eTS_FREE
@ G4eTS_OS
@ FatalException
virtual G4Plane3D GetTangentPlane(const G4ThreeVector &point) const =0
G4ErrorTargetType GetType() const
virtual G4eTSType GetTSType() const
G4Point3D GetPosition() const

Referenced by Propagate(), and PropagateOneStep().

◆ GetInitialTrajState()

const G4ErrorTrajState * G4ErrorPropagator::GetInitialTrajState ( ) const
inline

Definition at line 98 of file G4ErrorPropagator.hh.

99 { return theInitialTrajState; }

◆ GetStepLength()

G4double G4ErrorPropagator::GetStepLength ( ) const
inline

Definition at line 101 of file G4ErrorPropagator.hh.

102 { return theStepLength; }

◆ InitFreeTrajState()

G4ErrorFreeTrajState * G4ErrorPropagator::InitFreeTrajState ( G4ErrorTrajState currentTS)

Definition at line 452 of file G4ErrorPropagator.cc.

453{
454 G4ErrorFreeTrajState* currentTS_FREE = 0;
455
456 //----- Transform the TrajState to Free coordinates if it is OnSurface
457 //
458 if( currentTS->GetTSType() == G4eTS_OS )
459 {
461 static_cast<G4ErrorSurfaceTrajState*>(currentTS);
462 currentTS_FREE = new G4ErrorFreeTrajState( *tssd );
463 }
464 else if( currentTS->GetTSType() == G4eTS_FREE )
465 {
466 currentTS_FREE = static_cast<G4ErrorFreeTrajState*>(currentTS);
467 }
468 else
469 {
470 std::ostringstream message;
471 message << "Wrong trajectory state: " << currentTS->GetTSType();
472 G4Exception("G4ErrorPropagator::InitFreeTrajState()", "InvalidState",
473 FatalException, message);
474 }
475 return currentTS_FREE;
476}

Referenced by Propagate(), and PropagateOneStep().

◆ InitG4Track()

G4Track * G4ErrorPropagator::InitG4Track ( G4ErrorTrajState initialTS)

Definition at line 254 of file G4ErrorPropagator.cc.

255{
256 if( verbose >= 5 ) { G4cout << "InitG4Track " << G4endl; }
257
258 //----- Create Particle
259 //
260 const G4String partType = initialTS.GetParticleType();
262 G4ParticleDefinition* particle = particleTable->FindParticle(partType);
263 if( particle == 0)
264 {
265 std::ostringstream message;
266 message << "Particle type not defined: " << partType;
267 G4Exception("G4ErrorPropagator::InitG4Track()", "InvalidSetup",
268 FatalException, message);
269 }
270
271 G4DynamicParticle* DP =
272 new G4DynamicParticle(particle,initialTS.GetMomentum());
273
274 DP->SetPolarization(0.,0.,0.);
275
276 // Set Charge
277 //
278 if( particle->GetPDGCharge() < 0 )
279 {
280 DP->SetCharge(-eplus);
281 }
282 else
283 {
284 DP->SetCharge(eplus);
285 }
286
287 //----- Create Track
288 //
289 theG4Track = new G4Track(DP, 0., initialTS.GetPosition() );
290 theG4Track->SetParentID(0);
291
292#ifdef G4EVERBOSE
293 if(verbose >= 3)
294 {
295 G4cout << " G4ErrorPropagator new track of energy: "
296 << theG4Track->GetKineticEnergy() << G4endl;
297 }
298#endif
299
300 //---- Reproduce G4TrackingManager::ProcessOneTrack initialization
301 InvokePreUserTrackingAction( theG4Track );
302
303 if( fpSteppingManager == 0 )
304 {
305 G4Exception("G4ErrorPropagator::InitG4Track()", "InvalidSetup",
306 FatalException, "G4SteppingManager not initialized yet!");
307 }
308 else
309 {
310 fpSteppingManager->SetInitialStep(theG4Track);
311 }
312
313 // Give SteppingManger the maximum number of processes
314 //
315 fpSteppingManager->GetProcessNumber();
316
317 // Give track the pointer to the Step
318 //
319 theG4Track->SetStep(fpSteppingManager->GetStep());
320
321 // Inform beginning of tracking to physics processes
322 //
323 theG4Track->GetDefinition()->GetProcessManager()->StartTracking(theG4Track);
324
325 initialTS.SetG4Track( theG4Track );
326
327 return theG4Track;
328}
void SetPolarization(const G4ThreeVector &)
void SetCharge(G4double charge)
void InvokePreUserTrackingAction(G4Track *fpTrack)
void SetG4Track(G4Track *trk)
G4Vector3D GetMomentum() const
const G4String & GetParticleType() const
G4ProcessManager * GetProcessManager() const
G4double GetPDGCharge() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
void StartTracking(G4Track *aTrack=nullptr)
void SetInitialStep(G4Track *valueTrack)
G4Step * GetStep() const
void SetStep(const G4Step *aValue)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
void SetParentID(const G4int aValue)

Referenced by Propagate(), and PropagateOneStep().

◆ InvokePostUserTrackingAction()

void G4ErrorPropagator::InvokePostUserTrackingAction ( G4Track fpTrack)

Definition at line 607 of file G4ErrorPropagator.cc.

608{
609 const G4UserTrackingAction* fpUserTrackingAction =
611 if( fpUserTrackingAction != 0 )
612 {
613 const_cast<G4UserTrackingAction*>(fpUserTrackingAction)
614 ->PostUserTrackingAction((fpTrack) );
615 }
616}
G4UserTrackingAction * GetUserTrackingAction()

Referenced by Propagate().

◆ InvokePreUserTrackingAction()

void G4ErrorPropagator::InvokePreUserTrackingAction ( G4Track fpTrack)

Definition at line 594 of file G4ErrorPropagator.cc.

595{
596 const G4UserTrackingAction* fpUserTrackingAction =
598 if( fpUserTrackingAction != 0 )
599 {
600 const_cast<G4UserTrackingAction*>(fpUserTrackingAction)
601 ->PreUserTrackingAction((fpTrack) );
602 }
603}

Referenced by InitG4Track().

◆ MakeOneStep()

G4int G4ErrorPropagator::MakeOneStep ( G4ErrorFreeTrajState currentTS_FREE)

Definition at line 358 of file G4ErrorPropagator.cc.

359{
360 G4ErrorPropagatorData* g4edata =
362 G4int ierr = 0;
363
364 //---------- Track one step
365#ifdef G4EVERBOSE
366 if(verbose >= 2 )
367 {
368 G4cout << G4endl
369 << "@@@@@@@@@@@@@@@@@@@@@@@@@ NEW STEP " << G4endl;
370 }
371#endif
372
373 theG4Track->IncrementCurrentStepNumber();
374
375 fpSteppingManager->Stepping();
376
377 //---------- Check if Target has been reached (and then set G4ErrorState)
378
379 // G4ErrorPropagationNavigator limits the step if target is closer than
380 // boundary (but the winner process is always "Transportation": then
381 // error propagator will stop the track
382
383 if( theG4Track->GetStep()->GetPostStepPoint()
384 ->GetProcessDefinedStep()->GetProcessName() == "Transportation" )
385 {
386 if( g4edata->GetState()
388 { // target or step length reached
389
390#ifdef G4EVERBOSE
391 if(verbose >= 5 )
392 {
393 G4cout << " transportation determined by geant4e " << G4endl;
394 }
395#endif
397 }
398 else if( g4edata->GetTarget()->GetType() == G4ErrorTarget_GeomVolume )
399 {
401 (G4ErrorGeomVolumeTarget*)(g4edata->GetTarget());
402 if( static_cast<G4ErrorGeomVolumeTarget*>( target )
403 ->TargetReached( theG4Track->GetStep() ) )
404 {
406 }
407 }
408 }
409 else if( theG4Track->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()
410 ->GetProcessName() == "G4ErrorTrackLengthTarget" )
411 {
413 }
414
415 //---------- Propagate error
416
417#ifdef G4EVERBOSE
418 if(verbose >= 2 )
419 {
420 G4cout << " propagating error " << *currentTS_FREE << G4endl;
421 }
422#endif
423 const G4Track* cTrack = const_cast<G4Track*>(theG4Track);
424 ierr = currentTS_FREE->PropagateError( cTrack );
425
426#ifdef G4EVERBOSE
427 if(verbose >= 3 )
428 {
429 G4cout << " PropagateError returns " << ierr << G4endl;
430 }
431#endif
432
433 currentTS_FREE->Update( cTrack );
434
435 theStepLength += theG4Track->GetStepLength();
436
437 if(ierr != 0 )
438 {
439 std::ostringstream message;
440 message << "Error returned: " << ierr;
441 G4Exception("G4ErrorPropagator::MakeOneStep()",
442 "GEANT4e-Notification", JustWarning, message,
443 "Geant4 tracking will be stopped !");
444 }
445
446 return ierr;
447}
@ G4ErrorState_TargetCloserThanBoundary
@ G4ErrorTarget_GeomVolume
virtual G4int Update(const G4Track *aTrack)
virtual G4int PropagateError(const G4Track *aTrack)
virtual G4bool TargetReached(const G4Step *aStep)
const G4ErrorTarget * GetTarget(G4bool mustExist=false) const
void SetState(G4ErrorState sta)
const G4VProcess * GetProcessDefinedStep() const
G4StepPoint * GetPostStepPoint() const
G4StepStatus Stepping()
G4double GetStepLength() const
void IncrementCurrentStepNumber()
const G4Step * GetStep() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

Referenced by PropagateOneStep().

◆ Propagate()

G4int G4ErrorPropagator::Propagate ( G4ErrorTrajState currentTS,
const G4ErrorTarget target,
G4ErrorMode  mode = G4ErrorMode_PropForwards 
)

Definition at line 71 of file G4ErrorPropagator.cc.

73{
74 // to start ierror is set to 1 (= OK)
75 //
76 G4int ierr = 1;
77
78 G4ErrorPropagatorData* g4edata =
80
81 //----- Do not propagate zero or too low energy particles
82 //
83 if( currentTS->GetMomentum().mag() < 1.E-9*MeV )
84 {
85 std::ostringstream message;
86 message << "Energy too low to be propagated: "
87 << G4BestUnit(currentTS->GetMomentum().mag(),"Energy");
88 G4Exception("G4ErrorPropagator::Propagate()", "GEANT4e-Notification",
89 JustWarning, message);
90 return -3;
91 }
92
93 g4edata->SetMode( mode );
94
95#ifdef G4EVERBOSE
96 if( verbose >= 1 )
97 {
98 G4cout << " =====> starting GEANT4E tracking for "
99 << currentTS->GetParticleType()
100 << " Forwards= " << g4edata->GetMode() << G4endl;
101 }
102 if(verbose >= 1 )
103 {
104 G4cout << G4endl << "@@@@@@@@@@@@@@@@@@@@@@@@@ NEW STEP " << G4endl;
105 }
106
107 if( verbose >= 3 )
108 {
109 G4cout << " G4ErrorPropagator::Propagate initialTS ";
110 G4cout << *currentTS << G4endl;
111 target->Dump(G4String(" to target "));
112 }
113#endif
114
115 g4edata->SetTarget( target );
116
117 //----- Create a track
118 //
119 if( theG4Track != 0 ) { delete theG4Track; }
120 theG4Track = InitG4Track( *currentTS );
121
122 //----- Create a G4ErrorFreeTrajState
123 //
124 G4ErrorFreeTrajState* currentTS_FREE = InitFreeTrajState( currentTS );
125
126 //----- Track the particle
127 //
128 ierr = MakeSteps( currentTS_FREE );
129
130 //------ Tracking ended, check if target has been reached
131 // if target not found
132 //
133 if( g4edata->GetState() != G4ErrorState_StoppedAtTarget )
134 {
135 if( theG4Track->GetKineticEnergy() > 0. )
136 {
137 ierr = -ierr - 10;
138 }
139 else
140 {
141 ierr = -ierr - 20;
142 }
143 *currentTS = *currentTS_FREE;
144 if(verbose >= 0 )
145 {
146 std::ostringstream message;
147 message << "Particle does not reach target: " << *currentTS;
148 G4Exception("G4ErrorPropagator::Propagate()", "GEANT4e-Notification",
149 JustWarning, message);
150 }
151 }
152 else
153 {
154 GetFinalTrajState( currentTS, currentTS_FREE, target );
155 }
156
157#ifdef G4EVERBOSE
158 if( verbose >= 1 )
159 {
160 G4cout << " G4ErrorPropagator: propagation ended " << G4endl;
161 }
162 if( verbose >= 2 )
163 {
164 G4cout << " Current TrajState " << currentTS << G4endl;
165 }
166#endif
167
168 // Inform end of tracking to physics processes
169 //
170 theG4Track->GetDefinition()->GetProcessManager()->EndTracking();
171
172 InvokePostUserTrackingAction( theG4Track );
173
174 // delete currentTS_FREE;
175
176 return ierr;
177}
#define G4BestUnit(a, b)
void SetMode(G4ErrorMode mode)
void SetTarget(const G4ErrorTarget *target)
G4ErrorMode GetMode() const
G4Track * InitG4Track(G4ErrorTrajState &initialTS)
void InvokePostUserTrackingAction(G4Track *fpTrack)
void GetFinalTrajState(G4ErrorTrajState *currentTS, G4ErrorFreeTrajState *currentTS_FREE, const G4ErrorTarget *target)
G4ErrorFreeTrajState * InitFreeTrajState(G4ErrorTrajState *currentTS)
virtual void Dump(const G4String &msg) const =0

Referenced by G4ErrorPropagatorManager::Propagate().

◆ PropagateOneStep()

G4int G4ErrorPropagator::PropagateOneStep ( G4ErrorTrajState currentTS)

Definition at line 181 of file G4ErrorPropagator.cc.

182{
183 G4ErrorPropagatorData* g4edata =
185
186 if ( (g4edata->GetState() == G4ErrorState_PreInit)
189 {
190 std::ostringstream message;
191 message << "Called before initialization is done for this track!";
192 G4Exception("G4ErrorPropagator::PropagateOneStep()",
193 "InvalidCall", FatalException, message,
194 "Please call G4ErrorPropagatorManager::InitGeant4e().");
195 }
196
197 // to start ierror is set to 0 (= OK)
198 //
199 G4int ierr = 0;
200
201 //--- Do not propagate zero or too low energy particles
202 //
203 if( currentTS->GetMomentum().mag() < 1.E-9*MeV )
204 {
205 std::ostringstream message;
206 message << "Energy too low to be propagated: "
207 << G4BestUnit(currentTS->GetMomentum().mag(),"Energy");
208 G4Exception("G4ErrorPropagator::PropagateOneStep()",
209 "GEANT4e-Notification", JustWarning, message);
210 return -3;
211 }
212
213#ifdef G4EVERBOSE
214 if( verbose >= 1 )
215 {
216 G4cout << " =====> starting GEANT4E tracking for "
217 << currentTS->GetParticleType()
218 << " Forwards= " << g4edata->GetMode() << G4endl;
219 }
220 if( verbose >= 3 )
221 {
222 G4cout << " G4ErrorPropagator::Propagate initialTS ";
223 G4cout << *currentTS << G4endl;
224 }
225#endif
226
227 //----- If it is the first step, create a track
228 //
229 if( theStepN == 0 )
230 {
231 if( theG4Track != 0 ) { delete theG4Track; }
232 theG4Track = InitG4Track( *currentTS );
233 }
234 // set to 0 by the initialization in G4ErrorPropagatorManager
235 theStepN++;
236
237 //----- Create a G4ErrorFreeTrajState
238 //
239 G4ErrorFreeTrajState* currentTS_FREE = InitFreeTrajState( currentTS );
240
241 //----- Track the particle one step
242 //
243 ierr = MakeOneStep( currentTS_FREE );
244
245 //----- Get the state on target
246 //
247 GetFinalTrajState( currentTS, currentTS_FREE, g4edata->GetTarget() );
248
249 return ierr;
250}
@ G4State_GeomClosed
@ G4ErrorState_PreInit
G4int MakeOneStep(G4ErrorFreeTrajState *currentTS_FREE)
const G4ApplicationState & GetCurrentState() const
static G4StateManager * GetStateManager()

Referenced by G4ErrorPropagatorManager::PropagateOneStep().

◆ SetStepLength()

void G4ErrorPropagator::SetStepLength ( const G4double  sl)
inline

Definition at line 104 of file G4ErrorPropagator.hh.

105 { theStepLength = sl; }

◆ SetStepN()

void G4ErrorPropagator::SetStepN ( const G4int  sn)
inline

Definition at line 107 of file G4ErrorPropagator.hh.

108 { theStepN = sn; }

Referenced by G4ErrorPropagatorManager::InitTrackPropagation().


The documentation for this class was generated from the following files: