CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DummyChargeFlipProcess.cc
Go to the documentation of this file.
2 
3 #include <iostream>
4 #include "G4ParticleTable.hh"
5 #include "Randomize.hh"
6 #include "G4NeutronElasticXS.hh"
7 #include "G4Step.hh"
8 #include "G4TrackStatus.hh"
9 #include "G4Element.hh"
10 
11 using namespace CLHEP;
12 
13 DummyChargeFlipProcess::DummyChargeFlipProcess(const G4String& pname) : G4HadronicProcess(pname, fHadronic) {
14  AddDataSet(new G4NeutronElasticXS());
15  fPartChange = new G4ParticleChange();
16 }
17 
19 
20 G4bool DummyChargeFlipProcess::IsApplicable(const G4ParticleDefinition& aParticleType) {
21  return (aParticleType.GetParticleType() == "rhadron");
22 }
23 
24 G4VParticleChange* DummyChargeFlipProcess::PostStepDoIt(const G4Track& aTrack, const G4Step&) {
25  fPartChange->Initialize(aTrack);
26  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
27 
28  G4double ParentEnergy = aParticle->GetTotalEnergy();
29  const G4ThreeVector& ParentDirection(aParticle->GetMomentumDirection());
30 
31  G4double energyDeposit = 0.0;
32  G4double finalGlobalTime = aTrack.GetGlobalTime();
33 
34  G4int numberOfSecondaries = 1;
35  fPartChange->SetNumberOfSecondaries(numberOfSecondaries);
36  const G4TouchableHandle& thand = aTrack.GetTouchableHandle();
37 
38  // get current position of the track
39  aTrack.GetPosition();
40  // create a new track object
41  G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
42  float randomParticle = G4UniformRand();
43  G4ParticleDefinition* newType;
44  if (randomParticle < 0.333)
45  newType = particleTable->FindParticle(1009213);
46  else if (randomParticle > 0.667)
47  newType = particleTable->FindParticle(-1009213);
48  else
49  newType = particleTable->FindParticle(1009113);
50 
51  //G4cout << "RHADRON: New charge = " << newType->GetPDGCharge() << G4endl;
52 
53  G4DynamicParticle* newP = new G4DynamicParticle(newType, ParentDirection, ParentEnergy);
54  G4Track* secondary = new G4Track(newP, finalGlobalTime, aTrack.GetPosition());
55  // switch on good for tracking flag
56  secondary->SetGoodForTrackingFlag();
57  secondary->SetTouchableHandle(thand);
58  // add the secondary track in the List
59  fPartChange->AddSecondary(secondary);
60 
61  // Kill the parent particle
62  fPartChange->ProposeTrackStatus(fStopAndKill);
63  fPartChange->ProposeLocalEnergyDeposit(energyDeposit);
64  fPartChange->ProposeGlobalTime(finalGlobalTime);
65 
66  ClearNumberOfInteractionLengthLeft();
67 
68  return fPartChange;
69 }
G4ParticleChange * fPartChange
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
DummyChargeFlipProcess(const G4String &processName="Dummy")
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override