G4VParticleChange Class Reference

#include <G4VParticleChange.hh>

Inheritance diagram for G4VParticleChange:

G4FastStep G4ParticleChange G4ParticleChangeForDecay G4ParticleChangeForGamma G4ParticleChangeForLoss G4ParticleChangeForMSC G4ParticleChangeForTransport G4ParticleChangeForRadDecay

Public Member Functions

 G4VParticleChange ()
virtual ~G4VParticleChange ()
G4bool operator== (const G4VParticleChange &right) const
G4bool operator!= (const G4VParticleChange &right) const
virtual G4StepUpdateStepForAtRest (G4Step *Step)
virtual G4StepUpdateStepForAlongStep (G4Step *Step)
virtual G4StepUpdateStepForPostStep (G4Step *Step)
virtual void Initialize (const G4Track &)
G4double GetTrueStepLength () const
void ProposeTrueStepLength (G4double truePathLength)
G4double GetLocalEnergyDeposit () const
void ProposeLocalEnergyDeposit (G4double anEnergyPart)
G4double GetNonIonizingEnergyDeposit () const
void ProposeNonIonizingEnergyDeposit (G4double anEnergyPart)
G4TrackStatus GetTrackStatus () const
void ProposeTrackStatus (G4TrackStatus status)
G4SteppingControl GetSteppingControl () const
void ProposeSteppingControl (G4SteppingControl StepControlFlag)
G4bool GetFirstStepInVolume () const
G4bool GetLastStepInVolume () const
void ProposeFirstStepInVolume (G4bool flag)
void ProposeLastStepInVolume (G4bool flag)
void Clear ()
void SetNumberOfSecondaries (G4int totSecondaries)
G4int GetNumberOfSecondaries () const
G4TrackGetSecondary (G4int anIndex) const
void AddSecondary (G4Track *aSecondary)
G4double GetWeight () const
G4double GetParentWeight () const
void ProposeWeight (G4double finalWeight)
void ProposeParentWeight (G4double finalWeight)
void SetSecondaryWeightByProcess (G4bool)
G4bool IsSecondaryWeightSetByProcess () const
void SetParentWeightByProcess (G4bool)
G4bool IsParentWeightSetByProcess () const
virtual void DumpInfo () const
void SetVerboseLevel (G4int vLevel)
G4int GetVerboseLevel () const
virtual G4bool CheckIt (const G4Track &)
void ClearDebugFlag ()
void SetDebugFlag ()
G4bool GetDebugFlag () const

Protected Member Functions

 G4VParticleChange (const G4VParticleChange &right)
G4VParticleChangeoperator= (const G4VParticleChange &right)
G4StepUpdateStepInfo (G4Step *Step)
void InitializeTrueStepLength (const G4Track &)
void InitializeLocalEnergyDeposit (const G4Track &)
void InitializeSteppingControl (const G4Track &)
void InitializeParentWeight (const G4Track &)
void InitializeParentGlobalTime (const G4Track &)
void InitializeStatusChange (const G4Track &)
void InitializeSecondaries (const G4Track &)
void InitializeStepInVolumeFlags (const G4Track &)
G4bool CheckSecondary (G4Track &)
G4double GetAccuracyForWarning () const
G4double GetAccuracyForException () const

Protected Attributes

G4TrackFastVectortheListOfSecondaries
G4int theNumberOfSecondaries
G4int theSizeOftheListOfSecondaries
G4TrackStatus theStatusChange
G4SteppingControl theSteppingControlFlag
G4double theLocalEnergyDeposit
G4double theNonIonizingEnergyDeposit
G4double theTrueStepLength
G4bool theFirstStepInVolume
G4bool theLastStepInVolume
G4double theParentWeight
G4bool isParentWeightProposed
G4bool fSetSecondaryWeightByProcess
G4double theParentGlobalTime
G4int verboseLevel
G4bool debugFlag

Static Protected Attributes

static const G4double accuracyForWarning = 1.0e-9
static const G4double accuracyForException = 0.001

Detailed Description

Definition at line 94 of file G4VParticleChange.hh.


Constructor & Destructor Documentation

G4VParticleChange::G4VParticleChange (  ) 

Definition at line 48 of file G4VParticleChange.cc.

References debugFlag, and theListOfSecondaries.

00049   :theListOfSecondaries(0),
00050    theNumberOfSecondaries(0),
00051    theSizeOftheListOfSecondaries(G4TrackFastVectorSize),
00052    theStatusChange(fAlive),
00053    theSteppingControlFlag(NormalCondition),     
00054    theLocalEnergyDeposit(0.0),
00055    theNonIonizingEnergyDeposit(0.0),
00056    theTrueStepLength(0.0),
00057    theFirstStepInVolume(false),
00058    theLastStepInVolume(false),
00059    theParentWeight(1.0),
00060    isParentWeightProposed(false),
00061    fSetSecondaryWeightByProcess(false),
00062    theParentGlobalTime(0.0),
00063    verboseLevel(1), 
00064    debugFlag(false)
00065 {
00066 #ifdef G4VERBOSE
00067   // activate CHeckIt if in VERBOSE mode
00068   debugFlag = true;
00069 #endif
00070   theListOfSecondaries = new G4TrackFastVector();
00071 }

G4VParticleChange::~G4VParticleChange (  )  [virtual]

Definition at line 73 of file G4VParticleChange.cc.

References G4cout, theListOfSecondaries, theNumberOfSecondaries, and verboseLevel.

00073                                       {
00074   // check if tracks still exist in theListOfSecondaries
00075   if (theNumberOfSecondaries>0) {
00076 #ifdef G4VERBOSE
00077     if (verboseLevel>0) {
00078       G4cout << "G4VParticleChange::~G4VParticleChange() Warning  ";
00079       G4cout << "theListOfSecondaries is not empty ";
00080     }
00081 #endif
00082     for (G4int index= 0; index<theNumberOfSecondaries; index++){
00083       delete (*theListOfSecondaries)[index] ;
00084     }
00085   }
00086   delete theListOfSecondaries; 
00087 }

G4VParticleChange::G4VParticleChange ( const G4VParticleChange right  )  [protected]

Definition at line 89 of file G4VParticleChange.cc.

References G4FastVector< Type, N >::SetElement(), theListOfSecondaries, and theNumberOfSecondaries.


Member Function Documentation

void G4VParticleChange::AddSecondary ( G4Track aSecondary  ) 

Reimplemented in G4ParticleChange, and G4ParticleChangeForRadDecay.

Definition at line 171 of file G4VParticleChange.cc.

References CheckSecondary(), debugFlag, fSetSecondaryWeightByProcess, G4cout, G4endl, G4Exception(), JustWarning, G4FastVector< Type, N >::SetElement(), theListOfSecondaries, theNumberOfSecondaries, theParentWeight, theSizeOftheListOfSecondaries, and verboseLevel.

Referenced by G4ParticleChangeForLoss::AddSecondary(), G4ParticleChangeForGamma::AddSecondary(), G4ParticleChange::AddSecondary(), G4FastStep::CreateSecondaryTrack(), G4UnknownDecay::DecayIt(), G4Decay::DecayIt(), G4VEnergyLossProcess::PostStepDoIt(), and G4VEmProcess::PostStepDoIt().

00172 {
00173   if (debugFlag) CheckSecondary(*aTrack);
00174 
00175   // add a secondary after size check
00176   if (theSizeOftheListOfSecondaries > theNumberOfSecondaries) {
00177     // Set weight of secondary tracks
00178     if (!fSetSecondaryWeightByProcess) aTrack->SetWeight(theParentWeight);
00179     theListOfSecondaries->SetElement(theNumberOfSecondaries, aTrack);
00180     theNumberOfSecondaries++;
00181   } else {
00182     delete aTrack;
00183 
00184 #ifdef G4VERBOSE
00185     if (verboseLevel>0) {
00186       G4cout << "G4VParticleChange::AddSecondary() Warning  ";
00187       G4cout << "theListOfSecondaries is full !! " << G4endl;
00188       G4cout << " The track is deleted " << G4endl;
00189     }
00190 #endif
00191     G4Exception("G4VParticleChange::AddSecondary",
00192                 "TRACK101", JustWarning,
00193                 "Secondary Bug is full. The track is deleted"); 
00194   }
00195 }

G4bool G4VParticleChange::CheckIt ( const G4Track  )  [virtual]

Reimplemented in G4FastStep, G4ParticleChange, G4ParticleChangeForDecay, G4ParticleChangeForGamma, G4ParticleChangeForLoss, and G4ParticleChangeForMSC.

Definition at line 316 of file G4VParticleChange.cc.

References accuracyForException, accuracyForWarning, DumpInfo(), EventMustBeAborted, G4cout, G4endl, G4Exception(), theLocalEnergyDeposit, and theTrueStepLength.

Referenced by G4ParticleChangeForMSC::CheckIt(), G4ParticleChangeForLoss::CheckIt(), G4ParticleChangeForGamma::CheckIt(), G4ParticleChangeForDecay::CheckIt(), G4ParticleChange::CheckIt(), and G4FastStep::CheckIt().

00321 {
00322 
00323   G4bool    exitWithError = false;
00324   G4double  accuracy;
00325   static G4int nError = 0;
00326 #ifdef G4VERBOSE
00327   const  G4int maxError = 30;
00328 #endif
00329 
00330   // Energy deposit should not be negative
00331   G4bool itsOKforEnergy = true;
00332   accuracy = -1.0*theLocalEnergyDeposit/MeV;
00333   if (accuracy > accuracyForWarning) {
00334     itsOKforEnergy = false;
00335     nError += 1;
00336     exitWithError =  (accuracy > accuracyForException);
00337 #ifdef G4VERBOSE
00338     if (nError < maxError) {
00339       G4cout << "  G4VParticleChange::CheckIt    : ";
00340       G4cout << "the energy deposit  is negative  !!" 
00341              << "  Difference:  " << accuracy  << "[MeV] " <<G4endl;
00342       G4cout << aTrack.GetDefinition()->GetParticleName()
00343              << " E=" << aTrack.GetKineticEnergy()/MeV
00344              << " pos=" << aTrack.GetPosition().x()/m
00345              << ", " << aTrack.GetPosition().y()/m
00346              << ", " << aTrack.GetPosition().z()/m
00347              <<G4endl;
00348     }
00349 #endif
00350   }
00351  
00352   // true path length should not be negative
00353   G4bool itsOKforStepLength = true;
00354   accuracy = -1.0*theTrueStepLength/mm;
00355   if (accuracy > accuracyForWarning) {
00356     itsOKforStepLength = false;
00357     nError += 1;
00358     exitWithError =  (accuracy > accuracyForException);
00359 #ifdef G4VERBOSE
00360     if (nError < maxError) {
00361       G4cout << "  G4VParticleChange::CheckIt    : ";
00362       G4cout << "the true step length is negative  !!"
00363              << "  Difference:  " << accuracy  << "[MeV] " <<G4endl;
00364       G4cout << aTrack.GetDefinition()->GetParticleName()
00365              << " E=" << aTrack.GetKineticEnergy()/MeV
00366              << " pos=" << aTrack.GetPosition().x()/m
00367              << ", " << aTrack.GetPosition().y()/m
00368              << ", " << aTrack.GetPosition().z()/m
00369              <<G4endl;
00370     }
00371 #endif
00372 
00373   }
00374 #ifdef G4VERBOSE
00375   if (!itsOKforStepLength || !itsOKforEnergy ){
00376   // dump out information of this particle change
00377     DumpInfo();
00378   }
00379 #endif
00380 
00381   // Exit with error
00382   if (exitWithError) {
00383     G4Exception("G4VParticleChange::CheckIt",
00384                 "TRACK001", EventMustBeAborted,
00385                 "Step length and/or energy deposit was illegal");
00386   }
00387 
00388   // correction 
00389   if ( !itsOKforStepLength ) {
00390     theTrueStepLength =  (1.e-12)*mm;
00391   } 
00392   if ( !itsOKforEnergy ) {
00393     theLocalEnergyDeposit = 0.0;
00394   }
00395   return (itsOKforStepLength && itsOKforEnergy );
00396 }

G4bool G4VParticleChange::CheckSecondary ( G4Track  )  [protected]

Definition at line 398 of file G4VParticleChange.cc.

References accuracyForException, accuracyForWarning, EventMustBeAborted, G4cout, G4endl, G4Exception(), G4Track::GetDefinition(), G4Track::GetGlobalTime(), G4Track::GetKineticEnergy(), G4Track::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), ns, G4Track::SetKineticEnergy(), G4Track::SetMomentumDirection(), and theParentGlobalTime.

Referenced by AddSecondary().

00399 {
00400   G4bool    exitWithError = false;
00401   G4double  accuracy;
00402   static G4int nError = 0;
00403 #ifdef G4VERBOSE
00404   const  G4int maxError = 30;
00405 #endif
00406 
00407   // MomentumDirection should be unit vector
00408   G4bool itsOKforMomentum = true;  
00409   if (aTrack.GetKineticEnergy()>0.){
00410     accuracy = std::fabs((aTrack.GetMomentumDirection()).mag2()-1.0);
00411     if (accuracy > accuracyForWarning) {
00412       itsOKforMomentum = false;
00413       nError += 1;
00414       exitWithError = exitWithError || (accuracy > accuracyForException);
00415 #ifdef G4VERBOSE
00416       if (nError < maxError) {
00417         G4cout << " G4VParticleChange::CheckSecondary  :   ";
00418         G4cout << "the Momentum direction is not unit vector !! " 
00419                << "  Difference:  " << accuracy << G4endl;
00420         G4cout << aTrack.GetDefinition()->GetParticleName()
00421                << " E=" << aTrack.GetKineticEnergy()/MeV
00422                << " pos=" << aTrack.GetPosition().x()/m
00423                << ", " << aTrack.GetPosition().y()/m
00424                << ", " << aTrack.GetPosition().z()/m
00425                <<G4endl;
00426       }
00427 #endif
00428     }
00429   }
00430   
00431   // Kinetic Energy should not be negative
00432   G4bool itsOKforEnergy = true;
00433   accuracy = -1.0*(aTrack.GetKineticEnergy())/MeV;
00434   if (accuracy > accuracyForWarning) {
00435     itsOKforEnergy = false;
00436     nError += 1;
00437     exitWithError = exitWithError ||  (accuracy > accuracyForException);
00438 #ifdef G4VERBOSE
00439     if (nError < maxError) {
00440       G4cout << " G4VParticleChange::CheckSecondary  :   ";
00441       G4cout << "the kinetic energy is negative  !!" 
00442              << "  Difference:  " << accuracy  << "[MeV] " <<G4endl;
00443       G4cout << " G4VParticleChange::CheckSecondary  :   ";
00444       G4cout << "the global time of secondary is earlier than the parent  !!" 
00445              << "  Difference:  " << accuracy  << "[ns] " <<G4endl;
00446       G4cout << aTrack.GetDefinition()->GetParticleName()
00447              << " E=" << aTrack.GetKineticEnergy()/MeV
00448              << " pos=" << aTrack.GetPosition().x()/m
00449              << ", " << aTrack.GetPosition().y()/m
00450              << ", " << aTrack.GetPosition().z()/m
00451            <<G4endl;
00452     }
00453 #endif
00454   }
00455   // Check timing of secondaries
00456   G4bool itsOKforTiming = true;
00457 
00458   accuracy = (theParentGlobalTime - aTrack.GetGlobalTime())/ns;
00459   if (accuracy > accuracyForWarning){
00460     itsOKforTiming = false;
00461     nError += 1;
00462     exitWithError = (accuracy > accuracyForException);
00463 #ifdef G4VERBOSE
00464     if (nError < maxError) {
00465       G4cout << " G4VParticleChange::CheckSecondary  :   ";
00466       G4cout << "the global time of secondary goes back comapared to the parent  !!" 
00467              << "  Difference:  " << accuracy  << "[ns] " <<G4endl;
00468       G4cout << aTrack.GetDefinition()->GetParticleName()
00469              << " E=" << aTrack.GetKineticEnergy()/MeV
00470              << " pos=" << aTrack.GetPosition().x()/m
00471              << ", " << aTrack.GetPosition().y()/m
00472                << ", " << aTrack.GetPosition().z()/m
00473              << " time=" << aTrack.GetGlobalTime()/ns
00474              << " parent time=" << theParentGlobalTime/ns
00475              <<G4endl;
00476     }
00477 #endif
00478   }
00479 
00480   // Exit with error
00481   if (exitWithError) {
00482     G4Exception("G4VParticleChange::CheckSecondary",
00483                 "TRACK001", EventMustBeAborted,
00484                 "Secondary with illegal energy/momentum ");
00485   }
00486 
00487   G4bool itsOK = itsOKforMomentum && itsOKforEnergy && itsOKforTiming;
00488 
00489   //correction
00490   if (!itsOKforMomentum) {
00491     G4double vmag = (aTrack.GetMomentumDirection()).mag();
00492     aTrack.SetMomentumDirection((1./vmag)*aTrack.GetMomentumDirection());
00493   }
00494   if (!itsOKforEnergy) {
00495     aTrack.SetKineticEnergy(0.0);
00496   }
00497 
00498   return itsOK;
00499 }

void G4VParticleChange::Clear (  )  [inline]

Definition at line 180 of file G4VParticleChange.icc.

References theFirstStepInVolume, theLastStepInVolume, and theNumberOfSecondaries.

Referenced by G4QCaptureAtRest::AtRestDoIt(), G4ITStepProcessor::InvokeAlongStepDoItProcs(), G4ITStepProcessor::InvokeAtRestDoItProcs(), G4ITStepProcessor::InvokePSDIP(), G4WHadronElasticProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), and G4HadronElasticProcess::PostStepDoIt().

00181 {
00182   theNumberOfSecondaries = 0;
00183   theFirstStepInVolume = false;
00184   theLastStepInVolume = false;
00185 }

void G4VParticleChange::ClearDebugFlag (  )  [inline]

Definition at line 290 of file G4VParticleChange.icc.

References debugFlag.

Referenced by G4ErrorEnergyLoss::AlongStepDoIt().

00291 {
00292   debugFlag = false;
00293 }

void G4VParticleChange::DumpInfo (  )  const [virtual]

Reimplemented in G4FastStep, G4ParticleChange, G4ParticleChangeForDecay, G4ParticleChangeForGamma, G4ParticleChangeForLoss, G4ParticleChangeForMSC, and G4ParticleChangeForTransport.

Definition at line 252 of file G4VParticleChange.cc.

References fAlive, fKillTrackAndSecondaries, fPostponeToNextEvent, fStopAndKill, fStopButAlive, fSuspend, G4cout, G4endl, GetSecondary(), theFirstStepInVolume, theLastStepInVolume, theLocalEnergyDeposit, theNonIonizingEnergyDeposit, theNumberOfSecondaries, theStatusChange, theSteppingControlFlag, and theTrueStepLength.

Referenced by CheckIt(), G4ParticleChangeForMSC::DumpInfo(), G4ParticleChangeForLoss::DumpInfo(), G4ParticleChangeForGamma::DumpInfo(), G4ParticleChangeForDecay::DumpInfo(), G4ParticleChange::DumpInfo(), G4FastStep::DumpInfo(), and G4SteppingVerbose::VerboseParticleChange().

00253 {
00254 
00255 // Show header
00256   G4int olprc = G4cout.precision(3);
00257   G4cout << "      -----------------------------------------------" 
00258        << G4endl;
00259   G4cout << "        G4ParticleChange Information  " << std::setw(20) << G4endl;
00260   G4cout << "      -----------------------------------------------" 
00261        << G4endl;
00262 
00263   G4cout << "        # of 2ndaries       : " 
00264        << std::setw(20) << theNumberOfSecondaries
00265        << G4endl;
00266 
00267   if (theNumberOfSecondaries >0) {
00268     G4cout << "        Pointer to 2ndaries : " 
00269          << std::setw(20) << GetSecondary(0)
00270          << G4endl;
00271     G4cout << "        (Showed only 1st one)"
00272          << G4endl;
00273   }
00274   G4cout << "      -----------------------------------------------" 
00275        << G4endl;
00276 
00277   G4cout << "        Energy Deposit (MeV): " 
00278        << std::setw(20) << theLocalEnergyDeposit/MeV
00279        << G4endl;
00280 
00281   G4cout << "        Non-ionizing Energy Deposit (MeV): " 
00282        << std::setw(20) << theNonIonizingEnergyDeposit/MeV
00283        << G4endl;
00284 
00285   G4cout << "        Track Status        : " 
00286        << std::setw(20);
00287        if( theStatusChange == fAlive ){
00288          G4cout << " Alive";
00289        } else if( theStatusChange == fStopButAlive ){
00290            G4cout << " StopButAlive";
00291        } else if( theStatusChange == fStopAndKill ){
00292            G4cout << " StopAndKill";
00293        } else if( theStatusChange  == fKillTrackAndSecondaries ){
00294            G4cout << " KillTrackAndSecondaries";
00295        } else if( theStatusChange  == fSuspend ){
00296            G4cout << " Suspend";
00297        } else if( theStatusChange == fPostponeToNextEvent ){
00298            G4cout << " PostponeToNextEvent";
00299        }
00300        G4cout << G4endl;
00301   G4cout << "        True Path Length (mm) : " 
00302        << std::setw(20) << theTrueStepLength/mm
00303        << G4endl;
00304   G4cout << "        Stepping Control      : " 
00305        << std::setw(20) << theSteppingControlFlag
00306        << G4endl;   
00307   if (theFirstStepInVolume) {
00308     G4cout << "    First Step In the voulme  : "  << G4endl;
00309   }
00310   if (theLastStepInVolume) {
00311     G4cout << "    Last Step In the voulme  : "  << G4endl;
00312   }
00313   G4cout.precision(olprc);
00314 }

G4double G4VParticleChange::GetAccuracyForException (  )  const [protected]

Definition at line 507 of file G4VParticleChange.cc.

References accuracyForException.

Referenced by G4FastStep::CheckIt().

00508 {
00509   return accuracyForException;
00510 }

G4double G4VParticleChange::GetAccuracyForWarning (  )  const [protected]

Definition at line 502 of file G4VParticleChange.cc.

References accuracyForWarning.

Referenced by G4FastStep::CheckIt().

00503 {
00504   return accuracyForWarning;
00505 }

G4bool G4VParticleChange::GetDebugFlag (  )  const [inline]

Definition at line 302 of file G4VParticleChange.icc.

References debugFlag.

00303 {
00304   return debugFlag;
00305 }

G4bool G4VParticleChange::GetFirstStepInVolume (  )  const [inline]

Definition at line 71 of file G4VParticleChange.icc.

References theFirstStepInVolume.

00072 {
00073    return theFirstStepInVolume;
00074 }

G4bool G4VParticleChange::GetLastStepInVolume (  )  const [inline]

Definition at line 77 of file G4VParticleChange.icc.

References theLastStepInVolume.

00078 {
00079   return theLastStepInVolume;
00080 }

G4double G4VParticleChange::GetLocalEnergyDeposit (  )  const [inline]

Definition at line 99 of file G4VParticleChange.icc.

References theLocalEnergyDeposit.

Referenced by G4FastStep::GetTotalEnergyDeposited(), G4VEnergyLossProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4PolarizedComptonModel::SampleSecondaries(), and G4HeatedKleinNishinaCompton::SampleSecondaries().

00100 {
00101   return theLocalEnergyDeposit;
00102 }

G4double G4VParticleChange::GetNonIonizingEnergyDeposit (  )  const [inline]

Definition at line 111 of file G4VParticleChange.icc.

References theNonIonizingEnergyDeposit.

00112 {
00113   return theNonIonizingEnergyDeposit;
00114 }

G4int G4VParticleChange::GetNumberOfSecondaries (  )  const [inline]

Definition at line 41 of file G4VParticleChange.icc.

References theNumberOfSecondaries.

Referenced by G4HadronicProcess::CheckEnergyMomentumConservation(), G4ITStepProcessor::DealWithSecondaries(), G4FastStep::GetNumberOfSecondaryTracks(), G4Scintillation::PostStepDoIt(), G4OpWLS::PostStepDoIt(), and G4Cerenkov::PostStepDoIt().

00042 {
00043   return theNumberOfSecondaries;
00044 }

G4double G4VParticleChange::GetParentWeight (  )  const [inline]

Definition at line 148 of file G4VParticleChange.icc.

References theParentWeight.

Referenced by G4VEnergyLossProcess::AlongStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), and G4VEmProcess::PostStepDoIt().

00149 {
00150    return theParentWeight;
00151 }

G4Track * G4VParticleChange::GetSecondary ( G4int  anIndex  )  const [inline]

Definition at line 35 of file G4VParticleChange.icc.

References theListOfSecondaries.

Referenced by G4QCaptureAtRest::AtRestDoIt(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4ITStepProcessor::DealWithSecondaries(), DumpInfo(), and G4FastStep::GetSecondaryTrack().

00036 {
00037   return (*theListOfSecondaries)[anIndex];
00038 }

G4SteppingControl G4VParticleChange::GetSteppingControl (  )  const [inline]

Definition at line 59 of file G4VParticleChange.icc.

References theSteppingControlFlag.

00060 {
00061   return theSteppingControlFlag;    
00062 }

G4TrackStatus G4VParticleChange::GetTrackStatus (  )  const [inline]

Definition at line 53 of file G4VParticleChange.icc.

References theStatusChange.

Referenced by G4HadronicProcess::CheckEnergyMomentumConservation(), G4ParticleChange::CheckIt(), G4ITStepProcessor::InvokeAlongStepDoItProcs(), G4ITStepProcessor::InvokePSDIP(), G4VEnergyLossProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4QInelastic::PostStepDoIt(), and G4FastSimulationManagerProcess::PostStepDoIt().

00054 {
00055   return theStatusChange;
00056 }

G4double G4VParticleChange::GetTrueStepLength (  )  const [inline]

Definition at line 123 of file G4VParticleChange.icc.

References theTrueStepLength.

00124 {
00125   return theTrueStepLength;
00126 }

G4int G4VParticleChange::GetVerboseLevel (  )  const [inline]

Definition at line 142 of file G4VParticleChange.icc.

References verboseLevel.

00143 { 
00144   return verboseLevel; 
00145 }

G4double G4VParticleChange::GetWeight (  )  const [inline]

Definition at line 154 of file G4VParticleChange.icc.

References theParentWeight.

Referenced by G4AdjointhIonisationModel::SampleSecondaries(), G4AdjointComptonModel::SampleSecondaries(), and G4AdjointBremsstrahlungModel::SampleSecondaries().

00155 {
00156    return theParentWeight;
00157 }

void G4VParticleChange::Initialize ( const G4Track  )  [inline, virtual]

Reimplemented in G4ParticleChange, G4ParticleChangeForDecay, G4ParticleChangeForMSC, and G4ParticleChangeForTransport.

Definition at line 277 of file G4VParticleChange.icc.

References InitializeLocalEnergyDeposit(), InitializeParentGlobalTime(), InitializeParentWeight(), InitializeSecondaries(), InitializeStatusChange(), InitializeStepInVolumeFlags(), InitializeSteppingControl(), and InitializeTrueStepLength().

Referenced by G4WeightWindowProcess::AlongStepDoIt(), G4WeightCutOffProcess::AlongStepDoIt(), G4ScoreSplittingProcess::AlongStepDoIt(), G4ParallelWorldScoringProcess::AlongStepDoIt(), G4ParallelWorldProcess::AlongStepDoIt(), G4ImportanceProcess::AlongStepDoIt(), G4FastSimulationManagerProcess::AlongStepDoIt(), G4ScoreSplittingProcess::AtRestDoIt(), G4ParallelWorldScoringProcess::AtRestDoIt(), G4ParallelWorldProcess::AtRestDoIt(), G4ParticleChangeForDecay::Initialize(), G4ParticleChange::Initialize(), G4FastStep::Initialize(), G4VErrorLimitProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4ParallelWorldScoringProcess::PostStepDoIt(), G4ParallelWorldProcess::PostStepDoIt(), G4NeutronKiller::PostStepDoIt(), and G4ErrorTrackLengthTarget::PostStepDoIt().

00278 {
00279   InitializeStatusChange(track);
00280   InitializeLocalEnergyDeposit(track);
00281   InitializeSteppingControl(track);
00282   InitializeTrueStepLength(track);
00283   InitializeSecondaries(track);
00284   InitializeParentWeight(track);
00285   InitializeParentGlobalTime(track);
00286   InitializeStepInVolumeFlags(track);
00287 }

void G4VParticleChange::InitializeLocalEnergyDeposit ( const G4Track  )  [inline, protected]

Definition at line 165 of file G4VParticleChange.icc.

References theLocalEnergyDeposit, and theNonIonizingEnergyDeposit.

Referenced by Initialize().

00166 {  
00167   // clear theLocalEnergyDeposited   
00168   theLocalEnergyDeposit = 0.0;
00169   theNonIonizingEnergyDeposit = 0.0;
00170 }

void G4VParticleChange::InitializeParentGlobalTime ( const G4Track  )  [inline, protected]

Definition at line 204 of file G4VParticleChange.icc.

References G4StepPoint::GetGlobalTime(), G4Step::GetPreStepPoint(), G4Track::GetStep(), and theParentGlobalTime.

Referenced by Initialize().

00205 {
00206   // set the parent track's global time at the pre-step point
00207   theParentGlobalTime = track.GetStep()->GetPreStepPoint()->GetGlobalTime();
00208 } 

void G4VParticleChange::InitializeParentWeight ( const G4Track  )  [inline, protected]

Definition at line 197 of file G4VParticleChange.icc.

References G4Track::GetWeight(), isParentWeightProposed, and theParentWeight.

Referenced by Initialize().

00198 {
00199   // set the parent track's weight
00200   theParentWeight = track.GetWeight();
00201   isParentWeightProposed = false;
00202 }

void G4VParticleChange::InitializeSecondaries ( const G4Track  )  [inline, protected]

Definition at line 230 of file G4VParticleChange.icc.

References G4cerr, G4endl, theListOfSecondaries, theNumberOfSecondaries, and verboseLevel.

Referenced by Initialize(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForLoss::InitializeForPostStep(), and G4ParticleChangeForGamma::InitializeForPostStep().

00231 {
00232   // clear secondaries
00233   if (theNumberOfSecondaries>0) {
00234 #ifdef G4VERBOSE
00235     if (verboseLevel>0) {
00236       G4cerr << "G4VParticleChange::Initialize() Warning  ";
00237       G4cerr << "theListOfSecondaries is not empty " << G4endl;
00238       G4cerr << "All objects in theListOfSecondaries are destroyed!" << G4endl;
00239     }
00240 #endif
00241     for (G4int index= 0; index<theNumberOfSecondaries; index++){
00242       if ( (*theListOfSecondaries)[index] ){ 
00243          delete (*theListOfSecondaries)[index] ;
00244        }
00245     }
00246   }
00247   theNumberOfSecondaries = 0;
00248 }

void G4VParticleChange::InitializeStatusChange ( const G4Track  )  [inline, protected]

Definition at line 191 of file G4VParticleChange.icc.

References G4Track::GetTrackStatus(), and theStatusChange.

Referenced by Initialize(), and G4ParticleChangeForTransport::Initialize().

00192 {
00193   // set TrackStatus equal to the parent track's one
00194   theStatusChange = track.GetTrackStatus();
00195 }

void G4VParticleChange::InitializeStepInVolumeFlags ( const G4Track  )  [inline, protected]

Definition at line 223 of file G4VParticleChange.icc.

References G4Track::GetStep(), theFirstStepInVolume, and theLastStepInVolume.

Referenced by Initialize().

00224 {
00225    const G4Step* aStep = track.GetStep();
00226    theFirstStepInVolume = aStep-> IsFirstStepInVolume();
00227    theLastStepInVolume  = aStep-> IsLastStepInVolume();
00228 }

void G4VParticleChange::InitializeSteppingControl ( const G4Track  )  [inline, protected]

Definition at line 173 of file G4VParticleChange.icc.

References NormalCondition, and theSteppingControlFlag.

Referenced by Initialize(), and G4ParticleChangeForTransport::Initialize().

00174 {
00175   // SteppingControlFlag
00176   theSteppingControlFlag = NormalCondition;     
00177 }

void G4VParticleChange::InitializeTrueStepLength ( const G4Track  )  [inline, protected]

Definition at line 211 of file G4VParticleChange.icc.

References G4Track::GetStep(), G4Step::GetStepLength(), and theTrueStepLength.

Referenced by Initialize().

00212 {
00213   // Reset theTrueStepLength
00214   // !! Caution  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00215   theTrueStepLength = track.GetStep()->GetStepLength();
00216   // !!  TrueStepLength should be copied from G4Step not G4Track
00217   // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00218 }

G4bool G4VParticleChange::IsParentWeightSetByProcess (  )  const

Definition at line 521 of file G4VParticleChange.cc.

00522 {
00523   return  true;
00524 }

G4bool G4VParticleChange::IsSecondaryWeightSetByProcess (  )  const [inline]

Definition at line 314 of file G4VParticleChange.icc.

References fSetSecondaryWeightByProcess.

00315 {
00316   return fSetSecondaryWeightByProcess;
00317 }

G4bool G4VParticleChange::operator!= ( const G4VParticleChange right  )  const

Definition at line 166 of file G4VParticleChange.cc.

00167 {
00168    return (this != (G4VParticleChange *) &right);
00169 }

G4VParticleChange & G4VParticleChange::operator= ( const G4VParticleChange right  )  [protected]

Definition at line 116 of file G4VParticleChange.cc.

References debugFlag, fSetSecondaryWeightByProcess, G4cout, isParentWeightProposed, G4FastVector< Type, N >::SetElement(), theFirstStepInVolume, theLastStepInVolume, theListOfSecondaries, theLocalEnergyDeposit, theNonIonizingEnergyDeposit, theNumberOfSecondaries, theParentGlobalTime, theParentWeight, theStatusChange, theSteppingControlFlag, theTrueStepLength, and verboseLevel.

Referenced by G4FastStep::operator=().

00117 {
00118   if (this != &right){
00119     if (theNumberOfSecondaries>0) {
00120 #ifdef G4VERBOSE
00121       if (verboseLevel>0) {
00122         G4cout << "G4VParticleChange: assignment operator Warning  ";
00123         G4cout << "theListOfSecondaries is not empty ";
00124       }
00125 #endif
00126       for (G4int index = 0; index<theNumberOfSecondaries; index++){
00127         if ( (*theListOfSecondaries)[index] ) delete ((*theListOfSecondaries)[index]) ;
00128       }
00129     }
00130     delete theListOfSecondaries; 
00131       
00132     theListOfSecondaries =  new G4TrackFastVector();
00133     theNumberOfSecondaries = right.theNumberOfSecondaries;
00134    for (G4int index = 0; index<theNumberOfSecondaries; index++){
00135     G4Track* newTrack =  new G4Track(*((*right.theListOfSecondaries)[index] ));
00136     theListOfSecondaries->SetElement(index, newTrack);                      
00137   }
00138     theStatusChange = right.theStatusChange;
00139     theSteppingControlFlag = right.theSteppingControlFlag;
00140     theLocalEnergyDeposit = right.theLocalEnergyDeposit;
00141     theNonIonizingEnergyDeposit = right.theNonIonizingEnergyDeposit;
00142     theTrueStepLength = right.theTrueStepLength;
00143     
00144     theFirstStepInVolume = right.theFirstStepInVolume;
00145     theLastStepInVolume =  right.theLastStepInVolume;
00146 
00147     theParentWeight = right.theParentWeight;
00148     isParentWeightProposed = right.isParentWeightProposed;
00149     fSetSecondaryWeightByProcess = right.fSetSecondaryWeightByProcess;
00150 
00151     theParentGlobalTime = right.theParentGlobalTime;
00152 
00153     verboseLevel = right.verboseLevel;
00154     debugFlag = right.debugFlag;
00155     
00156   }
00157   return *this;
00158 }

G4bool G4VParticleChange::operator== ( const G4VParticleChange right  )  const

Definition at line 160 of file G4VParticleChange.cc.

00161 {
00162    return (this == (G4VParticleChange *) &right);
00163 }

void G4VParticleChange::ProposeFirstStepInVolume ( G4bool  flag  )  [inline]

Definition at line 83 of file G4VParticleChange.icc.

References theFirstStepInVolume.

00084 {
00085     theFirstStepInVolume = flag;
00086 }

void G4VParticleChange::ProposeLastStepInVolume ( G4bool  flag  )  [inline]

Definition at line 89 of file G4VParticleChange.icc.

References theLastStepInVolume.

Referenced by G4Transportation::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), and G4CoupledTransportation::PostStepDoIt().

00090 {
00091      theLastStepInVolume = flag;
00092 }

void G4VParticleChange::ProposeLocalEnergyDeposit ( G4double  anEnergyPart  )  [inline]

Definition at line 105 of file G4VParticleChange.icc.

References theLocalEnergyDeposit.

Referenced by G4VEnergyLossProcess::AlongStepDoIt(), G4NuclearStopping::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4QCaptureAtRest::AtRestDoIt(), G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4MuonMinusCaptureAtRest::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4UnknownDecay::DecayIt(), G4RadioactiveDecay::DecayIt(), G4DNAMolecularDecay::DecayIt(), G4Decay::DecayIt(), G4HadronicProcess::FillResult(), SpecialCuts::PostStepDoIt(), G4WHadronElasticProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4QSynchRad::PostStepDoIt(), G4QNGamma::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QIonIonElastic::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QElastic::PostStepDoIt(), G4QDiffraction::PostStepDoIt(), G4QCoherentChargeExchange::PostStepDoIt(), G4QAtomicElectronScattering::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4FastStep::ProposeTotalEnergyDeposited(), G4PolarizedComptonModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeIonisationModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeBremsstrahlungModel::SampleSecondaries(), G4PEEffectModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4MuElecInelasticModel::SampleSecondaries(), G4MuElecElasticModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LivermorePolarizedPhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermoreIonisationModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4HeatedKleinNishinaCompton::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4DNATransformElectronModel::SampleSecondaries(), G4DNAScreenedRutherfordElasticModel::SampleSecondaries(), G4DNASancheSolvatationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNAMillerGreenExcitationModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNAEmfietzoglouExcitationModel::SampleSecondaries(), G4DNADingfelderChargeIncreaseModel::SampleSecondaries(), G4DNADingfelderChargeDecreaseModel::SampleSecondaries(), G4DNAChampionElasticModel::SampleSecondaries(), G4DNABornIonisationModel::SampleSecondaries(), G4DNABornExcitationModel::SampleSecondaries(), and G4BoldyshevTripletModel::SampleSecondaries().

00106 {
00107   theLocalEnergyDeposit = anEnergyPart;
00108 }

void G4VParticleChange::ProposeNonIonizingEnergyDeposit ( G4double  anEnergyPart  )  [inline]

Definition at line 117 of file G4VParticleChange.icc.

References theNonIonizingEnergyDeposit.

Referenced by G4NuclearStopping::AlongStepDoIt(), G4WHadronElasticProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4IonCoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), and G4eCoulombScatteringModel::SampleSecondaries().

00118 {
00119   theNonIonizingEnergyDeposit = anEnergyPart;
00120 }

void G4VParticleChange::ProposeParentWeight ( G4double  finalWeight  )  [inline]

Definition at line 327 of file G4VParticleChange.icc.

References ProposeWeight().

Referenced by G4ContinuousGainOfEnergy::AlongStepDoIt(), G4AdjointAlongStepWeightCorrection::AlongStepDoIt(), G4VEmAdjointModel::CorrectPostStepWeight(), G4AdjointPhotoElectricModel::CorrectPostStepWeight(), G4AdjointIonIonisationModel::CorrectPostStepWeight(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), and G4AdjointeIonisationModel::SampleSecondaries().

00328 {
00329   ProposeWeight(w);
00330 }

void G4VParticleChange::ProposeSteppingControl ( G4SteppingControl  StepControlFlag  )  [inline]

Definition at line 65 of file G4VParticleChange.icc.

References theSteppingControlFlag.

Referenced by G4FastStep::ForceSteppingHitInvocation(), and G4ScoreSplittingProcess::PostStepDoIt().

00066 {
00067   theSteppingControlFlag = StepControlFlag;
00068 }

void G4VParticleChange::ProposeTrackStatus ( G4TrackStatus  status  )  [inline]

Definition at line 47 of file G4VParticleChange.icc.

References theStatusChange.

Referenced by G4Transportation::AlongStepDoIt(), G4ITTransportation::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4CoupledTransportation::AlongStepDoIt(), G4QCaptureAtRest::AtRestDoIt(), G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4MuonMinusCaptureAtRest::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4eplusPolarizedAnnihilation::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4UnknownDecay::DecayIt(), G4RadioactiveDecay::DecayIt(), G4DNAMolecularDecay::DecayIt(), G4Decay::DecayIt(), G4DNABrownianTransportation::Diffusion(), G4HadronicProcess::FillResult(), G4FastStep::KillPrimaryTrack(), G4WeightWindowProcess::KillTrack(), G4ImportanceProcess::KillTrack(), SpecialCuts::PostStepDoIt(), G4WHadronElasticProcess::PostStepDoIt(), G4WeightCutOffProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4QSynchRad::PostStepDoIt(), G4QNGamma::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QIonIonElastic::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QElastic::PostStepDoIt(), G4QDiffraction::PostStepDoIt(), G4QCoherentChargeExchange::PostStepDoIt(), G4QAtomicElectronScattering::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4NeutronKiller::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4FastSimulationManagerProcess::PostStepDoIt(), G4DNASecondOrderReaction::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4SeltzerBergerModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeAnnihilationModel::SampleSecondaries(), G4PEEffectModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4MuElecElasticModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LivermorePolarizedPhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermoreGammaConversionModelRC::SampleSecondaries(), G4LivermoreGammaConversionModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4HeatedKleinNishinaCompton::SampleSecondaries(), G4eeToHadronsMultiModel::SampleSecondaries(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4eBremsstrahlungModel::SampleSecondaries(), G4eBremParametrizedModel::SampleSecondaries(), G4DNATransformElectronModel::SampleSecondaries(), G4DNAScreenedRutherfordElasticModel::SampleSecondaries(), G4DNASancheSolvatationModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNADingfelderChargeIncreaseModel::SampleSecondaries(), G4DNADingfelderChargeDecreaseModel::SampleSecondaries(), G4DNAChampionElasticModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4AdjointIonIonisationModel::SampleSecondaries(), G4AdjointhIonisationModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), G4AdjointComptonModel::SampleSecondaries(), and G4AdjointBremsstrahlungModel::SampleSecondaries().

00048 {
00049   theStatusChange = aStatus;
00050 }

void G4VParticleChange::ProposeTrueStepLength ( G4double  truePathLength  )  [inline]

Definition at line 129 of file G4VParticleChange.icc.

References theTrueStepLength.

Referenced by G4VMultipleScattering::AlongStepDoIt(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), and G4FastStep::ProposePrimaryTrackPathLength().

00130 {
00131   theTrueStepLength = aLength;
00132 }

void G4VParticleChange::ProposeWeight ( G4double  finalWeight  )  [inline]

Definition at line 320 of file G4VParticleChange.icc.

References isParentWeightProposed, and theParentWeight.

Referenced by G4VEnergyLossProcess::AlongStepDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4SamplingPostStepAction::DoIt(), G4WHadronElasticProcess::PostStepDoIt(), G4WeightCutOffProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), and ProposeParentWeight().

00321 {
00322   theParentWeight = w;
00323   isParentWeightProposed = true;
00324 }

void G4VParticleChange::SetDebugFlag (  )  [inline]

Definition at line 296 of file G4VParticleChange.icc.

References debugFlag.

00297 {
00298   debugFlag = true;
00299 }

void G4VParticleChange::SetNumberOfSecondaries ( G4int  totSecondaries  )  [inline]

Definition at line 254 of file G4VParticleChange.icc.

References G4cerr, G4FastVector< Type, N >::Initialize(), theListOfSecondaries, theNumberOfSecondaries, theSizeOftheListOfSecondaries, and verboseLevel.

Referenced by G4ErrorEnergyLoss::AlongStepDoIt(), G4QCaptureAtRest::AtRestDoIt(), G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4MuonMinusCaptureAtRest::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4eplusPolarizedAnnihilation::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4UnknownDecay::DecayIt(), G4RadioactiveDecay::DecayIt(), G4DNAMolecularDecay::DecayIt(), G4Decay::DecayIt(), G4HadronicProcess::FillResult(), G4WHadronElasticProcess::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4QSynchRad::PostStepDoIt(), G4QNGamma::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QAtomicElectronScattering::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), and G4FastStep::SetNumberOfSecondaryTracks().

00255 {
00256   // check if tracks still exist in theListOfSecondaries
00257   if (theNumberOfSecondaries>0) {
00258 #ifdef G4VERBOSE
00259     if (verboseLevel>0) {
00260       G4cerr << "G4VParticleChange::SetNumberOfSecondaries() Warning  ";
00261       G4cerr << "theListOfSecondaries is not empty ";
00262     }
00263 #endif
00264     for (G4int index= 0; index<theNumberOfSecondaries; index++){
00265       if ( (*theListOfSecondaries)[index] ){
00266         delete (*theListOfSecondaries)[index] ;
00267       }
00268     }
00269   }
00270   theNumberOfSecondaries = 0;
00271   theSizeOftheListOfSecondaries = totSecondaries;
00272 
00273   // Initialize ListOfSecondaries
00274   theListOfSecondaries->Initialize(totSecondaries);
00275 }

void G4VParticleChange::SetParentWeightByProcess ( G4bool   ) 

Definition at line 516 of file G4VParticleChange.cc.

Referenced by G4ContinuousGainOfEnergy::AlongStepDoIt(), G4AdjointAlongStepWeightCorrection::AlongStepDoIt(), G4VEmAdjointModel::CorrectPostStepWeight(), G4AdjointPhotoElectricModel::CorrectPostStepWeight(), G4AdjointIonIonisationModel::CorrectPostStepWeight(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), and G4AdjointeIonisationModel::SampleSecondaries().

00517 {
00518 }

void G4VParticleChange::SetSecondaryWeightByProcess ( G4bool   )  [inline]

Definition at line 308 of file G4VParticleChange.icc.

References fSetSecondaryWeightByProcess.

Referenced by G4AdjointAlongStepWeightCorrection::AlongStepDoIt(), G4VEmAdjointModel::CorrectPostStepWeight(), G4AdjointPhotoElectricModel::CorrectPostStepWeight(), G4AdjointIonIonisationModel::CorrectPostStepWeight(), G4HadronicProcess::G4HadronicProcess(), G4VEnergyLossProcess::G4VEnergyLossProcess(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), and G4AdjointeIonisationModel::SampleSecondaries().

00309 {
00310   fSetSecondaryWeightByProcess = flag;
00311 }

void G4VParticleChange::SetVerboseLevel ( G4int  vLevel  )  [inline]

Definition at line 136 of file G4VParticleChange.icc.

References verboseLevel.

00137 { 
00138   verboseLevel = vLevel; 
00139 }

G4Step * G4VParticleChange::UpdateStepForAlongStep ( G4Step Step  )  [virtual]

Reimplemented in G4ParticleChange, G4ParticleChangeForLoss, G4ParticleChangeForMSC, and G4ParticleChangeForTransport.

Definition at line 228 of file G4VParticleChange.cc.

References G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4StepPoint::GetWeight(), isParentWeightProposed, G4StepPoint::SetWeight(), theParentWeight, and UpdateStepInfo().

Referenced by G4ITStepProcessor::InvokeAlongStepDoItProcs().

00229 {
00230   if (isParentWeightProposed ){ 
00231     G4double initialWeight = Step->GetPreStepPoint()->GetWeight();
00232     G4double currentWeight = Step->GetPostStepPoint()->GetWeight();
00233     G4double finalWeight   = (theParentWeight/initialWeight)*currentWeight;
00234     Step->GetPostStepPoint()->SetWeight( finalWeight );
00235   }   
00236   return UpdateStepInfo(Step);
00237 }

G4Step * G4VParticleChange::UpdateStepForAtRest ( G4Step Step  )  [virtual]

Reimplemented in G4FastStep, G4ParticleChange, G4ParticleChangeForDecay, G4ParticleChangeForGamma, and G4ParticleChangeForTransport.

Definition at line 219 of file G4VParticleChange.cc.

References G4Step::GetPostStepPoint(), isParentWeightProposed, G4StepPoint::SetWeight(), theParentWeight, and UpdateStepInfo().

Referenced by G4ITStepProcessor::InvokeAtRestDoItProcs().

00220 { 
00221   if (isParentWeightProposed ){
00222     Step->GetPostStepPoint()->SetWeight( theParentWeight );
00223   }
00224   return UpdateStepInfo(Step);
00225 }

G4Step * G4VParticleChange::UpdateStepForPostStep ( G4Step Step  )  [virtual]

Reimplemented in G4FastStep, G4ParticleChange, G4ParticleChangeForDecay, G4ParticleChangeForGamma, G4ParticleChangeForLoss, G4ParticleChangeForMSC, and G4ParticleChangeForTransport.

Definition at line 239 of file G4VParticleChange.cc.

References G4Step::GetPostStepPoint(), isParentWeightProposed, G4StepPoint::SetWeight(), theParentWeight, and UpdateStepInfo().

Referenced by G4ITStepProcessor::InvokePSDIP().

00240 {
00241   if (isParentWeightProposed ){ 
00242     Step->GetPostStepPoint()->SetWeight( theParentWeight );
00243   }
00244   return UpdateStepInfo(Step);
00245 }

G4Step * G4VParticleChange::UpdateStepInfo ( G4Step Step  )  [protected]

Reimplemented in G4ParticleChange.

Definition at line 202 of file G4VParticleChange.cc.

References G4Step::AddNonIonizingEnergyDeposit(), G4Step::AddTotalEnergyDeposit(), G4Step::ClearFirstStepFlag(), G4Step::ClearLastStepFlag(), G4Step::SetControlFlag(), G4Step::SetFirstStepFlag(), G4Step::SetLastStepFlag(), G4Step::SetStepLength(), theFirstStepInVolume, theLastStepInVolume, theLocalEnergyDeposit, theNonIonizingEnergyDeposit, theSteppingControlFlag, and theTrueStepLength.

Referenced by UpdateStepForAlongStep(), UpdateStepForAtRest(), G4ParticleChangeForDecay::UpdateStepForAtRest(), G4FastStep::UpdateStepForAtRest(), UpdateStepForPostStep(), G4ParticleChangeForDecay::UpdateStepForPostStep(), G4FastStep::UpdateStepForPostStep(), and G4ParticleChange::UpdateStepInfo().

00203 {
00204   // Update the G4Step specific attributes
00205   pStep->SetStepLength( theTrueStepLength );
00206   pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
00207   pStep->AddNonIonizingEnergyDeposit( theNonIonizingEnergyDeposit );
00208   pStep->SetControlFlag( theSteppingControlFlag );
00209 
00210   if (theFirstStepInVolume) {pStep->SetFirstStepFlag();}
00211   else                      {pStep->ClearFirstStepFlag();} 
00212   if (theLastStepInVolume)  {pStep->SetLastStepFlag();}
00213   else                      {pStep->ClearLastStepFlag();} 
00214 
00215   return pStep;
00216 }


Field Documentation

const G4double G4VParticleChange::accuracyForException = 0.001 [static, protected]

Definition at line 318 of file G4VParticleChange.hh.

Referenced by CheckIt(), G4ParticleChangeForMSC::CheckIt(), G4ParticleChangeForLoss::CheckIt(), G4ParticleChangeForGamma::CheckIt(), G4ParticleChangeForDecay::CheckIt(), G4ParticleChange::CheckIt(), CheckSecondary(), and GetAccuracyForException().

const G4double G4VParticleChange::accuracyForWarning = 1.0e-9 [static, protected]

Definition at line 317 of file G4VParticleChange.hh.

Referenced by CheckIt(), G4ParticleChangeForMSC::CheckIt(), G4ParticleChangeForLoss::CheckIt(), G4ParticleChangeForGamma::CheckIt(), G4ParticleChangeForDecay::CheckIt(), G4ParticleChange::CheckIt(), CheckSecondary(), and GetAccuracyForWarning().

G4bool G4VParticleChange::debugFlag [protected]

Definition at line 314 of file G4VParticleChange.hh.

Referenced by AddSecondary(), ClearDebugFlag(), G4ParticleChangeForGamma::G4ParticleChangeForGamma(), G4ParticleChangeForLoss::G4ParticleChangeForLoss(), G4VParticleChange(), GetDebugFlag(), operator=(), SetDebugFlag(), G4ParticleChangeForTransport::UpdateStepForAlongStep(), G4ParticleChange::UpdateStepForAlongStep(), G4ParticleChangeForDecay::UpdateStepForAtRest(), G4ParticleChange::UpdateStepForAtRest(), G4FastStep::UpdateStepForAtRest(), G4ParticleChange::UpdateStepForPostStep(), and G4FastStep::UpdateStepForPostStep().

G4bool G4VParticleChange::fSetSecondaryWeightByProcess [protected]

Definition at line 286 of file G4VParticleChange.hh.

Referenced by AddSecondary(), IsSecondaryWeightSetByProcess(), operator=(), G4ParticleChangeForLoss::operator=(), and SetSecondaryWeightByProcess().

G4bool G4VParticleChange::isParentWeightProposed [protected]

Definition at line 284 of file G4VParticleChange.hh.

Referenced by G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4ParticleChangeForGamma::InitializeForPostStep(), InitializeParentWeight(), operator=(), G4ParticleChangeForLoss::operator=(), ProposeWeight(), UpdateStepForAlongStep(), G4ParticleChangeForLoss::UpdateStepForAlongStep(), G4ParticleChange::UpdateStepForAlongStep(), UpdateStepForAtRest(), G4ParticleChangeForGamma::UpdateStepForAtRest(), G4ParticleChangeForDecay::UpdateStepForAtRest(), G4ParticleChange::UpdateStepForAtRest(), UpdateStepForPostStep(), G4ParticleChangeForLoss::UpdateStepForPostStep(), G4ParticleChangeForGamma::UpdateStepForPostStep(), G4ParticleChangeForDecay::UpdateStepForPostStep(), and G4ParticleChange::UpdateStepForPostStep().

G4bool G4VParticleChange::theFirstStepInVolume [protected]

Definition at line 278 of file G4VParticleChange.hh.

Referenced by Clear(), DumpInfo(), GetFirstStepInVolume(), InitializeStepInVolumeFlags(), operator=(), ProposeFirstStepInVolume(), and UpdateStepInfo().

G4bool G4VParticleChange::theLastStepInVolume [protected]

Definition at line 279 of file G4VParticleChange.hh.

Referenced by Clear(), DumpInfo(), GetLastStepInVolume(), InitializeStepInVolumeFlags(), operator=(), ProposeLastStepInVolume(), and UpdateStepInfo().

G4TrackFastVector* G4VParticleChange::theListOfSecondaries [protected]

Definition at line 245 of file G4VParticleChange.hh.

Referenced by AddSecondary(), G4ParticleChangeForRadDecay::AddSecondary(), G4VParticleChange(), GetSecondary(), InitializeSecondaries(), operator=(), G4ParticleChangeForTransport::operator=(), G4ParticleChangeForMSC::operator=(), G4ParticleChangeForLoss::operator=(), G4ParticleChangeForGamma::operator=(), G4ParticleChangeForDecay::operator=(), G4ParticleChange::operator=(), G4FastStep::operator=(), SetNumberOfSecondaries(), and ~G4VParticleChange().

G4double G4VParticleChange::theLocalEnergyDeposit [protected]

Definition at line 260 of file G4VParticleChange.hh.

Referenced by CheckIt(), DumpInfo(), GetLocalEnergyDeposit(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4ParticleChangeForGamma::InitializeForPostStep(), InitializeLocalEnergyDeposit(), operator=(), G4ParticleChangeForTransport::operator=(), G4ParticleChangeForMSC::operator=(), G4ParticleChangeForLoss::operator=(), G4ParticleChangeForGamma::operator=(), G4ParticleChangeForDecay::operator=(), G4ParticleChange::operator=(), G4FastStep::operator=(), ProposeLocalEnergyDeposit(), G4ParticleChangeForLoss::UpdateStepForAlongStep(), G4ParticleChangeForGamma::UpdateStepForAtRest(), G4ParticleChangeForLoss::UpdateStepForPostStep(), G4ParticleChangeForGamma::UpdateStepForPostStep(), and UpdateStepInfo().

G4double G4VParticleChange::theNonIonizingEnergyDeposit [protected]

Definition at line 269 of file G4VParticleChange.hh.

Referenced by DumpInfo(), GetNonIonizingEnergyDeposit(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4ParticleChangeForGamma::InitializeForPostStep(), InitializeLocalEnergyDeposit(), operator=(), ProposeNonIonizingEnergyDeposit(), G4ParticleChangeForLoss::UpdateStepForAlongStep(), G4ParticleChangeForLoss::UpdateStepForPostStep(), G4ParticleChangeForGamma::UpdateStepForPostStep(), and UpdateStepInfo().

G4int G4VParticleChange::theNumberOfSecondaries [protected]

Definition at line 248 of file G4VParticleChange.hh.

Referenced by AddSecondary(), G4ParticleChangeForRadDecay::AddSecondary(), Clear(), DumpInfo(), G4VParticleChange(), GetNumberOfSecondaries(), InitializeSecondaries(), operator=(), G4ParticleChangeForTransport::operator=(), G4ParticleChangeForMSC::operator=(), G4ParticleChangeForLoss::operator=(), G4ParticleChangeForGamma::operator=(), G4ParticleChangeForDecay::operator=(), G4ParticleChange::operator=(), G4FastStep::operator=(), SetNumberOfSecondaries(), and ~G4VParticleChange().

G4double G4VParticleChange::theParentGlobalTime [protected]

Definition at line 289 of file G4VParticleChange.hh.

Referenced by CheckSecondary(), InitializeParentGlobalTime(), and operator=().

G4double G4VParticleChange::theParentWeight [protected]

Definition at line 282 of file G4VParticleChange.hh.

Referenced by AddSecondary(), GetParentWeight(), GetWeight(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4ParticleChangeForGamma::InitializeForPostStep(), InitializeParentWeight(), operator=(), G4ParticleChangeForLoss::operator=(), G4ParticleChangeForGamma::operator=(), ProposeWeight(), UpdateStepForAlongStep(), G4ParticleChangeForLoss::UpdateStepForAlongStep(), G4ParticleChange::UpdateStepForAlongStep(), UpdateStepForAtRest(), G4ParticleChangeForGamma::UpdateStepForAtRest(), G4ParticleChangeForDecay::UpdateStepForAtRest(), G4ParticleChange::UpdateStepForAtRest(), UpdateStepForPostStep(), G4ParticleChangeForLoss::UpdateStepForPostStep(), G4ParticleChangeForGamma::UpdateStepForPostStep(), G4ParticleChangeForDecay::UpdateStepForPostStep(), and G4ParticleChange::UpdateStepForPostStep().

G4int G4VParticleChange::theSizeOftheListOfSecondaries [protected]

Definition at line 251 of file G4VParticleChange.hh.

Referenced by AddSecondary(), G4ParticleChangeForRadDecay::AddSecondary(), G4ParticleChangeForTransport::operator=(), G4ParticleChangeForMSC::operator=(), G4FastStep::operator=(), and SetNumberOfSecondaries().

G4TrackStatus G4VParticleChange::theStatusChange [protected]

Definition at line 254 of file G4VParticleChange.hh.

Referenced by DumpInfo(), GetTrackStatus(), G4ParticleChangeForMSC::Initialize(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4ParticleChangeForGamma::InitializeForPostStep(), InitializeStatusChange(), operator=(), G4ParticleChangeForTransport::operator=(), G4ParticleChangeForMSC::operator=(), G4ParticleChangeForLoss::operator=(), G4ParticleChangeForGamma::operator=(), G4ParticleChangeForDecay::operator=(), G4ParticleChange::operator=(), G4FastStep::operator=(), ProposeTrackStatus(), and G4ParticleChangeForMSC::UpdateStepForAlongStep().

G4SteppingControl G4VParticleChange::theSteppingControlFlag [protected]

Definition at line 257 of file G4VParticleChange.hh.

Referenced by DumpInfo(), G4ParticleChangeForGamma::G4ParticleChangeForGamma(), G4ParticleChangeForLoss::G4ParticleChangeForLoss(), GetSteppingControl(), G4FastStep::Initialize(), InitializeSteppingControl(), operator=(), G4ParticleChangeForTransport::operator=(), G4ParticleChangeForMSC::operator=(), G4ParticleChangeForLoss::operator=(), G4ParticleChangeForGamma::operator=(), G4ParticleChangeForDecay::operator=(), G4ParticleChange::operator=(), G4FastStep::operator=(), ProposeSteppingControl(), G4ParticleChangeForTransport::UpdateStepForAlongStep(), and UpdateStepInfo().

G4double G4VParticleChange::theTrueStepLength [protected]

Definition at line 274 of file G4VParticleChange.hh.

Referenced by CheckIt(), DumpInfo(), GetTrueStepLength(), InitializeTrueStepLength(), operator=(), G4ParticleChangeForTransport::operator=(), G4ParticleChangeForMSC::operator=(), G4ParticleChangeForDecay::operator=(), G4ParticleChange::operator=(), G4FastStep::operator=(), ProposeTrueStepLength(), G4ParticleChangeForMSC::UpdateStepForAlongStep(), and UpdateStepInfo().

G4int G4VParticleChange::verboseLevel [protected]

Definition at line 293 of file G4VParticleChange.hh.

Referenced by AddSecondary(), G4ParticleChangeForRadDecay::AddSecondary(), G4FastStep::G4FastStep(), G4ParticleChange::G4ParticleChange(), G4ParticleChangeForDecay::G4ParticleChangeForDecay(), G4ParticleChangeForGamma::G4ParticleChangeForGamma(), G4ParticleChangeForLoss::G4ParticleChangeForLoss(), G4ParticleChangeForMSC::G4ParticleChangeForMSC(), G4ParticleChangeForTransport::G4ParticleChangeForTransport(), GetVerboseLevel(), InitializeSecondaries(), operator=(), G4ParticleChangeForTransport::operator=(), G4ParticleChangeForMSC::operator=(), G4ParticleChangeForLoss::operator=(), G4ParticleChangeForGamma::operator=(), G4ParticleChangeForDecay::operator=(), G4ParticleChange::operator=(), SetNumberOfSecondaries(), SetVerboseLevel(), G4ParticleChangeForTransport::UpdateStepForAtRest(), G4FastStep::~G4FastStep(), G4ParticleChange::~G4ParticleChange(), G4ParticleChangeForDecay::~G4ParticleChangeForDecay(), G4ParticleChangeForGamma::~G4ParticleChangeForGamma(), G4ParticleChangeForLoss::~G4ParticleChangeForLoss(), G4ParticleChangeForMSC::~G4ParticleChangeForMSC(), G4ParticleChangeForTransport::~G4ParticleChangeForTransport(), and ~G4VParticleChange().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:52 2013 for Geant4 by  doxygen 1.4.7