#include <G4ErrorPropagator.hh>
Definition at line 54 of file G4ErrorPropagator.hh.
◆ G4ErrorPropagator()
G4ErrorPropagator::G4ErrorPropagator |
( |
| ) |
|
Definition at line 57 of file G4ErrorPropagator.cc.
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
68 thePropIsInitialized = false;
69}
G4DLLIMPORT std::ostream G4cout
static G4EventManager * GetEventManager()
G4TrackingManager * GetTrackingManager() const
G4SteppingManager * GetSteppingManager() const
◆ ~G4ErrorPropagator()
G4ErrorPropagator::~G4ErrorPropagator |
( |
| ) |
|
|
inline |
◆ CheckIfLastStep()
Definition at line 520 of file G4ErrorPropagator.cc.
521{
523 G4bool lastG4eStep =
false;
526
527#ifdef G4EVERBOSE
528 if( verbose >= 4 )
529 {
530 G4cout <<
" G4ErrorPropagator::CheckIfLastStep G4ErrorState= "
532 }
533#endif
534
535
536
537
539 {
540 lastG4eStep = true;
541#ifdef G4EVERBOSE
542 if(verbose >= 5 )
543 {
544 G4cout <<
" G4ErrorPropagator::CheckIfLastStep " << lastG4eStep
546 }
547#endif
548 }
550 {
551
552
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()",
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 }
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()",
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}
@ G4ErrorState_StoppedAtTarget
G4DLLIMPORT std::ostream G4cerr
G4ErrorState GetState() const
static G4ErrorPropagatorData * GetErrorPropagatorData()
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetNextVolume() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
◆ GetFinalTrajState()
Definition at line 474 of file G4ErrorPropagator.cc.
477{
480
481#ifdef G4EVERBOSE
482 if(verbose >= 1 )
483 {
484 G4cout <<
" G4ErrorPropagator::Propagate: final state "
487 }
488#endif
489
492 {
493 currentTS = currentTS_FREE;
494 }
496 {
498 {
499 G4Exception(
"G4ErrorPropagator:GetFinalTrajState()",
501 "Using a G4ErrorSurfaceTrajState with wrong target");
502 }
508#ifdef G4EVERBOSE
509 if(verbose >= 1 )
510 {
511 G4cout << currentTS <<
" returning tssd " << *currentTS <<
G4endl;
512 }
513#endif
514 delete currentTS_FREE;
515 }
516}
virtual G4Plane3D GetTangentPlane(const G4ThreeVector &point) const =0
G4ErrorTargetType GetType() const
virtual G4eTSType GetTSType() const
G4Point3D GetPosition() const
Referenced by Propagate(), and PropagateOneStep().
◆ GetInitialTrajState()
◆ GetStepLength()
G4double G4ErrorPropagator::GetStepLength |
( |
| ) |
const |
|
inline |
◆ InitFreeTrajState()
◆ InitG4Track()
Definition at line 250 of file G4ErrorPropagator.cc.
251{
252 if( verbose >= 5 ) {
G4cout <<
"InitG4Track " <<
G4endl; }
253
254
255
259 if( particle == 0)
260 {
261 G4cerr <<
"ERROR - G4ErrorPropagator::InitG4Track()" <<
G4endl
262 <<
" Particle type not defined " + partType <<
G4endl;
263 G4Exception(
"G4ErrorPropagator::InitG4Track()",
"InvalidSetup",
265 }
266
269
271
272
273
275 {
277 }
278 else
279 {
281 }
282
283
284
287
288#ifdef G4EVERBOSE
289 if(verbose >= 3)
290 {
291 G4cout <<
" G4ErrorPropagator new track of energy: "
293 }
294#endif
295
296
298
299 if( fpSteppingManager == 0 )
300 {
301 G4Exception(
"G4ErrorPropagator::InitG4Track()",
"InvalidSetup",
303 }
304 else
305 {
307 }
308
309
310
312
313
314
316
317
318
320
322
323 return theG4Track;
324}
void SetCharge(G4double charge)
void SetPolarization(G4double polX, G4double polY, G4double polZ)
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=0)
void SetInitialStep(G4Track *valueTrack)
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 622 of file G4ErrorPropagator.cc.
623{
626 if( fpUserTrackingAction != 0 )
627 {
629 ->PostUserTrackingAction((fpTrack) );
630 }
631}
G4UserTrackingAction * GetUserTrackingAction()
Referenced by Propagate().
◆ InvokePreUserTrackingAction()
void G4ErrorPropagator::InvokePreUserTrackingAction |
( |
G4Track * |
fpTrack | ) |
|
Definition at line 609 of file G4ErrorPropagator.cc.
610{
613 if( fpUserTrackingAction != 0 )
614 {
616 ->PreUserTrackingAction((fpTrack) );
617 }
618}
Referenced by InitG4Track().
◆ MakeOneStep()
Definition at line 354 of file G4ErrorPropagator.cc.
355{
359
360
361#ifdef G4EVERBOSE
362 if(verbose >= 2 )
363 {
365 <<
"@@@@@@@@@@@@@@@@@@@@@@@@@ NEW STEP " <<
G4endl;
366 }
367#endif
368
370
372
373
374
375
376
377
378
381 {
384 {
385
386#ifdef G4EVERBOSE
387 if(verbose >= 5 )
388 {
389 G4cout <<
" transportation determined by geant4e " <<
G4endl;
390 }
391#endif
393 }
395 {
400 {
402 }
403 }
404 }
407 {
409 }
410
411
412
413#ifdef G4EVERBOSE
414 if(verbose >= 2 )
415 {
416 G4cout <<
" propagating error " << *currentTS_FREE <<
G4endl;
417 }
418#endif
421
422#ifdef G4EVERBOSE
423 if(verbose >= 3 )
424 {
426 }
427#endif
428
429 currentTS_FREE->
Update( cTrack );
430
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}
@ G4ErrorState_TargetCloserThanBoundary
@ G4ErrorTarget_GeomVolume
virtual G4int Update(const G4Track *aTrack)
virtual G4int PropagateError(const G4Track *aTrack)
const G4ErrorTarget * GetTarget(G4bool mustExist=0) const
void SetState(G4ErrorState sta)
const G4VProcess * GetProcessDefinedStep() const
G4StepPoint * GetPostStepPoint() const
G4double GetStepLength() const
void IncrementCurrentStepNumber()
const G4Step * GetStep() const
const G4String & GetProcessName() const
Referenced by PropagateOneStep().
◆ Propagate()
Definition at line 73 of file G4ErrorPropagator.cc.
75{
76
77
79
82
83
84
86 {
87 G4cerr <<
"ERROR - G4ErrorPropagator::Propagate()" <<
G4endl
88 << " Energy too low to be propagated: "
90 return -3;
91 }
92
94
95#ifdef G4EVERBOSE
96 if( verbose >= 1 )
97 {
98 G4cout <<
" =====> starting GEANT4E tracking for "
101 }
102 if(verbose >= 1 )
103 {
105 }
106
107 if( verbose >= 3 )
108 {
109 G4cout <<
" G4ErrorPropagator::Propagate initialTS ";
112 }
113#endif
114
116
117
118
119 if( theG4Track != 0 ) { delete theG4Track; }
121
122
123
125
126
127
128 ierr = MakeSteps( currentTS_FREE );
129
130
131
132
134 {
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
149 }
150 }
151 else
152 {
154 }
155
156#ifdef G4EVERBOSE
157 if( verbose >= 1 )
158 {
159 G4cout <<
" G4ErrorPropagator: propagation ended " <<
G4endl;
160 }
161 if( verbose >= 2 )
162 {
164 }
165#endif
166
167
168
170
172
173
174
175 return ierr;
176}
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
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()
Definition at line 180 of file G4ErrorPropagator.cc.
181{
184
188 {
189 G4cout <<
"ERROR - G4ErrorPropagator::PropagateOneStep()" <<
G4endl
190 << " Called before initialization is done for this track."
192 << " Please call G4ErrorPropagatorManager::InitGeant4e()."
194 G4Exception(
"G4ErrorPropagator::PropagateOneStep()",
196 "Called before initialization is done for this track!");
197 }
198
199
200
202
203
204
206 {
207 G4cerr <<
"ERROR - G4ErrorPropagator::PropagateOneStep()" <<
G4endl
208 << " Energy too low to be propagated: "
210 return -3;
211 }
212
213#ifdef G4EVERBOSE
214 if( verbose >= 1 )
215 {
216 G4cout <<
" =====> starting GEANT4E tracking for "
219 }
220 if( verbose >= 3 )
221 {
222 G4cout <<
" G4ErrorPropagator::Propagate initialTS ";
224 }
225#endif
226
227
228
229 if( theStepN == 0 ) { theG4Track =
InitG4Track( *currentTS ); }
230
231 theStepN++;
232
233
234
236
237
238
240
241
242
244
245 return ierr;
246}
G4int MakeOneStep(G4ErrorFreeTrajState *currentTS_FREE)
G4ApplicationState GetCurrentState() const
static G4StateManager * GetStateManager()
Referenced by G4ErrorPropagatorManager::PropagateOneStep().
◆ SetStepLength()
◆ SetStepN()
void G4ErrorPropagator::SetStepN |
( |
const G4int |
sn | ) |
|
|
inline |
The documentation for this class was generated from the following files: