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 
15  : thePhotonEnergy(0),
16  thenParticleInGun(0),
17  thenParticle(0),
18  theLaserBeamsInTEC1(),
19  theLaserBeamsInTEC2(),
20  theLaserBeamsInTECTIBTOBTEC()
21 {
22  // {{{ LaserPrimaryGeneratorAction constructor
23 
24  // get the PhotonEnergy from the parameter set
25  thePhotonEnergy = theConf.getUntrackedParameter<double>("PhotonEnergy",1.15) * eV;
26 
27  // number of particles in the Laser beam
28  thenParticleInGun = theConf.getUntrackedParameter<int>("NumberOfPhotonsInParticleGun",1);
29 
30  // number of particles in one beam. ATTENTION: each beam contains nParticleInGun with the same
31  // startpoint and direction. nParticle gives the number of particles in the beam with a different
32  // startpoint. They are used to simulate the gaussian beamprofile of the Laser Beams.
33  thenParticle = theConf.getUntrackedParameter<int>("NumberOfPhotonsInEachBeam",1);
34 
35  // create a messenger for this class
36 // theGunMessenger = new LaserPrimaryGeneratorMessenger(this);
37 
38  // create the beams in the right endcap
40 
41  // create the beams in the left endcap
43 
44  // create the beams to connect the TECs with TOB and TIB
46  // }}}
47 }
48 
50 {
51  // {{{ LaserPrimaryGeneratorAction destructor
52 
53  if ( theLaserBeamsInTEC1 != 0 ) { delete theLaserBeamsInTEC1; }
54  if ( theLaserBeamsInTEC2 != 0 ) { delete theLaserBeamsInTEC2; }
56  // }}}
57 }
58 
60 {
61  // {{{ GeneratePrimaries (G4Event * myEvent)
62 
63  // this function is called at the beginning of an Event in LaserAlignment::upDate(const BeginOfEvent * myEvent)
64  LogDebug("LaserPrimaryGeneratorAction") << "<LaserPrimaryGeneratorAction::GeneratePrimaries(G4Event*)>: create a new Laser Event";
65 
66  // shoot in the right endcap
68 
69  // shoot in the left endcap
71 
72  // shoot in the barrel
74 
75  // loop over all the generated vertices, get the primaries and set the user information
76  int theID = 0;
77 
78  for (int i = 1; i < myEvent->GetNumberOfPrimaryVertex(); i++)
79  {
80  G4PrimaryVertex * theVertex = myEvent->GetPrimaryVertex(i);
81 
82  for (int j = 0; j < theVertex->GetNumberOfParticle(); j++)
83  {
84  G4PrimaryParticle * thePrimary = theVertex->GetPrimary(j);
85 
86  setGeneratorId(thePrimary, theID);
87  theID++;
88  }
89  }
90  // }}}
91 }
92 
93 void LaserPrimaryGeneratorAction::setGeneratorId(G4PrimaryParticle * aParticle, int ID) const
94 {
95  // {{{ SetGeneratorId(G4PrimaryParticle * aParticle, int ID) const
96 
97  /* *********************************************************************** */
98  /* OSCAR expacts each G4PrimaryParticle to have some User Information *
99  * therefore this function have been implemented */
100  /* *********************************************************************** */
101 
102  aParticle->SetUserInformation(new GenParticleInfo(ID));
103  // }}}
104 }
#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