CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LaserPrimaryGeneratorAction.cc
Go to the documentation of this file.
1 
11 
13 
14 #include "G4SystemOfUnits.hh"
15 
17  : thePhotonEnergy(0),
18  thenParticleInGun(0),
19  thenParticle(0),
20  theLaserBeamsInTEC1(),
21  theLaserBeamsInTEC2(),
22  theLaserBeamsInTECTIBTOBTEC()
23 {
24  // {{{ LaserPrimaryGeneratorAction constructor
25 
26  // get the PhotonEnergy from the parameter set
27  thePhotonEnergy = theConf.getUntrackedParameter<double>("PhotonEnergy",1.15) * eV;
28 
29  // number of particles in the Laser beam
30  thenParticleInGun = theConf.getUntrackedParameter<int>("NumberOfPhotonsInParticleGun",1);
31 
32  // number of particles in one beam. ATTENTION: each beam contains nParticleInGun with the same
33  // startpoint and direction. nParticle gives the number of particles in the beam with a different
34  // startpoint. They are used to simulate the gaussian beamprofile of the Laser Beams.
35  thenParticle = theConf.getUntrackedParameter<int>("NumberOfPhotonsInEachBeam",1);
36 
37  // create a messenger for this class
38 // theGunMessenger = new LaserPrimaryGeneratorMessenger(this);
39 
40  // create the beams in the right endcap
42 
43  // create the beams in the left endcap
45 
46  // create the beams to connect the TECs with TOB and TIB
48  // }}}
49 }
50 
52 {
53  // {{{ LaserPrimaryGeneratorAction destructor
54 
55  if ( theLaserBeamsInTEC1 != 0 ) { delete theLaserBeamsInTEC1; }
56  if ( theLaserBeamsInTEC2 != 0 ) { delete theLaserBeamsInTEC2; }
58  // }}}
59 }
60 
62 {
63  // {{{ GeneratePrimaries (G4Event * myEvent)
64 
65  // this function is called at the beginning of an Event in LaserAlignment::upDate(const BeginOfEvent * myEvent)
66  LogDebug("LaserPrimaryGeneratorAction") << "<LaserPrimaryGeneratorAction::GeneratePrimaries(G4Event*)>: create a new Laser Event";
67 
68  // shoot in the right endcap
70 
71  // shoot in the left endcap
73 
74  // shoot in the barrel
76 
77  // loop over all the generated vertices, get the primaries and set the user information
78  int theID = 0;
79 
80  for (int i = 1; i < myEvent->GetNumberOfPrimaryVertex(); i++)
81  {
82  G4PrimaryVertex * theVertex = myEvent->GetPrimaryVertex(i);
83 
84  for (int j = 0; j < theVertex->GetNumberOfParticle(); j++)
85  {
86  G4PrimaryParticle * thePrimary = theVertex->GetPrimary(j);
87 
88  setGeneratorId(thePrimary, theID);
89  theID++;
90  }
91  }
92  // }}}
93 }
94 
95 void LaserPrimaryGeneratorAction::setGeneratorId(G4PrimaryParticle * aParticle, int ID) const
96 {
97  // {{{ SetGeneratorId(G4PrimaryParticle * aParticle, int ID) const
98 
99  /* *********************************************************************** */
100  /* OSCAR expacts each G4PrimaryParticle to have some User Information *
101  * therefore this function have been implemented */
102  /* *********************************************************************** */
103 
104  aParticle->SetUserInformation(new GenParticleInfo(ID));
105  // }}}
106 }
#define LogDebug(id)
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void GeneratePrimaries(G4Event *myEvent)
shoot optical photons into the detector at the beginning of an event
uint32_t ID
Definition: Definitions.h:26
void GeneratePrimaries(G4Event *myEvent)
shoot optical photons into the detector at the beginning of an event
LaserPrimaryGeneratorAction(edm::ParameterSet const &theConf)
constructor
void GeneratePrimaries(G4Event *myEvent)
shoot optical photons into the detector at the beginning of an event
int j
Definition: DBlmapReader.cc:9
void setGeneratorId(G4PrimaryParticle *aParticle, int ID) const
set Id of the optical photons
void GeneratePrimaries(G4Event *myEvent)
call the corresponding GeneratePrimaries routines for both TEC&#39;s and the Barrel