CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
11  G4HadronicProcess(processName)
12 {
13  AddDataSet(new G4HadronElasticDataSet);
14 }
15 
17 {
18 }
19 
20 
22 BuildPhysicsTable(const G4ParticleDefinition& aParticleType)
23 {
24  if (!G4HadronicProcess::GetCrossSectionDataStore()) {
25  return;
26  }
27  G4HadronicProcess::GetCrossSectionDataStore()->BuildPhysicsTable(aParticleType);
28 }
29 
30 
33  const G4DynamicParticle* /*aParticle*/, const G4Element* element, G4double /*aTemp*/)
34 {
35 
36  return 30*millibarn*element->GetN();
37 }
38 
39 
40 G4bool
42 IsApplicable(const G4ParticleDefinition& aParticleType)
43 {
44  if(aParticleType.GetParticleType() == "rhadron")
45  return true;
46  else
47  return false;
48 }
49 
50 
51 void
53 DumpPhysicsTable(const G4ParticleDefinition& /*aParticleType*/)
54 {
55 
56 }
57 
58 
60  const G4Track &aTrack, const G4Step &/*aStep*/)
61 {
62  G4ParticleChange * pc = new G4ParticleChange();
63  pc->Initialize(aTrack);
64  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
65  //G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
66 
67  G4double ParentEnergy = aParticle->GetTotalEnergy();
68  G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
69 
70  G4double energyDeposit = 0.0;
71  G4double finalGlobalTime = aTrack.GetGlobalTime();
72 
73  G4int numberOfSecondaries = 1;
74  pc->SetNumberOfSecondaries(numberOfSecondaries);
75  const G4TouchableHandle thand = aTrack.GetTouchableHandle();
76 
77  // get current position of the track
78  aTrack.GetPosition();
79  // create a new track object
80  G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
81  float randomParticle = G4UniformRand();
82  G4ParticleDefinition * newType;
83  if(randomParticle < 0.333)
84  newType=particleTable->FindParticle(1009213);
85  else if(randomParticle > 0.667)
86  newType=particleTable->FindParticle(-1009213);
87  else
88  newType=particleTable->FindParticle(1009113);
89 
90  G4cout << "RHADRON: New charge = " << newType->GetPDGCharge() << G4endl;
91 
92  G4DynamicParticle * newP = new G4DynamicParticle(newType,ParentDirection,ParentEnergy);
93  G4Track* secondary = new G4Track( newP ,
94  finalGlobalTime ,
95  aTrack.GetPosition()
96  );
97  // switch on good for tracking flag
98  secondary->SetGoodForTrackingFlag();
99  secondary->SetTouchableHandle(thand);
100  // add the secondary track in the List
101  pc->AddSecondary(secondary);
102 
103  // Kill the parent particle
104  pc->ProposeTrackStatus( fStopAndKill ) ;
105  pc->ProposeLocalEnergyDeposit(energyDeposit);
106  pc->ProposeGlobalTime( finalGlobalTime );
107  // Clear NumberOfInteractionLengthLeft
108  ClearNumberOfInteractionLengthLeft();
109 
110  return pc;
111 }
112 
void DumpPhysicsTable(const G4ParticleDefinition &aParticleType)
G4double GetMicroscopicCrossSection(const G4DynamicParticle *aParticle, const G4Element *anElement, G4double aTemp)
DummyChargeFlipProcess(const G4String &processName="Dummy")
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)