/scratch0/jsnuveri/BDSIM/BDSIMgit/bdsim/src/BDSLaserCompton.cc

00001 /* BDSIM code.    Version 1.0
00002    Author: Grahame A. Blair, Royal Holloway, Univ. of London.
00003    Last modified 24.7.2002
00004    Copyright (c) 2002 by G.A.Blair.  ALL RIGHTS RESERVED. 
00005 */
00006 //      ------------ BDSLaserCompton physics process --------
00007 //                     by Grahame Blair, 18 October 2001
00008 #include "BDSGlobalConstants.hh" 
00009 
00010 #include "BDSLaserCompton.hh"
00011 #include "G4ios.hh"
00012 #include "G4Gamma.hh"
00013 #include "CLHEP/Units/PhysicalConstants.h"
00014 
00015 BDSLaserCompton::BDSLaserCompton(const G4String& processName)
00016   :  G4VDiscreteProcess(processName),//isInitialised(false),
00017      itsLaserEnergy(0.0)
00018 {
00019   itsLaserWavelength=BDSGlobalConstants::Instance()->GetLaserwireWavelength();
00020   itsLaserDirection=BDSGlobalConstants::Instance()->GetLaserwireDir();
00021 
00022 
00023   //    if(itsLaserWavelength<=0.)
00024   //     {G4Exception("BDSLaserCompton: Invalid Wavelength");}
00025   // itsLaserEnergy=CLHEP::twopi*CLHEP::hbarc/itsLaserWavelength;
00026  // point laserwire in x:     P_x        Py Pz   E
00027  //G4LorentzVector Laser4mom(itsLaserEnergy,0,0,itsLaserEnergy);
00028  //itsComptonEngine=new BDSComptonEngine(Laser4mom);
00029   itsComptonEngine=new BDSComptonEngine();
00030 } 
00031 
00032  
00033 BDSLaserCompton::~BDSLaserCompton()
00034 {
00035   delete itsComptonEngine;
00036 }
00037 
00038 
00039 G4VParticleChange* BDSLaserCompton::PostStepDoIt(const G4Track& trackData,
00040                                                  const G4Step& stepData)
00041 {
00042  
00043  
00044  aParticleChange.Initialize(trackData);
00045  
00046  // ensure that Laserwire can only occur once in an event
00047  G4cout << "FireLaserCompton == " << FireLaserCompton << G4endl;
00048  if(!FireLaserCompton){
00049    return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
00050  }
00051  G4Material* aMaterial=trackData.GetMaterial() ;
00052  
00053  if(aMaterial==BDSMaterials::Instance()->GetMaterial("LaserVac"))
00054    {
00055      itsLaserWavelength=BDSGlobalConstants::Instance()->GetLaserwireWavelength();
00056      itsLaserDirection=BDSGlobalConstants::Instance()->GetLaserwireDir();
00057      
00058      //G4cout << "&&&&&" << itsLaserDirection << "&&&&&\n";
00059      if(itsLaserWavelength<=0.)
00060        {G4Exception("BDSLaserCompton::PostStepDoIt - Invalid Wavelength", "-1", FatalException, "");}
00061      itsLaserEnergy=CLHEP::twopi*CLHEP::hbarc/itsLaserWavelength;
00062      // point laserwire in x:     P_x        Py Pz   E
00063      G4LorentzVector Laser4mom(itsLaserEnergy*itsLaserDirection.unit(),itsLaserEnergy);
00064      
00065      const G4DynamicParticle* aDynamicParticle=trackData.GetDynamicParticle();
00066      
00067      itsComptonEngine->
00068        SetIncomingElectron4Vec(aDynamicParticle->Get4Momentum());
00069      itsComptonEngine->SetIncomingPhoton4Vec(Laser4mom);
00070      
00071      itsComptonEngine->PerformCompton();
00072      
00073      if(BDSGlobalConstants::Instance()->GetLaserwireTrackPhotons())
00074        {
00075 
00076          // create G4DynamicParticle object for the Gamma 
00077          G4LorentzVector ScatGam=itsComptonEngine->GetScatteredGamma();
00078          //      G4cout<<" Gamma Energy="<<ScatGam.e()/GeV<<" GeV"<<G4endl;
00079          G4DynamicParticle* aGamma= 
00080            new G4DynamicParticle (G4Gamma::Gamma(), 
00081                                   ScatGam.vect().unit(),// direction 
00082                                   ScatGam.e());
00083          
00084          aParticleChange.SetNumberOfSecondaries(1);
00085          aParticleChange.AddSecondary(aGamma); 
00086          if(!BDSGlobalConstants::Instance()->GetLaserwireTrackElectrons())
00087            {
00088              aParticleChange.ProposeEnergy( 0. );
00089              aParticleChange.ProposeLocalEnergyDeposit (0.);
00090              aParticleChange.ProposeTrackStatus(fStopAndKill);
00091            }
00092        }
00093      else
00094        {
00095          aParticleChange.SetNumberOfSecondaries(0);
00096          aParticleChange.ProposeLocalEnergyDeposit (0.);
00097        }
00098      //
00099      // Update the incident particle 
00100      //
00101 
00102     
00103      G4double NewKinEnergy=
00104        itsComptonEngine->GetScatteredElectron().e()-CLHEP::electron_mass_c2;
00105      
00106      //  G4double NewKinEnergy=0; // tmp to track photon only
00107      
00108      G4LorentzVector ScatEl=itsComptonEngine->GetScatteredElectron();
00109      
00110      
00111      if (NewKinEnergy > 0.)
00112        {
00113          aParticleChange.ProposeMomentumDirection(ScatEl.vect().unit());
00114          aParticleChange.ProposeEnergy(NewKinEnergy);
00115          aParticleChange.ProposeLocalEnergyDeposit (0.); 
00116        } 
00117      else
00118        { 
00119          aParticleChange.ProposeEnergy( 0. );
00120          aParticleChange.ProposeLocalEnergyDeposit (0.);
00121          G4double charge= aDynamicParticle->GetCharge();
00122          if (charge<0.) aParticleChange.ProposeTrackStatus(fStopAndKill);
00123          else       aParticleChange.ProposeTrackStatus(fStopButAlive);
00124        }    
00125    }
00126  
00127  FireLaserCompton=false;
00128  
00129  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
00130 }

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7