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

00001 /* BDSIM code.    Version 1.0
00002    Author: Grahame A. Blair, Royal Holloway, Univ. of London.
00003    Last modified 28.12.2004
00004    Copyright (c) 2004 by G.A.Blair.  ALL RIGHTS RESERVED. 
00005 
00006    Modified 22.03.05 by J.C.Carter, Royal Holloway, Univ. of London.
00007    Added GABs SynchGen code
00008 */
00009 
00010 //==================================================================
00011 //==================================================================
00012 #include "BDSGlobalConstants.hh"
00013 #include "BDSPrimaryGeneratorAction.hh"
00014 #include "BDSBunch.hh"
00015 #include "BDSParticle.hh"
00016 #include "BDSDebug.hh"
00017 
00018 #include "G4Event.hh"
00019 #include "G4ParticleGun.hh"
00020 #include "G4ParticleDefinition.hh"
00021 
00022 BDSPrimaryGeneratorAction::BDSPrimaryGeneratorAction(BDSBunch* bdsBunchIn):
00023   weight(1),
00024   bdsBunch(bdsBunchIn)
00025 {
00026   particleGun  = new G4ParticleGun(1); // 1-particle gun
00027 
00028 #ifdef BDSDEBUG
00029   G4cout << "BDSPrimaryGeneratorAction.cc: Primary particle is " << BDSGlobalConstants::Instance()->GetParticleDefinition()->GetParticleName() << G4endl;
00030 #endif
00031   
00032 #ifdef BDSDEBUG
00033   G4cout << "Setting momentum..." << G4endl;
00034 #endif
00035   particleGun->SetParticleMomentumDirection(G4ThreeVector(0.,0.,1.));
00036   particleGun->SetParticlePosition(G4ThreeVector(0.*CLHEP::cm,0.*CLHEP::cm,0.*CLHEP::cm));
00037   particleGun->SetParticleEnergy(BDSGlobalConstants::Instance()->GetBeamKineticEnergy());
00038   particleGun->SetParticleTime(0);
00039 }
00040 
00041 //===================================================
00042 
00043 BDSPrimaryGeneratorAction::~BDSPrimaryGeneratorAction()
00044 {
00045   delete particleGun;
00046 }
00047 
00048 //===================================================
00049 
00050 void BDSPrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
00051 {
00052   //this function is called at the begining of event
00053 
00054   G4double x0=0.0, y0=0.0, z0=0.0, xp=0.0, yp=0.0, zp=0.0, t=0.0, E=0.0;
00055 
00056   particleGun->SetParticleDefinition(BDSGlobalConstants::Instance()->GetParticleDefinition());
00057   bdsBunch->GetNextParticle(x0,y0,z0,xp,yp,zp,t,E,weight); // get next starting point
00058   if(E==0) G4cout << __METHOD_NAME__ << "Particle energy is 0! This will not be tracked." << G4endl;
00059 #ifdef BDSDEBUG 
00060   G4cout << __METHOD_NAME__ 
00061          << x0 << " " << y0 << " " << z0 << " " 
00062          << xp << " " << yp << " " << zp << " " 
00063          << t  << " " << E  << " " << weight << G4endl;
00064 #endif
00065 
00066   G4ThreeVector LocalPos;
00067   G4ThreeVector LocalMomDir;
00068   
00069   G4ThreeVector PartMomDir(xp,yp,zp);
00070   G4ThreeVector PartPosition(x0,y0,z0);
00071   
00072   particleGun->SetParticlePosition(PartPosition);
00073   particleGun->SetParticleEnergy(E);
00074   particleGun->SetParticleMomentumDirection(PartMomDir);
00075   particleGun->SetParticleTime(t);
00076   
00077   particleGun->GeneratePrimaryVertex(anEvent);
00078 
00079   //Set the weight
00080 #ifdef BDSDEBUG
00081   G4cout << "BDSPrimaryGeneratorAction: setting weight = " << weight << G4endl;
00082 #endif
00083   anEvent->GetPrimaryVertex()->SetWeight(weight);
00084   
00085   if(BDSGlobalConstants::Instance()->holdingQueue.size()!=0){
00086     
00087     //    anEvent->    GetTrack()->SetTrackID(outputParticle.trackID);
00088     //    anEvent->    GetTrack()->SetParentID(outputParticle.parentID);
00089     
00090     BDSGlobalConstants::Instance()->holdingQueue.pop_front();
00091     BDSGlobalConstants::Instance()->outputQueue.pop_front();
00092   }
00093 
00094   G4double totalE = E+particleGun->GetParticleDefinition()->GetPDGMass();
00095 #ifdef BDSDEBUG
00096   G4cout
00097     << "BDSPrimaryGeneratorAction: " << G4endl
00098     << "  position= " << particleGun->GetParticlePosition()/CLHEP::m<<" m"<<G4endl
00099     << "  kinetic energy= " << E/CLHEP::GeV << " GeV" << G4endl
00100     << "  total energy= " << totalE/CLHEP::GeV << " GeV" << G4endl
00101     << "  momentum direction= " << PartMomDir << G4endl
00102     << "  weight= " << anEvent->GetPrimaryVertex()->GetWeight() << G4endl;
00103 #endif
00104 
00105   // save initial values outside scope for entry into the samplers:
00106   BDSParticle initialPoint(x0,y0,z0,xp,yp,zp,totalE,t,weight);
00107   BDSGlobalConstants::Instance()->SetInitialPoint(initialPoint);
00108 }
00109 
00110 //===================================================

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7