CMS 3D CMS Logo

DummyChargeFlipProcess.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include "G4ParticleTable.hh"
00003 #include "Randomize.hh"
00004 
00005 #include "SimG4Core/CustomPhysics/interface/DummyChargeFlipProcess.h"
00006 
00007 DummyChargeFlipProcess::
00008 DummyChargeFlipProcess(const G4String& processName) : 
00009       G4HadronicProcess(processName)
00010 {
00011   AddDataSet(new G4HadronElasticDataSet);
00012 }
00013 
00014 DummyChargeFlipProcess::~DummyChargeFlipProcess()
00015 {
00016 }
00017  
00018 
00019 void DummyChargeFlipProcess::
00020 BuildPhysicsTable(const G4ParticleDefinition& aParticleType)
00021 {
00022    if (!G4HadronicProcess::GetCrossSectionDataStore()) {
00023       return;
00024    }
00025    G4HadronicProcess::GetCrossSectionDataStore()->BuildPhysicsTable(aParticleType);
00026 }
00027 
00028 
00029 G4double DummyChargeFlipProcess::
00030 GetMicroscopicCrossSection(
00031       const G4DynamicParticle* /*aParticle*/, const G4Element* element, G4double /*aTemp*/)
00032 {
00033 
00034    return 30*millibarn*element->GetN(); 
00035 }
00036 
00037 
00038 G4bool
00039 DummyChargeFlipProcess::
00040 IsApplicable(const G4ParticleDefinition& aParticleType)
00041 {
00042    if(aParticleType.GetParticleType()  == "rhadron")
00043     return true;
00044   else
00045     return false;
00046 }
00047 
00048 
00049 void 
00050 DummyChargeFlipProcess::
00051 DumpPhysicsTable(const G4ParticleDefinition& /*aParticleType*/)
00052 {
00053 
00054 }
00055 
00056 
00057 G4VParticleChange *DummyChargeFlipProcess::PostStepDoIt(
00058   const G4Track &aTrack, const G4Step &/*aStep*/)
00059 {
00060   SetDispatch(this);
00061   G4ParticleChange * pc = new G4ParticleChange();
00062   pc->Initialize(aTrack);
00063   const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
00064   G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
00065 
00066   G4double   ParentEnergy  = aParticle->GetTotalEnergy();
00067   G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
00068 
00069   G4double energyDeposit = 0.0;
00070   G4double finalGlobalTime = aTrack.GetGlobalTime();
00071 
00072   G4int numberOfSecondaries = 1;
00073   pc->SetNumberOfSecondaries(numberOfSecondaries);
00074   const G4TouchableHandle thand = aTrack.GetTouchableHandle();
00075 
00076      // get current position of the track
00077      aTrack.GetPosition();
00078      // create a new track object
00079       G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
00080       float randomParticle = G4UniformRand();
00081       G4ParticleDefinition * newType = aParticleDef;
00082       if(randomParticle < 0.333)
00083         newType=particleTable->FindParticle(1009213);
00084       else if(randomParticle > 0.667)
00085         newType=particleTable->FindParticle(-1009213);
00086       else
00087         newType=particleTable->FindParticle(1009113);
00088      
00089       std::cout << "RHADRON: New charge = " << newType->GetPDGCharge() << std::endl;
00090        
00091      G4DynamicParticle * newP =  new G4DynamicParticle(newType,ParentDirection,ParentEnergy);
00092      G4Track* secondary = new G4Track( newP ,
00093                                       finalGlobalTime ,
00094                                        aTrack.GetPosition()
00095      );
00096      // switch on good for tracking flag
00097      secondary->SetGoodForTrackingFlag();
00098      secondary->SetTouchableHandle(thand);
00099      // add the secondary track in the List
00100      pc->AddSecondary(secondary);
00101 
00102   // Kill the parent particle
00103   pc->ProposeTrackStatus( fStopAndKill ) ;
00104   pc->ProposeLocalEnergyDeposit(energyDeposit); 
00105   pc->ProposeGlobalTime( finalGlobalTime );
00106   // Clear NumberOfInteractionLengthLeft
00107   ClearNumberOfInteractionLengthLeft();
00108 
00109 
00110 
00111    return pc;
00112 
00113 }
00114 

Generated on Tue Jun 9 17:47:00 2009 for CMSSW by  doxygen 1.5.4