CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
LaserPrimaryGeneratorAction Class Reference

#include <LaserPrimaryGeneratorAction.h>

Inheritance diagram for LaserPrimaryGeneratorAction:

Public Member Functions

void GeneratePrimaries (G4Event *myEvent) override
 
 LaserPrimaryGeneratorAction (edm::ParameterSet const &theConf)
 constructor More...
 
void setGeneratorId (G4PrimaryParticle *aParticle, int ID) const
 set Id of the optical photons More...
 
 ~LaserPrimaryGeneratorAction () override
 destructor More...
 

Private Attributes

LaserBeamsTEC1theLaserBeamsInTEC1
 
LaserBeamsTEC2theLaserBeamsInTEC2
 
LaserBeamsBarreltheLaserBeamsInTECTIBTOBTEC
 
G4int thenParticle
 
G4int thenParticleInGun
 
G4double thePhotonEnergy
 

Detailed Description

Primary Generator Action for the Laser Events

Date
2007/03/20 12:00:59
Revision
1.2
Author
Maarten Thomas

Definition at line 26 of file LaserPrimaryGeneratorAction.h.

Constructor & Destructor Documentation

LaserPrimaryGeneratorAction::LaserPrimaryGeneratorAction ( edm::ParameterSet const &  theConf)

constructor

Definition at line 16 of file LaserPrimaryGeneratorAction.cc.

References edm::ParameterSet::getUntrackedParameter(), theLaserBeamsInTEC1, theLaserBeamsInTEC2, theLaserBeamsInTECTIBTOBTEC, thenParticle, thenParticleInGun, and thePhotonEnergy.

17  : thePhotonEnergy(0),
19  thenParticle(0),
23  // {{{ LaserPrimaryGeneratorAction constructor
24 
25  // get the PhotonEnergy from the parameter set
26  thePhotonEnergy = theConf.getUntrackedParameter<double>("PhotonEnergy", 1.15) * eV;
27 
28  // number of particles in the Laser beam
29  thenParticleInGun = theConf.getUntrackedParameter<int>("NumberOfPhotonsInParticleGun", 1);
30 
31  // number of particles in one beam. ATTENTION: each beam contains
32  // nParticleInGun with the same startpoint and direction. nParticle gives the
33  // number of particles in the beam with a different startpoint. They are used
34  // 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 }
LaserPrimaryGeneratorAction::~LaserPrimaryGeneratorAction ( )
override

destructor

Definition at line 51 of file LaserPrimaryGeneratorAction.cc.

References theLaserBeamsInTEC1, theLaserBeamsInTEC2, and theLaserBeamsInTECTIBTOBTEC.

51  {
52  // {{{ LaserPrimaryGeneratorAction destructor
53 
54  if (theLaserBeamsInTEC1 != nullptr) {
55  delete theLaserBeamsInTEC1;
56  }
57  if (theLaserBeamsInTEC2 != nullptr) {
58  delete theLaserBeamsInTEC2;
59  }
60  if (theLaserBeamsInTECTIBTOBTEC != nullptr) {
62  }
63  // }}}
64 }

Member Function Documentation

void LaserPrimaryGeneratorAction::GeneratePrimaries ( G4Event *  myEvent)
override

call the corresponding GeneratePrimaries routines for both TEC's and the Barrel

Definition at line 66 of file LaserPrimaryGeneratorAction.cc.

References LaserBeamsTEC1::GeneratePrimaries(), LaserBeamsTEC2::GeneratePrimaries(), LaserBeamsBarrel::GeneratePrimaries(), mps_fire::i, dqmiolumiharvest::j, LogDebug, setGeneratorId(), theLaserBeamsInTEC1, theLaserBeamsInTEC2, and theLaserBeamsInTECTIBTOBTEC.

Referenced by LaserAlignmentSimulation::update().

66  {
67  // {{{ GeneratePrimaries (G4Event * myEvent)
68 
69  // this function is called at the beginning of an Event in
70  // LaserAlignment::upDate(const BeginOfEvent * myEvent)
71  LogDebug("LaserPrimaryGeneratorAction") << "<LaserPrimaryGeneratorAction::GeneratePrimaries(G4Event*)>: create a "
72  "new Laser Event";
73 
74  // shoot in the right endcap
76 
77  // shoot in the left endcap
79 
80  // shoot in the barrel
82 
83  // loop over all the generated vertices, get the primaries and set the user
84  // information
85  int theID = 0;
86 
87  for (int i = 1; i < myEvent->GetNumberOfPrimaryVertex(); i++) {
88  G4PrimaryVertex *theVertex = myEvent->GetPrimaryVertex(i);
89 
90  for (int j = 0; j < theVertex->GetNumberOfParticle(); j++) {
91  G4PrimaryParticle *thePrimary = theVertex->GetPrimary(j);
92 
93  setGeneratorId(thePrimary, theID);
94  theID++;
95  }
96  }
97  // }}}
98 }
#define LogDebug(id)
void GeneratePrimaries(G4Event *myEvent) override
shoot optical photons into the detector at the beginning of an event
void setGeneratorId(G4PrimaryParticle *aParticle, int ID) const
set Id of the optical photons
void GeneratePrimaries(G4Event *myEvent) override
shoot optical photons into the detector at the beginning of an event
void GeneratePrimaries(G4Event *myEvent) override
shoot optical photons into the detector at the beginning of an event
void LaserPrimaryGeneratorAction::setGeneratorId ( G4PrimaryParticle *  aParticle,
int  ID 
) const

set Id of the optical photons

Definition at line 100 of file LaserPrimaryGeneratorAction.cc.

Referenced by GeneratePrimaries().

100  {
101  // {{{ SetGeneratorId(G4PrimaryParticle * aParticle, int ID) const
102 
103  /* *********************************************************************** */
104  /* OSCAR expacts each G4PrimaryParticle to have some User Information *
105  * therefore this function have been implemented */
106  /* *********************************************************************** */
107 
108  aParticle->SetUserInformation(new GenParticleInfo(ID));
109  // }}}
110 }
uint32_t ID
Definition: Definitions.h:24

Member Data Documentation

LaserBeamsTEC1* LaserPrimaryGeneratorAction::theLaserBeamsInTEC1
private
LaserBeamsTEC2* LaserPrimaryGeneratorAction::theLaserBeamsInTEC2
private
LaserBeamsBarrel* LaserPrimaryGeneratorAction::theLaserBeamsInTECTIBTOBTEC
private
G4int LaserPrimaryGeneratorAction::thenParticle
private

Definition at line 43 of file LaserPrimaryGeneratorAction.h.

Referenced by LaserPrimaryGeneratorAction().

G4int LaserPrimaryGeneratorAction::thenParticleInGun
private

Definition at line 42 of file LaserPrimaryGeneratorAction.h.

Referenced by LaserPrimaryGeneratorAction().

G4double LaserPrimaryGeneratorAction::thePhotonEnergy
private

Definition at line 41 of file LaserPrimaryGeneratorAction.h.

Referenced by LaserPrimaryGeneratorAction().