G4ParticleChangeForDecay Class Reference

#include <G4ParticleChangeForDecay.hh>

Inheritance diagram for G4ParticleChangeForDecay:

G4VParticleChange G4ParticleChangeForRadDecay

Public Member Functions

 G4ParticleChangeForDecay ()
virtual ~G4ParticleChangeForDecay ()
G4bool operator== (const G4ParticleChangeForDecay &right) const
G4bool operator!= (const G4ParticleChangeForDecay &right) const
virtual G4StepUpdateStepForAtRest (G4Step *Step)
virtual G4StepUpdateStepForPostStep (G4Step *Step)
virtual void Initialize (const G4Track &)
void ProposeGlobalTime (G4double t)
void ProposeLocalTime (G4double t)
G4double GetGlobalTime (G4double timeDelay=0.0) const
G4double GetLocalTime (G4double timeDelay=0.0) const
const G4ThreeVectorGetPolarization () const
void ProposePolarization (G4double Px, G4double Py, G4double Pz)
void ProposePolarization (const G4ThreeVector &finalPoralization)
virtual void DumpInfo () const
virtual G4bool CheckIt (const G4Track &)

Protected Member Functions

 G4ParticleChangeForDecay (const G4ParticleChangeForDecay &right)
G4ParticleChangeForDecayoperator= (const G4ParticleChangeForDecay &right)

Protected Attributes

G4double theGlobalTime0
G4double theLocalTime0
G4double theTimeChange
G4ThreeVector thePolarizationChange

Detailed Description

Definition at line 54 of file G4ParticleChangeForDecay.hh.


Constructor & Destructor Documentation

G4ParticleChangeForDecay::G4ParticleChangeForDecay (  ) 

Definition at line 48 of file G4ParticleChangeForDecay.cc.

References G4cout, G4endl, and G4VParticleChange::verboseLevel.

00049   : G4VParticleChange(), 
00050     theGlobalTime0(0.),theLocalTime0(0.),theTimeChange(0.)
00051 {
00052 #ifdef G4VERBOSE
00053   if (verboseLevel>2) {
00054     G4cout << "G4ParticleChangeForDecay::G4ParticleChangeForDecay() " << G4endl;
00055   }
00056 #endif
00057 }

G4ParticleChangeForDecay::~G4ParticleChangeForDecay (  )  [virtual]

Definition at line 59 of file G4ParticleChangeForDecay.cc.

References G4cout, G4endl, and G4VParticleChange::verboseLevel.

00060 {
00061 #ifdef G4VERBOSE
00062   if (verboseLevel>2) {
00063     G4cout << "G4ParticleChangeForDecay::~G4ParticleChangeForDecay() " << G4endl;
00064   }
00065 #endif
00066 }

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

Definition at line 68 of file G4ParticleChangeForDecay.cc.

References theGlobalTime0, theLocalTime0, thePolarizationChange, and theTimeChange.

00068                                                                                        : G4VParticleChange(right)
00069 {
00070   theGlobalTime0 = right.theGlobalTime0;
00071   theLocalTime0 = right.theLocalTime0;
00072   theTimeChange = right.theTimeChange;
00073   thePolarizationChange = right.thePolarizationChange;
00074 }


Member Function Documentation

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

Reimplemented from G4VParticleChange.

Definition at line 193 of file G4ParticleChangeForDecay.cc.

References G4VParticleChange::accuracyForException, G4VParticleChange::accuracyForWarning, G4VParticleChange::CheckIt(), DumpInfo(), EventMustBeAborted, G4cout, G4endl, G4Exception(), G4Track::GetDefinition(), G4Track::GetKineticEnergy(), G4Track::GetLocalTime(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), ns, theLocalTime0, and theTimeChange.

Referenced by UpdateStepForAtRest().

00194 {
00195   G4bool    exitWithError = false;
00196 
00197   G4double  accuracy;
00198 
00199   // local time should not go back
00200   G4bool itsOK =true;
00201   accuracy = -1.0*(theTimeChange - theLocalTime0)/ns;
00202   if (accuracy > accuracyForWarning) {
00203     itsOK = false;
00204     exitWithError = (accuracy > accuracyForException);
00205 #ifdef G4VERBOSE
00206     G4cout << "  G4ParticleChangeForDecay::CheckIt    : ";
00207     G4cout << "the local time goes back  !!" 
00208            << "  Difference:  " << accuracy  << "[ns] " <<G4endl;
00209     G4cout << aTrack.GetDefinition()->GetParticleName()
00210            << " E=" << aTrack.GetKineticEnergy()/MeV
00211            << " pos=" << aTrack.GetPosition().x()/m
00212            << ", " << aTrack.GetPosition().y()/m
00213            << ", " << aTrack.GetPosition().z()/m
00214            <<G4endl;
00215 #endif
00216   }
00217 
00218   // dump out information of this particle change
00219 #ifdef G4VERBOSE
00220   if (!itsOK) DumpInfo();
00221 #endif
00222 
00223   // Exit with error
00224   if (exitWithError) {
00225     G4Exception("G4ParticleChangeForDecay::CheckIt",
00226                 "TRACK005",EventMustBeAborted,
00227                 "time was  illegal");
00228   } 
00229 
00230   // correction
00231   if (!itsOK) {
00232     theTimeChange = aTrack.GetLocalTime();
00233   }
00234 
00235   itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
00236   return itsOK;
00237 }

void G4ParticleChangeForDecay::DumpInfo (  )  const [virtual]

Reimplemented from G4VParticleChange.

Definition at line 182 of file G4ParticleChangeForDecay.cc.

References G4VParticleChange::DumpInfo(), G4cout, G4endl, ns, and theTimeChange.

Referenced by CheckIt().

00183 {
00184 // Show header
00185   G4VParticleChange::DumpInfo();
00186 
00187   G4int oldprc = G4cout.precision(3);
00188   G4cout << "        Time (ns)           : " 
00189          << std::setw(20) << theTimeChange/ns << G4endl;
00190   G4cout.precision(oldprc);
00191 }

G4double G4ParticleChangeForDecay::GetGlobalTime ( G4double  timeDelay = 0.0  )  const [inline]

Definition at line 132 of file G4ParticleChangeForDecay.hh.

References theGlobalTime0, theLocalTime0, and theTimeChange.

Referenced by UpdateStepForAtRest().

00133 {
00134   //  Convert the time delay to the global time.
00135   return theGlobalTime0 + (theTimeChange-theLocalTime0) + timeDelay;
00136 }

G4double G4ParticleChangeForDecay::GetLocalTime ( G4double  timeDelay = 0.0  )  const [inline]

Definition at line 145 of file G4ParticleChangeForDecay.hh.

References theTimeChange.

00146 {
00147   //  Convert the time delay to the local time.
00148   return theTimeChange + timeDelay;
00149 }

const G4ThreeVector * G4ParticleChangeForDecay::GetPolarization (  )  const [inline]

Definition at line 152 of file G4ParticleChangeForDecay.hh.

References thePolarizationChange.

00153 {
00154   return &thePolarizationChange;
00155 }

void G4ParticleChangeForDecay::Initialize ( const G4Track  )  [virtual]

Reimplemented from G4VParticleChange.

Definition at line 124 of file G4ParticleChangeForDecay.cc.

References G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4Track::GetLocalTime(), G4DynamicParticle::GetPolarization(), G4VParticleChange::Initialize(), theGlobalTime0, theLocalTime0, thePolarizationChange, and theTimeChange.

Referenced by G4UnknownDecay::DecayIt(), G4RadioactiveDecay::DecayIt(), and G4Decay::DecayIt().

00125 {
00126   // use base class's method at first
00127   G4VParticleChange::Initialize(track);
00128 
00129   const G4DynamicParticle*  pParticle = track.GetDynamicParticle();
00130 
00131   // set TimeChange equal to local time of the parent track
00132   theTimeChange                = track.GetLocalTime();
00133 
00134   // set initial Local/Global time of the parent track
00135   theLocalTime0           = track.GetLocalTime();
00136   theGlobalTime0          = track.GetGlobalTime();
00137 
00138   // set the Polarization equal to those of the parent track
00139   thePolarizationChange  = pParticle->GetPolarization();
00140 }

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

Definition at line 116 of file G4ParticleChangeForDecay.cc.

00117 {
00118    return ((G4VParticleChange *)this != (G4VParticleChange *) &right);
00119 }

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

Definition at line 77 of file G4ParticleChangeForDecay.cc.

References G4cout, G4FastVector< Type, N >::SetElement(), theGlobalTime0, G4VParticleChange::theListOfSecondaries, G4VParticleChange::theLocalEnergyDeposit, theLocalTime0, G4VParticleChange::theNumberOfSecondaries, thePolarizationChange, G4VParticleChange::theStatusChange, G4VParticleChange::theSteppingControlFlag, theTimeChange, G4VParticleChange::theTrueStepLength, and G4VParticleChange::verboseLevel.

00078 {
00079   if (this != &right){
00080     if (theNumberOfSecondaries>0) {
00081 #ifdef G4VERBOSE
00082       if (verboseLevel>0) {
00083         G4cout << "G4ParticleChangeForDecay: assignment operator Warning  ";
00084         G4cout << "theListOfSecondaries is not empty ";
00085        }
00086 #endif
00087       for (G4int index= 0; index<theNumberOfSecondaries; index++){
00088         if ( (*theListOfSecondaries)[index] ) delete (*theListOfSecondaries)[index] ;
00089       }
00090     }
00091     delete theListOfSecondaries; 
00092     theListOfSecondaries =  new G4TrackFastVector();
00093     theNumberOfSecondaries = right.theNumberOfSecondaries;
00094     for (G4int index = 0; index<theNumberOfSecondaries; index++){
00095       G4Track* newTrack =  new G4Track(*((*right.theListOfSecondaries)[index] ));
00096       theListOfSecondaries->SetElement(index, newTrack);                            }
00097 
00098     theStatusChange = right.theStatusChange;
00099     theTrueStepLength = right.theTrueStepLength;
00100     theLocalEnergyDeposit = right.theLocalEnergyDeposit;
00101     theSteppingControlFlag = right.theSteppingControlFlag;
00102     
00103     theGlobalTime0 = right.theGlobalTime0;
00104     theLocalTime0 = right.theLocalTime0;
00105     theTimeChange = right.theTimeChange;
00106     thePolarizationChange = right.thePolarizationChange;
00107   }
00108   return *this;
00109 }

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

Definition at line 111 of file G4ParticleChangeForDecay.cc.

00112 {
00113    return ((G4VParticleChange *)this == (G4VParticleChange *) &right);
00114 }

void G4ParticleChangeForDecay::ProposeGlobalTime ( G4double  t  )  [inline]

Definition at line 126 of file G4ParticleChangeForDecay.hh.

References theGlobalTime0, theLocalTime0, and theTimeChange.

Referenced by G4UnknownDecay::DecayIt().

00127 {
00128   theTimeChange = (t-theGlobalTime0) + theLocalTime0;
00129 }

void G4ParticleChangeForDecay::ProposeLocalTime ( G4double  t  )  [inline]

Definition at line 139 of file G4ParticleChangeForDecay.hh.

References theTimeChange.

Referenced by G4RadioactiveDecay::DecayIt(), and G4Decay::DecayIt().

00140 {
00141   theTimeChange = t;
00142 }

void G4ParticleChangeForDecay::ProposePolarization ( const G4ThreeVector finalPoralization  )  [inline]

Definition at line 158 of file G4ParticleChangeForDecay.hh.

References thePolarizationChange.

00159 {
00160   thePolarizationChange = finalPoralization;
00161 }

void G4ParticleChangeForDecay::ProposePolarization ( G4double  Px,
G4double  Py,
G4double  Pz 
) [inline]

Definition at line 164 of file G4ParticleChangeForDecay.hh.

References thePolarizationChange.

Referenced by G4DecayWithSpin::DecayIt().

00168 {
00169   thePolarizationChange.setX(Px);
00170   thePolarizationChange.setY(Py);
00171   thePolarizationChange.setZ(Pz);
00172 }

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

Reimplemented from G4VParticleChange.

Definition at line 157 of file G4ParticleChangeForDecay.cc.

References G4StepPoint::AddProperTime(), CheckIt(), G4VParticleChange::debugFlag, GetGlobalTime(), G4Step::GetPostStepPoint(), G4Step::GetTrack(), G4VParticleChange::isParentWeightProposed, G4StepPoint::SetGlobalTime(), G4StepPoint::SetLocalTime(), G4StepPoint::SetPolarization(), G4StepPoint::SetWeight(), theLocalTime0, G4VParticleChange::theParentWeight, thePolarizationChange, theTimeChange, and G4VParticleChange::UpdateStepInfo().

00158 { 
00159   // A physics process always calculates the final state of the particle
00160 
00161   G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint(); 
00162 
00163   // update polarization
00164   pPostStepPoint->SetPolarization( thePolarizationChange );
00165  
00166   // update time
00167   pPostStepPoint->SetGlobalTime( GetGlobalTime() );
00168   pPostStepPoint->SetLocalTime(  theTimeChange ); 
00169   pPostStepPoint->AddProperTime (theTimeChange-theLocalTime0);
00170 
00171 #ifdef G4VERBOSE
00172   G4Track*     aTrack  = pStep->GetTrack();
00173   if (debugFlag) CheckIt(*aTrack);
00174 #endif
00175 
00176   if (isParentWeightProposed )pPostStepPoint->SetWeight( theParentWeight );
00177 
00178   //  Update the G4Step specific attributes 
00179   return UpdateStepInfo(pStep);
00180 }

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

Reimplemented from G4VParticleChange.

Definition at line 146 of file G4ParticleChangeForDecay.cc.

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

00147 { 
00148   if (isParentWeightProposed ){
00149     pStep->GetPostStepPoint()->SetWeight( theParentWeight );
00150   }
00151 
00152   //  Update the G4Step specific attributes 
00153   return UpdateStepInfo(pStep);
00154 }


Field Documentation

G4double G4ParticleChangeForDecay::theGlobalTime0 [protected]

Definition at line 109 of file G4ParticleChangeForDecay.hh.

Referenced by G4ParticleChangeForDecay(), GetGlobalTime(), Initialize(), operator=(), and ProposeGlobalTime().

G4double G4ParticleChangeForDecay::theLocalTime0 [protected]

Definition at line 111 of file G4ParticleChangeForDecay.hh.

Referenced by CheckIt(), G4ParticleChangeForDecay(), GetGlobalTime(), Initialize(), operator=(), ProposeGlobalTime(), and UpdateStepForAtRest().

G4ThreeVector G4ParticleChangeForDecay::thePolarizationChange [protected]

Definition at line 117 of file G4ParticleChangeForDecay.hh.

Referenced by G4ParticleChangeForDecay(), GetPolarization(), Initialize(), operator=(), ProposePolarization(), and UpdateStepForAtRest().

G4double G4ParticleChangeForDecay::theTimeChange [protected]

Definition at line 114 of file G4ParticleChangeForDecay.hh.

Referenced by CheckIt(), DumpInfo(), G4ParticleChangeForDecay(), GetGlobalTime(), GetLocalTime(), Initialize(), operator=(), ProposeGlobalTime(), ProposeLocalTime(), and UpdateStepForAtRest().


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