CMS 3D CMS Logo

DummyChargeFlipProcess.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include "G4ParticleTable.hh"
3 #include "Randomize.hh"
4 
6 
7 using namespace CLHEP;
8 
9 DummyChargeFlipProcess::DummyChargeFlipProcess(const G4String& processName) : G4HadronicProcess(processName) {
10  AddDataSet(new G4HadronElasticDataSet);
11 }
12 
14 
15 void DummyChargeFlipProcess::BuildPhysicsTable(const G4ParticleDefinition& aParticleType) {
16  if (!G4HadronicProcess::GetCrossSectionDataStore()) {
17  return;
18  }
19  G4HadronicProcess::GetCrossSectionDataStore()->BuildPhysicsTable(aParticleType);
20 }
21 
22 G4double DummyChargeFlipProcess::GetMicroscopicCrossSection(const G4DynamicParticle* /*aParticle*/,
23  const G4Element* element,
24  G4double /*aTemp*/) {
25  return 30 * millibarn * element->GetN();
26 }
27 
28 G4bool DummyChargeFlipProcess::IsApplicable(const G4ParticleDefinition& aParticleType) {
29  if (aParticleType.GetParticleType() == "rhadron")
30  return true;
31  else
32  return false;
33 }
34 
35 void DummyChargeFlipProcess::DumpPhysicsTable(const G4ParticleDefinition& /*aParticleType*/) {}
36 
37 G4VParticleChange* DummyChargeFlipProcess::PostStepDoIt(const G4Track& aTrack, const G4Step& /*aStep*/) {
38  G4ParticleChange* pc = new G4ParticleChange();
39  pc->Initialize(aTrack);
40  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
41  //G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
42 
43  G4double ParentEnergy = aParticle->GetTotalEnergy();
44  const G4ThreeVector& ParentDirection(aParticle->GetMomentumDirection());
45 
46  G4double energyDeposit = 0.0;
47  G4double finalGlobalTime = aTrack.GetGlobalTime();
48 
49  G4int numberOfSecondaries = 1;
50  pc->SetNumberOfSecondaries(numberOfSecondaries);
51  const G4TouchableHandle& thand = aTrack.GetTouchableHandle();
52 
53  // get current position of the track
54  aTrack.GetPosition();
55  // create a new track object
56  G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
57  float randomParticle = G4UniformRand();
58  G4ParticleDefinition* newType;
59  if (randomParticle < 0.333)
60  newType = particleTable->FindParticle(1009213);
61  else if (randomParticle > 0.667)
62  newType = particleTable->FindParticle(-1009213);
63  else
64  newType = particleTable->FindParticle(1009113);
65 
66  G4cout << "RHADRON: New charge = " << newType->GetPDGCharge() << G4endl;
67 
68  G4DynamicParticle* newP = new G4DynamicParticle(newType, ParentDirection, ParentEnergy);
69  G4Track* secondary = new G4Track(newP, finalGlobalTime, aTrack.GetPosition());
70  // switch on good for tracking flag
71  secondary->SetGoodForTrackingFlag();
72  secondary->SetTouchableHandle(thand);
73  // add the secondary track in the List
74  pc->AddSecondary(secondary);
75 
76  // Kill the parent particle
77  pc->ProposeTrackStatus(fStopAndKill);
78  pc->ProposeLocalEnergyDeposit(energyDeposit);
79  pc->ProposeGlobalTime(finalGlobalTime);
80  // Clear NumberOfInteractionLengthLeft
81  ClearNumberOfInteractionLengthLeft();
82 
83  return pc;
84 }
void DumpPhysicsTable(const G4ParticleDefinition &aParticleType)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4double GetMicroscopicCrossSection(const G4DynamicParticle *aParticle, const G4Element *anElement, G4double aTemp)
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
DummyChargeFlipProcess(const G4String &processName="Dummy")
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override