Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:32:09

0001 /** \file LaserPrimaryGeneratorAction.cc
0002  *
0003  *
0004  *  $Date: 2007/03/20 12:01:01 $
0005  *  $Revision: 1.2 $
0006  *  \author Maarten Thomas
0007  */
0008 
0009 #include "Alignment/LaserAlignmentSimulation/interface/LaserPrimaryGeneratorAction.h"
0010 #include "SimG4Core/Notification/interface/GenParticleInfo.h"
0011 
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 
0014 #include "G4SystemOfUnits.hh"
0015 
0016 LaserPrimaryGeneratorAction::LaserPrimaryGeneratorAction(edm::ParameterSet const &theConf)
0017     : thePhotonEnergy(0),
0018       thenParticleInGun(0),
0019       thenParticle(0),
0020       theLaserBeamsInTEC1(),
0021       theLaserBeamsInTEC2(),
0022       theLaserBeamsInTECTIBTOBTEC() {
0023   // {{{ LaserPrimaryGeneratorAction constructor
0024 
0025   // get the PhotonEnergy from the parameter set
0026   thePhotonEnergy = theConf.getUntrackedParameter<double>("PhotonEnergy", 1.15) * eV;
0027 
0028   // number of particles in the Laser beam
0029   thenParticleInGun = theConf.getUntrackedParameter<int>("NumberOfPhotonsInParticleGun", 1);
0030 
0031   // number of particles in one beam. ATTENTION: each beam contains
0032   // nParticleInGun with the same startpoint and direction. nParticle gives the
0033   // number of particles in the beam with a different startpoint. They are used
0034   // to simulate the gaussian beamprofile of the Laser Beams.
0035   thenParticle = theConf.getUntrackedParameter<int>("NumberOfPhotonsInEachBeam", 1);
0036 
0037   // create a messenger for this class
0038   //   theGunMessenger = new LaserPrimaryGeneratorMessenger(this);
0039 
0040   // create the beams in the right endcap
0041   theLaserBeamsInTEC1 = new LaserBeamsTEC1(thenParticleInGun, thenParticle, thePhotonEnergy);
0042 
0043   // create the beams in the left endcap
0044   theLaserBeamsInTEC2 = new LaserBeamsTEC2(thenParticleInGun, thenParticle, thePhotonEnergy);
0045 
0046   // create the beams to connect the TECs with TOB and TIB
0047   theLaserBeamsInTECTIBTOBTEC = new LaserBeamsBarrel(thenParticleInGun, thenParticle, thePhotonEnergy);
0048   // }}}
0049 }
0050 
0051 LaserPrimaryGeneratorAction::~LaserPrimaryGeneratorAction() {
0052   // {{{ LaserPrimaryGeneratorAction destructor
0053 
0054   if (theLaserBeamsInTEC1 != nullptr) {
0055     delete theLaserBeamsInTEC1;
0056   }
0057   if (theLaserBeamsInTEC2 != nullptr) {
0058     delete theLaserBeamsInTEC2;
0059   }
0060   if (theLaserBeamsInTECTIBTOBTEC != nullptr) {
0061     delete theLaserBeamsInTECTIBTOBTEC;
0062   }
0063   // }}}
0064 }
0065 
0066 void LaserPrimaryGeneratorAction::GeneratePrimaries(G4Event *myEvent) {
0067   // {{{ GeneratePrimaries (G4Event * myEvent)
0068 
0069   // this function is called at the beginning of an Event in
0070   // LaserAlignment::upDate(const BeginOfEvent * myEvent)
0071   LogDebug("LaserPrimaryGeneratorAction") << "<LaserPrimaryGeneratorAction::GeneratePrimaries(G4Event*)>: create a "
0072                                              "new Laser Event";
0073 
0074   // shoot in the right endcap
0075   theLaserBeamsInTEC1->GeneratePrimaries(myEvent);
0076 
0077   // shoot in the left endcap
0078   theLaserBeamsInTEC2->GeneratePrimaries(myEvent);
0079 
0080   // shoot in the barrel
0081   theLaserBeamsInTECTIBTOBTEC->GeneratePrimaries(myEvent);
0082 
0083   // loop over all the generated vertices, get the primaries and set the user
0084   // information
0085   int theID = 0;
0086 
0087   for (int i = 1; i < myEvent->GetNumberOfPrimaryVertex(); i++) {
0088     G4PrimaryVertex *theVertex = myEvent->GetPrimaryVertex(i);
0089 
0090     for (int j = 0; j < theVertex->GetNumberOfParticle(); j++) {
0091       G4PrimaryParticle *thePrimary = theVertex->GetPrimary(j);
0092 
0093       setGeneratorId(thePrimary, theID);
0094       theID++;
0095     }
0096   }
0097   // }}}
0098 }
0099 
0100 void LaserPrimaryGeneratorAction::setGeneratorId(G4PrimaryParticle *aParticle, int ID) const {
0101   // {{{ SetGeneratorId(G4PrimaryParticle * aParticle, int ID) const
0102 
0103   /* *********************************************************************** */
0104   /*   OSCAR expacts each G4PrimaryParticle to have some User Information    *
0105    *   therefore this function have been implemented                         */
0106   /* *********************************************************************** */
0107 
0108   aParticle->SetUserInformation(new GenParticleInfo(ID));
0109   // }}}
0110 }