Geant4-11
G4FastStep.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//
27//
28//---------------------------------------------------------------
29//
30// G4FastStep.cc
31//
32// Description:
33// Encapsulates a G4ParticleChange and insure friendly interface
34// methods to manage the primary/secondaries final state for
35// Fast Simulation Models.
36//
37// History:
38// Oct 97: Verderi && MoraDeFreitas - First Implementation.
39// Apr 98: MoraDeFreitas - G4FastStep becomes the G4ParticleChange
40// for the Fast Simulation Process.
41//
42//---------------------------------------------------------------
43
44#include "G4FastStep.hh"
45
46#include "G4SystemOfUnits.hh"
47#include "G4UnitsTable.hh"
48#include "G4Track.hh"
49#include "G4Step.hh"
50#include "G4TrackFastVector.hh"
51#include "G4DynamicParticle.hh"
52
53void G4FastStep::Initialize(const G4FastTrack& fastTrack)
54{
55 // keeps the fastTrack reference
56 fFastTrack=&fastTrack;
57
58 // currentTrack will be used to Initialize the other data members
59 const G4Track& currentTrack = *(fFastTrack->GetPrimaryTrack());
60
61 // use base class's method at first
63
64 // set Energy/Momentum etc. equal to those of the parent particle
65 const G4DynamicParticle* pParticle = currentTrack.GetDynamicParticle();
66 theEnergyChange = pParticle->GetKineticEnergy();
70
71 // set Position/Time etc. equal to those of the parent track
72 thePositionChange = currentTrack.GetPosition();
73 theTimeChange = currentTrack.GetGlobalTime();
74
75 // switch off stepping hit invokation by default:
77
78 // event biasing weigth:
79 theWeightChange = currentTrack.GetWeight();
80}
81
82//----------------------------------------
83// -- Set the StopAndKilled signal
84// -- and put kinetic energy to 0.0. in the
85// -- G4ParticleChange.
86//----------------------------------------
88{
91}
92
93//--------------------
94//
95//--------------------
96void
99 G4bool localCoordinates)
100{
101 // Compute the position coordinate in global
102 // reference system if needed ...
103 G4ThreeVector globalPosition = position;
104 if (localCoordinates)
105 globalPosition = fFastTrack->GetInverseAffineTransformation()->
106 TransformPoint(position);
107 // ...and feed the globalPosition:
108 thePositionChange = globalPosition;
109}
110
111void
114 G4bool localCoordinates)
115{
117}
118
119//--------------------
120//
121//--------------------
122void
125 G4bool localCoordinates)
126{
127 // Compute the momentum in global reference
128 // system if needed ...
129 G4ThreeVector globalMomentum = momentum;
130 if (localCoordinates)
131 globalMomentum = fFastTrack->GetInverseAffineTransformation()->
132 TransformAxis(momentum);
133 // ...and feed the globalMomentum (ensuring unitarity)
134 SetMomentumChange(globalMomentum.unit());
135}
136
137void
140 G4bool localCoordinates)
141{
142 ProposePrimaryTrackFinalMomentumDirection(momentum, localCoordinates);
143}
144
145//--------------------
146//
147//--------------------
148void
151 const G4ThreeVector &direction,
152 G4bool localCoordinates)
153{
154 // Compute global direction if needed...
155 G4ThreeVector globalDirection = direction;
156 if (localCoordinates)
157 globalDirection =fFastTrack->GetInverseAffineTransformation()->
158 TransformAxis(direction);
159 // ...and feed the globalMomentum (ensuring unitarity)
160 SetMomentumChange(globalDirection.unit());
162}
163
164void
167 const G4ThreeVector &direction,
168 G4bool localCoordinates)
169{
170 ProposePrimaryTrackFinalKineticEnergyAndDirection(kineticEnergy, direction, localCoordinates);
171}
172
173//--------------------
174//
175//--------------------
176void
179 G4bool localCoordinates)
180{
181 // Compute polarization in global system if needed:
182 G4ThreeVector globalPolarization(polarization);
183 if (localCoordinates)
184 globalPolarization = fFastTrack->GetInverseAffineTransformation()->
185 TransformAxis(globalPolarization);
186 // Feed the particle globalPolarization:
187 thePolarizationChange = globalPolarization;
188}
189
190void
193 G4bool localCoordinates)
194{
195 ProposePrimaryTrackFinalPolarization(polarization, localCoordinates);
196}
197
198//--------------------
199//
200//--------------------
203 G4ThreeVector polarization,
205 G4double time,
206 G4bool localCoordinates )
207{
208 G4DynamicParticle dummyDynamics(dynamics);
209
210 // ------------------------------------------
211 // Add the polarization to the dummyDynamics:
212 // ------------------------------------------
213 dummyDynamics.SetPolarization(polarization.x(),
214 polarization.y(),
215 polarization.z());
216
217 return CreateSecondaryTrack(dummyDynamics, position, time, localCoordinates);
218}
219
220//--------------------
221//
222//--------------------
226 G4double time,
227 G4bool localCoordinates )
228{
229 // ----------------------------------------
230 // Quantities in global coordinates system.
231 //
232 // The allocated globalDynamics is deleted
233 // by the destructor of the G4Track.
234 // ----------------------------------------
235 G4DynamicParticle* globalDynamics =
236 new G4DynamicParticle(dynamics);
237 G4ThreeVector globalPosition(position);
238
239 // -----------------------------------
240 // Convert to global system if needed:
241 // -----------------------------------
242 if (localCoordinates)
243 {
244 // -- Momentum Direction:
245 globalDynamics->SetMomentumDirection(fFastTrack->
246 GetInverseAffineTransformation()->
247 TransformAxis(globalDynamics->
248 GetMomentumDirection()));
249 // -- Polarization:
250 G4ThreeVector globalPolarization;
251 globalPolarization = fFastTrack->GetInverseAffineTransformation()->
252 TransformAxis(globalDynamics->GetPolarization());
253 globalDynamics->SetPolarization(
254 globalPolarization.x(),
255 globalPolarization.y(),
256 globalPolarization.z()
257 );
258
259 // -- Position:
260 globalPosition = fFastTrack->GetInverseAffineTransformation()->
261 TransformPoint(globalPosition);
262 }
263
264 //-------------------------------------
265 // Create the G4Track of the secondary:
266 //-------------------------------------
267 G4Track* secondary = new G4Track(
268 globalDynamics,
269 time,
270 globalPosition
271 );
272
273 //-------------------------------
274 // and feed the changes:
275 //-------------------------------
276 AddSecondary(secondary);
277
278 //--------------------------------------
279 // returns the pointer on the secondary:
280 //--------------------------------------
281 return secondary;
282}
283
284// G4FastStep should never be Initialized in this way
285// but we must define it to avoid warnings.
287{
288 G4ExceptionDescription tellWhatIsWrong;
289 tellWhatIsWrong << "G4FastStep can be initialised only through G4FastTrack."
290 << G4endl;
291 G4Exception("G4FastStep::Initialize(const G4Track&)",
292 "FastSim005",
294 tellWhatIsWrong);
295}
296
299{
300 if (verboseLevel>2)
301 {
302 G4cerr << "G4FastStep::G4FastStep()" << G4endl;
303 }
304}
305
307{
308 if (verboseLevel>2)
309 {
310 G4cerr << "G4FastStep::~G4FastStep()" << G4endl;
311 }
312}
313
314// copy and assignment operators are implemented as "shallow copy"
317{
318 *this = right;
319}
320
321
323{
324 if (this != &right)
325 {
341 fFastTrack = right.fFastTrack;
342 }
343 return *this;
344}
345
346
347
348
349
351{
352 return ((G4VParticleChange *)this == (G4VParticleChange *) &right);
353}
354
356{
357 return ((G4VParticleChange *)this != (G4VParticleChange *) &right);
358}
359
360//----------------------------------------------------------------
361// methods for updating G4Step
362//
363
365{
366 // A physics process always calculates the final state of the particle
367
368 // Take note that the return type of GetMomentumChange is a
369 // pointer to G4ParticleMometum. Also it is a normalized
370 // momentum vector.
371
372 // G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
373 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
374 G4Track* aTrack = pStep->GetTrack();
375 // G4double mass = aTrack->GetDynamicParticle()->GetMass();
376
377 // update kinetic energy and momentum direction
379 pPostStepPoint->SetKineticEnergy( theEnergyChange );
380
381 // update polarization
382 pPostStepPoint->SetPolarization( thePolarizationChange );
383
384 // update position and time
385 pPostStepPoint->SetPosition( thePositionChange );
386 pPostStepPoint->SetGlobalTime( theTimeChange );
387 pPostStepPoint->AddLocalTime( theTimeChange
388 - aTrack->GetGlobalTime());
389 pPostStepPoint->SetProperTime( theProperTimeChange );
390
391 // update weight
392 pPostStepPoint->SetWeight( theWeightChange );
393
394 if (debugFlag) CheckIt(*aTrack);
395
396
397 // Update the G4Step specific attributes
398 return UpdateStepInfo(pStep);
399}
400
402{
403 // A physics process always calculates the final state of the particle
404
405 // G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
406 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
407 G4Track* aTrack = pStep->GetTrack();
408 // G4double mass = aTrack->GetDynamicParticle()->GetMass();
409
410 // update kinetic energy and momentum direction
412 pPostStepPoint->SetKineticEnergy( theEnergyChange );
413
414 // update polarization
415 pPostStepPoint->SetPolarization( thePolarizationChange );
416
417 // update position and time
418 pPostStepPoint->SetPosition( thePositionChange );
419 pPostStepPoint->SetGlobalTime( theTimeChange );
420 pPostStepPoint->AddLocalTime( theTimeChange
421 - aTrack->GetGlobalTime());
422 pPostStepPoint->SetProperTime( theProperTimeChange );
423
424 // update weight
425 pPostStepPoint->SetWeight( theWeightChange );
426
427 if (debugFlag) CheckIt(*aTrack);
428
429 // Update the G4Step specific attributes
430 return UpdateStepInfo(pStep);
431}
432
433//----------------------------------------------------------------
434// methods for printing messages
435//
436
438{
439// use base-class DumpInfo
441
442 G4cout << " Position - x (mm) : " << G4BestUnit( thePositionChange.x(), "Length" ) << G4endl;
443 G4cout << " Position - y (mm) : " << G4BestUnit( thePositionChange.y(), "Length" ) << G4endl;
444 G4cout << " Position - z (mm) : " << G4BestUnit( thePositionChange.z(), "Length" ) << G4endl;
445 G4cout << " Time (ns) : " << G4BestUnit( theTimeChange, "Time" ) << G4endl;
446 G4cout << " Proper Time (ns) : " << G4BestUnit( theProperTimeChange, "Time" ) << G4endl;
447 G4int olprc = G4cout.precision(3);
448 G4cout << " Momentum Direct - x : " << std::setw(20) << theMomentumChange.x() << G4endl;
449 G4cout << " Momentum Direct - y : " << std::setw(20) << theMomentumChange.y() << G4endl;
450 G4cout << " Momentum Direct - z : " << std::setw(20) << theMomentumChange.z() << G4endl;
451 G4cout.precision(olprc);
452 G4cout << " Kinetic Energy (MeV): " << G4BestUnit( theEnergyChange, "Energy" ) << G4endl;
453 G4cout.precision(3);
454 G4cout << " Polarization - x : " << std::setw(20) << thePolarizationChange.x() << G4endl;
455 G4cout << " Polarization - y : " << std::setw(20) << thePolarizationChange.y() << G4endl;
456 G4cout << " Polarization - z : " << std::setw(20) << thePolarizationChange.z() << G4endl;
457 G4cout.precision(olprc);
458}
459
461{
462 //
463 // In the G4FastStep::CheckIt
464 // We only check a bit
465 //
466 // If the user violates the energy,
467 // We don't care, we agree.
468 //
469 // But for theMomentumDirectionChange,
470 // We do pay attention.
471 // And if too large is its range,
472 // We issue an Exception.
473 //
474 //
475 // It means, the G4FastStep::CheckIt issues an exception only for the
476 // theMomentumDirectionChange which should be an unit vector
477 // and it corrects it because it could cause problems for the ulterior
478 // tracking.For the rest, only warning are issued.
479
480 G4bool itsOK = true;
481 G4bool exitWithError = false;
482 G4double accuracy;
483
484 // Energy should not be larger than the initial value
485 accuracy = ( theEnergyChange - aTrack.GetKineticEnergy())/MeV;
486 if (accuracy > GetAccuracyForWarning())
487 {
489 ed << "The energy becomes larger than the initial value, difference = " << accuracy << " MeV" << G4endl;
490 G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
491 "FastSim006",
492 JustWarning, ed);
493 itsOK = false;
494 if (accuracy > GetAccuracyForException()) {exitWithError = true;}
495 }
496
497 G4bool itsOKforMomentum = true;
498 if ( theEnergyChange >0.)
499 {
500 accuracy = std::abs(theMomentumChange.mag2()-1.0);
501 if (accuracy > GetAccuracyForWarning())
502 {
504 ed << "The Momentum Change is not a unit vector, difference = " << accuracy << G4endl;
505 G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
506 "FastSim007",
507 JustWarning, ed);
508 itsOK = itsOKforMomentum = false;
509 if (accuracy > GetAccuracyForException()) {exitWithError = true;}
510 }
511 }
512
513 accuracy = (aTrack.GetGlobalTime()- theTimeChange)/ns;
514 if (accuracy > GetAccuracyForWarning())
515 {
517 ed << "The global time is getting backward, difference = " << accuracy << " ns" << G4endl;
518 G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
519 "FastSim008",
520 JustWarning, ed);
521 itsOK = false;
522 }
523
524 accuracy = (aTrack.GetProperTime() - theProperTimeChange )/ns;
525 if (accuracy > GetAccuracyForWarning())
526 {
528 ed << "The proper time is getting backward, difference = " << accuracy << " ns" << G4endl;
529 G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
530 "FastSim009",
531 JustWarning, ed);
532 itsOK = false;
533 }
534
535 if (!itsOK)
536 {
537 G4cout << "ERROR - G4FastStep::CheckIt() " << G4endl;
538 G4cout << " Pointer : " << this << G4endl ;
539 DumpInfo();
540 }
541
542 // Exit with error
543 if (exitWithError)
544 {
546 ed << "An inaccuracy in G4FastStep is beyond tolerance." << G4endl;
547 G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
548 "FastSim010",
549 FatalException, ed);
550 }
551
552 //correction for Momentum only.
553 if (!itsOKforMomentum) {
556 }
557
558 itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
559 return itsOK;
560}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static constexpr double MeV
Definition: G4SIunits.hh:200
@ AvoidHitInvocation
#define G4BestUnit(a, b)
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double mag2() const
double y() const
double mag() const
void SetPolarization(const G4ThreeVector &)
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetProperTime() const
const G4ThreeVector & GetPolarization() const
G4FastStep & operator=(const G4FastStep &right)
Definition: G4FastStep.cc:322
void SetPrimaryTrackFinalKineticEnergyAndDirection(G4double, const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:166
G4Step * UpdateStepForPostStep(G4Step *Step)
Definition: G4FastStep.cc:364
G4ThreeVector thePolarizationChange
Definition: G4FastStep.hh:354
virtual ~G4FastStep()
Definition: G4FastStep.cc:306
void KillPrimaryTrack()
Definition: G4FastStep.cc:87
void SetPrimaryTrackFinalPosition(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:113
G4bool operator==(const G4FastStep &right) const
Definition: G4FastStep.cc:350
G4bool operator!=(const G4FastStep &right) const
Definition: G4FastStep.cc:355
G4double theEnergyChange
Definition: G4FastStep.hh:357
void SetPrimaryTrackFinalPolarization(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:192
void DumpInfo() const
Definition: G4FastStep.cc:437
G4double theProperTimeChange
Definition: G4FastStep.hh:366
void ProposePrimaryTrackFinalKineticEnergyAndDirection(G4double, const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:150
void SetPrimaryTrackFinalMomentum(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:139
G4bool CheckIt(const G4Track &)
Definition: G4FastStep.cc:460
void ProposePrimaryTrackFinalPosition(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:98
G4double theWeightChange
Definition: G4FastStep.hh:372
G4ParticleMomentum theMomentumChange
Definition: G4FastStep.hh:351
void SetPrimaryTrackFinalKineticEnergy(G4double)
G4ThreeVector thePositionChange
Definition: G4FastStep.hh:360
void SetMomentumChange(G4double Px, G4double Py, G4double Pz)
const G4FastTrack * fFastTrack
Definition: G4FastStep.hh:369
G4double theTimeChange
Definition: G4FastStep.hh:363
void ProposePrimaryTrackFinalPolarization(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:178
void ProposePrimaryTrackFinalMomentumDirection(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:124
G4Step * UpdateStepForAtRest(G4Step *Step)
Definition: G4FastStep.cc:401
G4Track * CreateSecondaryTrack(const G4DynamicParticle &, G4ThreeVector, G4ThreeVector, G4double, G4bool localCoordinates=true)
Definition: G4FastStep.cc:202
void Initialize(const G4FastTrack &)
Definition: G4FastStep.cc:53
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:206
const G4AffineTransform * GetInverseAffineTransformation() const
Definition: G4FastTrack.hh:236
void SetKineticEnergy(const G4double aValue)
void SetWeight(G4double aValue)
void SetProperTime(const G4double aValue)
void SetGlobalTime(const G4double aValue)
void SetPosition(const G4ThreeVector &aValue)
void AddLocalTime(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPostStepPoint() const
G4double GetWeight() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4double GetProperTime() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
virtual G4bool CheckIt(const G4Track &)
void ProposeTrackStatus(G4TrackStatus status)
G4TrackStatus theStatusChange
G4TrackFastVector * theListOfSecondaries
G4SteppingControl theSteppingControlFlag
virtual void Initialize(const G4Track &)
G4double GetAccuracyForException() const
void AddSecondary(G4Track *aSecondary)
G4double GetAccuracyForWarning() const
virtual void DumpInfo() const
G4VParticleChange & operator=(const G4VParticleChange &right)
G4Step * UpdateStepInfo(G4Step *Step)
#define position
Definition: xmlparse.cc:622
#define ns
Definition: xmlparse.cc:614