Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ErrorPropagator.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$
27//
28// ------------------------------------------------------------
29// GEANT 4 class implementation file
30// ------------------------------------------------------------
31//
32
33#include "G4ErrorPropagator.hh"
39
40#include "G4SystemOfUnits.hh"
41#include "G4DynamicParticle.hh"
42#include "G4Track.hh"
43#include "G4SteppingManager.hh"
44#include "G4EventManager.hh"
45#include "G4TrackingManager.hh"
46#include "G4ParticleTable.hh"
47#include "G4StateManager.hh"
48
49#include "G4VPhysicalVolume.hh"
51#include "G4UnitsTable.hh"
52
53#include <vector>
54
55
56//---------------------------------------------------------------------------
58 : theStepLength(0.), theInitialTrajState(0),
59 theFinalTrajState(0), theStepN(0), theG4Track(0)
60{
62#ifdef G4EVERBOSE
63 if(verbose >= 5) { G4cout << "G4ErrorPropagator " << this << G4endl; }
64#endif
65
66 fpSteppingManager = G4EventManager::GetEventManager()
68 thePropIsInitialized = false;
69}
70
71
72//-----------------------------------------------------------------------
74 const G4ErrorTarget* target, G4ErrorMode mode )
75{
76 // to start ierror is set to 1 (= OK)
77 //
78 G4int ierr = 1;
79
80 G4ErrorPropagatorData* g4edata =
82
83 //----- Do not propagate zero or too low energy particles
84 //
85 if( currentTS->GetMomentum().mag() < 1.E-9*MeV )
86 {
87 G4cerr << "ERROR - G4ErrorPropagator::Propagate()" << G4endl
88 << " Energy too low to be propagated: "
89 << G4BestUnit(currentTS->GetMomentum().mag(),"Energy") << G4endl;
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 G4cerr << "ERROR - G4ErrorPropagator::Propagate()" << G4endl
147 << " Particle does not reach target: " << *currentTS
148 << G4endl;
149 }
150 }
151 else
152 {
153 GetFinalTrajState( currentTS, currentTS_FREE, target );
154 }
155
156#ifdef G4EVERBOSE
157 if( verbose >= 1 )
158 {
159 G4cout << " G4ErrorPropagator: propagation ended " << G4endl;
160 }
161 if( verbose >= 2 )
162 {
163 G4cout << " Current TrajState " << currentTS << G4endl;
164 }
165#endif
166
167 // Inform end of tracking to physics processes
168 //
169 theG4Track->GetDefinition()->GetProcessManager()->EndTracking();
170
171 InvokePostUserTrackingAction( theG4Track );
172
173 // delete currentTS_FREE;
174
175 return ierr;
176}
177
178
179//-----------------------------------------------------------------------
181{
182 G4ErrorPropagatorData* g4edata =
184
185 if ( (g4edata->GetState() == G4ErrorState_PreInit)
188 {
189 G4cout << "ERROR - G4ErrorPropagator::PropagateOneStep()" << G4endl
190 << " Called before initialization is done for this track."
191 << G4endl
192 << " Please call G4ErrorPropagatorManager::InitGeant4e()."
193 << G4endl;
194 G4Exception("G4ErrorPropagator::PropagateOneStep()",
195 "InvalidCall", FatalException,
196 "Called before initialization is done for this track!");
197 }
198
199 // to start ierror is set to 0 (= OK)
200 //
201 G4int ierr = 0;
202
203 //--- Do not propagate zero or too low energy particles
204 //
205 if( currentTS->GetMomentum().mag() < 1.E-9*MeV )
206 {
207 G4cerr << "ERROR - G4ErrorPropagator::PropagateOneStep()" << G4endl
208 << " Energy too low to be propagated: "
209 << G4BestUnit(currentTS->GetMomentum().mag(),"Energy") << G4endl;
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 ) { theG4Track = InitG4Track( *currentTS ); }
230 // set to 0 by the initialization in G4ErrorPropagatorManager
231 theStepN++;
232
233 //----- Create a G4ErrorFreeTrajState
234 //
235 G4ErrorFreeTrajState* currentTS_FREE = InitFreeTrajState( currentTS );
236
237 //----- Track the particle one step
238 //
239 ierr = MakeOneStep( currentTS_FREE );
240
241 //----- Get the state on target
242 //
243 GetFinalTrajState( currentTS, currentTS_FREE, g4edata->GetTarget() );
244
245 return ierr;
246}
247
248
249//---------------------------------------------------------------------------
251{
252 if( verbose >= 5 ) { G4cout << "InitG4Track " << G4endl; }
253
254 //----- Create Particle
255 //
256 const G4String partType = initialTS.GetParticleType();
258 G4ParticleDefinition* particle = particleTable->FindParticle(partType);
259 if( particle == 0)
260 {
261 G4cerr << "ERROR - G4ErrorPropagator::InitG4Track()" << G4endl
262 << " Particle type not defined " + partType << G4endl;
263 G4Exception( "G4ErrorPropagator::InitG4Track()", "InvalidSetup",
264 FatalException, "Particle type not defined !" );
265 }
266
267 G4DynamicParticle* DP =
268 new G4DynamicParticle(particle,initialTS.GetMomentum());
269
270 DP->SetPolarization(0.,0.,0.);
271
272 // Set Charge
273 //
274 if( particle->GetPDGCharge() < 0 )
275 {
276 DP->SetCharge(-eplus);
277 }
278 else
279 {
280 DP->SetCharge(eplus);
281 }
282
283 //----- Create Track
284 //
285 theG4Track = new G4Track(DP, 0., initialTS.GetPosition() );
286 theG4Track->SetParentID(0);
287
288#ifdef G4EVERBOSE
289 if(verbose >= 3)
290 {
291 G4cout << " G4ErrorPropagator new track of energy: "
292 << theG4Track->GetKineticEnergy() << G4endl;
293 }
294#endif
295
296 //---- Reproduce G4TrackingManager::ProcessOneTrack initialization
297 InvokePreUserTrackingAction( theG4Track );
298
299 if( fpSteppingManager == 0 )
300 {
301 G4Exception("G4ErrorPropagator::InitG4Track()", "InvalidSetup",
302 FatalException, "G4SteppingManager not initialized yet!");
303 }
304 else
305 {
306 fpSteppingManager->SetInitialStep(theG4Track);
307 }
308
309 // Give SteppingManger the maximum number of processes
310 //
311 fpSteppingManager->GetProcessNumber();
312
313 // Give track the pointer to the Step
314 //
315 theG4Track->SetStep(fpSteppingManager->GetStep());
316
317 // Inform beginning of tracking to physics processes
318 //
319 theG4Track->GetDefinition()->GetProcessManager()->StartTracking(theG4Track);
320
321 initialTS.SetG4Track( theG4Track );
322
323 return theG4Track;
324}
325
326
327//-----------------------------------------------------------------------
328G4int G4ErrorPropagator::MakeSteps( G4ErrorFreeTrajState* currentTS_FREE )
329{
330 G4int ierr = 0;
331
332 //----- Track the particle Step-by-Step while it is alive
333 //
334 theStepLength = 0.;
335
336 while( (theG4Track->GetTrackStatus() == fAlive) ||
337 (theG4Track->GetTrackStatus() == fStopButAlive) )
338 {
339 ierr = MakeOneStep( currentTS_FREE );
340 if( ierr != 0 ) { break; }
341
342 //----- Check if last step for error propagation
343 //
344 if( CheckIfLastStep( theG4Track ) )
345 {
346 break;
347 }
348 }
349 return ierr;
350}
351
352
353//-----------------------------------------------------------------------
355{
356 G4ErrorPropagatorData* g4edata =
358 G4int ierr = 0;
359
360 //---------- Track one step
361#ifdef G4EVERBOSE
362 if(verbose >= 2 )
363 {
364 G4cout << G4endl
365 << "@@@@@@@@@@@@@@@@@@@@@@@@@ NEW STEP " << G4endl;
366 }
367#endif
368
369 theG4Track->IncrementCurrentStepNumber();
370
371 fpSteppingManager->Stepping();
372
373 //---------- Check if Target has been reached (and then set G4ErrorState)
374
375 // G4ErrorPropagationNavigator limits the step if target is closer than
376 // boundary (but the winner process is always "Transportation": then
377 // error propagator will stop the track
378
379 if( theG4Track->GetStep()->GetPostStepPoint()
380 ->GetProcessDefinedStep()->GetProcessName() == "Transportation" )
381 {
382 if( g4edata->GetState()
384 { // target or step length reached
385
386#ifdef G4EVERBOSE
387 if(verbose >= 5 )
388 {
389 G4cout << " transportation determined by geant4e " << G4endl;
390 }
391#endif
393 }
394 else if( g4edata->GetTarget()->GetType() == G4ErrorTarget_GeomVolume )
395 {
397 (G4ErrorGeomVolumeTarget*)(g4edata->GetTarget());
398 if( static_cast<G4ErrorGeomVolumeTarget*>( target )
399 ->TargetReached( theG4Track->GetStep() ) )
400 {
402 }
403 }
404 }
405 else if( theG4Track->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()
406 ->GetProcessName() == "G4ErrorTrackLengthTarget" )
407 {
409 }
410
411 //---------- Propagate error
412
413#ifdef G4EVERBOSE
414 if(verbose >= 2 )
415 {
416 G4cout << " propagating error " << *currentTS_FREE << G4endl;
417 }
418#endif
419 const G4Track* cTrack = const_cast<G4Track*>(theG4Track);
420 ierr = currentTS_FREE->PropagateError( cTrack );
421
422#ifdef G4EVERBOSE
423 if(verbose >= 3 )
424 {
425 G4cout << " PropagateError returns " << ierr << G4endl;
426 }
427#endif
428
429 currentTS_FREE->Update( cTrack );
430
431 theStepLength += theG4Track->GetStepLength();
432
433 if(ierr != 0 )
434 {
435 G4cerr << "ERROR - G4ErrorPropagator:MakeOneStep()" << G4endl
436 << " Error returned: " << ierr << G4endl
437 << " Geant4 tracking will be stopped !" << G4endl;
438 }
439
440 return ierr;
441}
442
443
444//-----------------------------------------------------------------------
447{
448 G4ErrorFreeTrajState* currentTS_FREE = 0;
449
450 //----- Transform the TrajState to Free coordinates if it is OnSurface
451 //
452 if( currentTS->GetTSType() == G4eTS_OS )
453 {
455 static_cast<G4ErrorSurfaceTrajState*>(currentTS);
456 currentTS_FREE = new G4ErrorFreeTrajState( *tssd );
457 }
458 else if( currentTS->GetTSType() == G4eTS_FREE )
459 {
460 currentTS_FREE = static_cast<G4ErrorFreeTrajState*>(currentTS);
461 }
462 else
463 {
464 G4cerr << "ERROR - G4ErrorPropagator::InitFreeTrajState()" << G4endl
465 << "WRONG TrajState " + currentTS->GetTSType() << G4endl;
466 G4Exception("G4ErrorPropagator::InitFreeTrajState()", "InvalidState",
467 FatalException, "WRONG trajectory state !");
468 }
469 return currentTS_FREE;
470}
471
472
473//-----------------------------------------------------------------------
475 G4ErrorFreeTrajState* currentTS_FREE,
476 const G4ErrorTarget* target )
477{
478 G4ErrorPropagatorData* g4edata =
480
481#ifdef G4EVERBOSE
482 if(verbose >= 1 )
483 {
484 G4cout << " G4ErrorPropagator::Propagate: final state "
485 << G4int(g4edata->GetState()) << " TSType "
486 << currentTS->GetTSType() << G4endl;
487 }
488#endif
489
490 if( (currentTS->GetTSType() == G4eTS_FREE) ||
491 (g4edata->GetState() != G4ErrorState_StoppedAtTarget) )
492 {
493 currentTS = currentTS_FREE;
494 }
495 else if( currentTS->GetTSType() == G4eTS_OS )
496 {
497 if( target->GetType() == G4ErrorTarget_TrkL )
498 {
499 G4Exception("G4ErrorPropagator:GetFinalTrajState()",
500 "InvalidSetup", FatalException,
501 "Using a G4ErrorSurfaceTrajState with wrong target");
502 }
503 const G4ErrorTanPlaneTarget* targetWTP =
504 static_cast<const G4ErrorTanPlaneTarget*>(target);
505 *currentTS = G4ErrorSurfaceTrajState(
506 *(static_cast<G4ErrorFreeTrajState*>(currentTS_FREE)),
507 targetWTP->GetTangentPlane( currentTS_FREE->GetPosition() ) );
508#ifdef G4EVERBOSE
509 if(verbose >= 1 )
510 {
511 G4cout << currentTS << " returning tssd " << *currentTS << G4endl;
512 }
513#endif
514 delete currentTS_FREE;
515 }
516}
517
518
519//-------------------------------------------------------------------------
521{
522 G4bool exception = 0;
523 G4bool lastG4eStep = false;
524 G4ErrorPropagatorData* g4edata =
526
527#ifdef G4EVERBOSE
528 if( verbose >= 4 )
529 {
530 G4cout << " G4ErrorPropagator::CheckIfLastStep G4ErrorState= "
531 << G4int(g4edata->GetState()) << G4endl;
532 }
533#endif
534
535 //----- Check if this is the last step: track has reached the target
536 // or the end of world
537 //
539 {
540 lastG4eStep = true;
541#ifdef G4EVERBOSE
542 if(verbose >= 5 )
543 {
544 G4cout << " G4ErrorPropagator::CheckIfLastStep " << lastG4eStep
545 << " " << G4int(g4edata->GetState()) << G4endl;
546 }
547#endif
548 }
549 else if( aTrack->GetNextVolume() == 0 )
550 {
551 //----- If particle is out of world, without finding the G4ErrorTarget,
552 // give a n error/warning
553 //
554 lastG4eStep = true;
555 if( exception )
556 {
557 G4cerr << "ERROR - G4ErrorPropagator::CheckIfLastStep()" << G4endl
558 << " Track extrapolated until end of World" << G4endl
559 << " without finding the defined target " << G4endl;
560 G4Exception("G4ErrorPropagator::CheckIfLastStep()",
561 "InvalidSetup", FatalException,
562 "Track extrapolated without finding the defined target.");
563 }
564 else
565 {
566 if( verbose >= 1 )
567 {
568 G4cerr << "ERROR - G4ErrorPropagator::CheckIfLastStep()" << G4endl
569 << " Track extrapolated until end of World" << G4endl
570 << " without finding the defined target " << G4endl;
571 }
572 }
573 } //----- not last step from G4e, but track is stopped (energy exhausted)
574 else if( aTrack->GetTrackStatus() == fStopAndKill )
575 {
576 if( exception )
577 {
578 G4cerr << "ERROR - G4ErrorPropagator::CheckIfLastStep()" << G4endl
579 << " Track extrapolated until energy is exhausted" << G4endl
580 << " without finding the defined target !" << G4endl;
581 G4Exception("G4ErrorPropagator::CheckIfLastStep()",
582 "InvalidSetup", FatalException,
583 "Track extrapolated without finding the defined target.");
584 }
585 else
586 {
587 if( verbose >= 1 )
588 {
589 G4cerr << "ERROR - G4ErrorPropagator::CheckIfLastStep()" << G4endl
590 << " Track extrapolated until energy is exhausted" << G4endl
591 << " without finding the defined target !" << G4endl;
592 }
593 lastG4eStep = 1;
594 }
595 }
596
597#ifdef G4EVERBOSE
598 if( verbose >= 5 )
599 {
600 G4cout << " return CheckIfLastStep " << lastG4eStep << G4endl;
601 }
602#endif
603
604 return lastG4eStep;
605}
606
607
608//---------------------------------------------------------------------------
610{
611 const G4UserTrackingAction* fpUserTrackingAction =
613 if( fpUserTrackingAction != 0 )
614 {
615 const_cast<G4UserTrackingAction*>(fpUserTrackingAction)
616 ->PreUserTrackingAction((fpTrack) );
617 }
618}
619
620
621//---------------------------------------------------------------------------
623{
624 const G4UserTrackingAction* fpUserTrackingAction =
626 if( fpUserTrackingAction != 0 )
627 {
628 const_cast<G4UserTrackingAction*>(fpUserTrackingAction)
629 ->PostUserTrackingAction((fpTrack) );
630 }
631}
@ G4State_GeomClosed
@ G4ErrorState_TargetCloserThanBoundary
@ G4ErrorState_PreInit
@ G4ErrorState_StoppedAtTarget
@ G4ErrorTarget_GeomVolume
@ G4ErrorTarget_TrkL
@ G4eTS_FREE
@ G4eTS_OS
@ FatalException
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
@ fAlive
@ fStopAndKill
@ fStopButAlive
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4DLLIMPORT std::ostream G4cout
void SetCharge(G4double charge)
void SetPolarization(G4double polX, G4double polY, G4double polZ)
virtual G4int Update(const G4Track *aTrack)
virtual G4int PropagateError(const G4Track *aTrack)
virtual G4bool TargetReached(const G4Step *aStep)
const G4ErrorTarget * GetTarget(G4bool mustExist=0) const
G4ErrorState GetState() const
static G4ErrorPropagatorData * GetErrorPropagatorData()
void SetMode(G4ErrorMode mode)
void SetState(G4ErrorState sta)
void SetTarget(const G4ErrorTarget *target)
G4ErrorMode GetMode() const
G4Track * InitG4Track(G4ErrorTrajState &initialTS)
void InvokePostUserTrackingAction(G4Track *fpTrack)
void InvokePreUserTrackingAction(G4Track *fpTrack)
void GetFinalTrajState(G4ErrorTrajState *currentTS, G4ErrorFreeTrajState *currentTS_FREE, const G4ErrorTarget *target)
G4int PropagateOneStep(G4ErrorTrajState *currentTS)
G4bool CheckIfLastStep(G4Track *aTrack)
G4ErrorFreeTrajState * InitFreeTrajState(G4ErrorTrajState *currentTS)
G4int MakeOneStep(G4ErrorFreeTrajState *currentTS_FREE)
G4int Propagate(G4ErrorTrajState *currentTS, const G4ErrorTarget *target, G4ErrorMode mode=G4ErrorMode_PropForwards)
virtual G4Plane3D GetTangentPlane(const G4ThreeVector &point) const =0
virtual void Dump(const G4String &msg) const =0
G4ErrorTargetType GetType() const
void SetG4Track(G4Track *trk)
G4Vector3D GetMomentum() const
const G4String & GetParticleType() const
virtual G4eTSType GetTSType() const
G4Point3D GetPosition() const
G4UserTrackingAction * GetUserTrackingAction()
static G4EventManager * GetEventManager()
G4TrackingManager * GetTrackingManager() const
G4ProcessManager * GetProcessManager() const
G4double GetPDGCharge() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
void StartTracking(G4Track *aTrack=0)
G4ApplicationState GetCurrentState() const
static G4StateManager * GetStateManager()
const G4VProcess * GetProcessDefinedStep() const
G4StepPoint * GetPostStepPoint() const
G4StepStatus Stepping()
void SetInitialStep(G4Track *valueTrack)
G4Step * GetStep() const
G4TrackStatus GetTrackStatus() const
void SetStep(const G4Step *aValue)
G4VPhysicalVolume * GetNextVolume() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetStepLength() const
void IncrementCurrentStepNumber()
const G4Step * GetStep() const
void SetParentID(const G4int aValue)
G4SteppingManager * GetSteppingManager() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41