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 "G4UnitsTable.hh"
00013 
00014 #if G4VERSION_NUMBER > 899
00015 BDSLaserCompton::BDSLaserCompton(const G4String& processName)
00016   :  G4VDiscreteProcess(processName),//isInitialised(false),
00017      itsLaserEnergy(0.0)
00018 #else
00019 BDSLaserCompton::BDSLaserCompton(const G4String& processName)
00020      :  G4VeEnergyLoss(processName), itsLaserEnergy(0.0)
00021 #endif
00022 {
00023   itsLaserWavelength=BDSGlobalConstants::Instance()->GetLaserwireWavelength();
00024   itsLaserDirection=BDSGlobalConstants::Instance()->GetLaserwireDir();
00025 
00026 
00027   //    if(itsLaserWavelength<=0.)
00028   //     {G4Exception("BDSLaserCompton: Invalid Wavelength");}
00029   // itsLaserEnergy=twopi*hbarc/itsLaserWavelength;
00030  // point laserwire in x:     P_x        Py Pz   E
00031  //G4LorentzVector Laser4mom(itsLaserEnergy,0,0,itsLaserEnergy);
00032  //itsComptonEngine=new BDSComptonEngine(Laser4mom);
00033   itsComptonEngine=new BDSComptonEngine();
00034 } 
00035 
00036  
00037 BDSLaserCompton::~BDSLaserCompton()
00038 {
00039   delete itsComptonEngine;
00040 }
00041 
00042 
00043 G4VParticleChange* BDSLaserCompton::PostStepDoIt(const G4Track& trackData,
00044                                                  const G4Step& stepData)
00045 {
00046  
00047  
00048  aParticleChange.Initialize(trackData);
00049  
00050  // ensure that Laserwire can only occur once in an event
00051  G4cout << "FireLaserCompton == " << FireLaserCompton << G4endl;
00052  if(!FireLaserCompton){
00053          #if G4VERSION_NUMBER > 899
00054    return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
00055    #else
00056          return G4VContinuousDiscreteProcess::PostStepDoIt(trackData,stepData);
00057          #endif
00058  }
00059  G4Material* aMaterial=trackData.GetMaterial() ;
00060  
00061  if(aMaterial==BDSMaterials::Instance()->GetMaterial("LaserVac"))
00062    {
00063      itsLaserWavelength=BDSGlobalConstants::Instance()->GetLaserwireWavelength();
00064      itsLaserDirection=BDSGlobalConstants::Instance()->GetLaserwireDir();
00065      
00066      //G4cout << "&&&&&" << itsLaserDirection << "&&&&&\n";
00067      if(itsLaserWavelength<=0.)
00068        {G4Exception("BDSLaserCompton::PostStepDoIt - Invalid Wavelength", "-1", FatalException, "");}
00069      itsLaserEnergy=twopi*hbarc/itsLaserWavelength;
00070      // point laserwire in x:     P_x        Py Pz   E
00071      G4LorentzVector Laser4mom(itsLaserEnergy*itsLaserDirection.unit(),itsLaserEnergy);
00072      
00073      const G4DynamicParticle* aDynamicParticle=trackData.GetDynamicParticle();
00074      
00075      itsComptonEngine->
00076        SetIncomingElectron4Vec(aDynamicParticle->Get4Momentum());
00077      itsComptonEngine->SetIncomingPhoton4Vec(Laser4mom);
00078      
00079      itsComptonEngine->PerformCompton();
00080      
00081      if(BDSGlobalConstants::Instance()->GetLaserwireTrackPhotons())
00082        {
00083 
00084          // create G4DynamicParticle object for the Gamma 
00085          G4LorentzVector ScatGam=itsComptonEngine->GetScatteredGamma();
00086          //      G4cout<<" Gamma Energy="<<ScatGam.e()/GeV<<" GeV"<<G4endl;
00087          G4DynamicParticle* aGamma= 
00088            new G4DynamicParticle (G4Gamma::Gamma(), 
00089                                   ScatGam.vect().unit(),// direction 
00090                                   ScatGam.e());
00091          
00092          aParticleChange.SetNumberOfSecondaries(1);
00093          aParticleChange.AddSecondary(aGamma); 
00094          if(!BDSGlobalConstants::Instance()->GetLaserwireTrackElectrons())
00095                    {
00096 #if G4VERSION_NUMBER > 699
00097              aParticleChange.ProposeEnergy( 0. );
00098              aParticleChange.ProposeLocalEnergyDeposit (0.);
00099              aParticleChange.ProposeTrackStatus(fStopAndKill);
00100 #else
00101              aParticleChange.SetEnergyChange( 0. );
00102              aParticleChange.SetLocalEnergyDeposit (0.);
00103              aParticleChange.SetStatusChange(fStopAndKill);
00104 #endif
00105            }
00106        }
00107      else
00108        {
00109 #if G4VERSION_NUMBER > 699
00110          aParticleChange.SetNumberOfSecondaries(0);
00111          aParticleChange.ProposeLocalEnergyDeposit (0.);
00112 #else
00113          aParticleChange.SetNumberOfSecondaries(0);
00114          aParticleChange.SetLocalEnergyDeposit (0.);
00115 #endif
00116 
00117        }
00118      //
00119      // Update the incident particle 
00120      //
00121 
00122     
00123      G4double NewKinEnergy=
00124        itsComptonEngine->GetScatteredElectron().e()-electron_mass_c2;
00125      
00126      //  G4double NewKinEnergy=0; // tmp to track photon only
00127      
00128      G4LorentzVector ScatEl=itsComptonEngine->GetScatteredElectron();
00129      
00130      
00131      if (NewKinEnergy > 0.)
00132        {
00133 #if G4VERSION_NUMBER > 699
00134          aParticleChange.ProposeMomentumDirection(ScatEl.vect().unit());
00135          aParticleChange.ProposeEnergy(NewKinEnergy);
00136          aParticleChange.ProposeLocalEnergyDeposit (0.); 
00137 #else
00138          aParticleChange.SetMomentumChange(ScatEl.vect().unit());
00139          aParticleChange.SetEnergyChange(NewKinEnergy);
00140          aParticleChange.SetLocalEnergyDeposit (0.); 
00141 #endif
00142        } 
00143      else
00144        { 
00145 
00146 #if G4VERSION_NUMBER > 699
00147          aParticleChange.ProposeEnergy( 0. );
00148          aParticleChange.ProposeLocalEnergyDeposit (0.);
00149          G4double charge= aDynamicParticle->GetCharge();
00150          if (charge<0.) aParticleChange.ProposeTrackStatus(fStopAndKill);
00151          else       aParticleChange.ProposeTrackStatus(fStopButAlive);
00152 #else
00153          aParticleChange.SetEnergyChange( 0. );
00154          aParticleChange.SetLocalEnergyDeposit (0.);
00155          G4double charge= aDynamicParticle->GetCharge();
00156          if (charge<0.) aParticleChange.SetStatusChange(fStopAndKill);
00157          else       aParticleChange.SetStatusChange(fStopButAlive);
00158 #endif
00159 
00160        }    
00161      
00162    }
00163  
00164  FireLaserCompton=false;
00165  
00166  #if G4VERSION_NUMBER > 899
00167  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
00168  #else
00169  return G4VContinuousDiscreteProcess::PostStepDoIt(trackData,stepData);
00170  #endif
00171 }
00172 
00173 #if G4VERSION_NUMBER > 899
00174 /*
00175 #include "G4LossTableManager.hh"
00176 #include "G4eBremsstrahlungModel.hh"
00177 #include "G4UniversalFluctuation.hh"
00178 
00179 void BDSLaserCompton::InitialiseEnergyLossProcess(const G4ParticleDefinition* p, const G4ParticleDefinition*)
00180 {
00181   if(!isInitialised) {
00182     particle = p;
00183     SetSecondaryParticle(G4Gamma::Gamma());
00184     SetIonisation(false);
00185     if (!EmModel()) SetEmModel(new G4eBremsstrahlungModel());
00186     EmModel()->SetLowEnergyLimit (100*eV);
00187     EmModel()->SetHighEnergyLimit(100*TeV);
00188     if (!FluctModel()) SetFluctModel(new G4UniversalFluctuation());
00189 
00190     AddEmModel(1, EmModel(), FluctModel());
00191     isInitialised = true;
00192   }
00193   G4LossTableManager* man = G4LossTableManager::Instance();
00194   dynamic_cast<G4eBremsstrahlungModel*>(EmModel())
00195     ->SetEnergyThreshold(man->BremsstrahlungTh());
00196   dynamic_cast<G4eBremsstrahlungModel*>(EmModel())
00197     ->SetLPMflag(man->LPMFlag());
00198 }
00199 
00200 void BDSLaserCompton::PrintInfo()
00201 {
00202 
00203   if(EmModel()) {
00204     G4cout << "      Total cross sections and sampling from "
00205            << EmModel()->GetName() << " model"
00206            << " (based on the EEDL data library) "
00207            << "\n      Good description from 1 KeV to 100 GeV, "
00208            << "log scale extrapolation above 100 GeV."
00209            << " LPM flag "
00210            << dynamic_cast<G4eBremsstrahlungModel*>(EmModel())->LPMflag()
00211            << G4endl;
00212     G4double eth = dynamic_cast<G4eBremsstrahlungModel*>(EmModel())->EnergyThreshold();
00213     if(eth < DBL_MIN)
00214       G4cout << "      HighEnergyThreshold(GeV)= " << eth/GeV
00215              << G4endl;
00216   }
00217 
00218 }
00219 */
00220 #endif

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7