Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-10 02:20:15

0001 /** \file LaserBeamsTEC1.cc
0002  *
0003  *
0004  *  $Date: 2010/09/09 18:22:48 $
0005  *  $Revision: 1.7 $
0006  *  \author Maarten Thomas
0007  */
0008 
0009 #include "Alignment/LaserAlignmentSimulation/interface/LaserBeamsTEC1.h"
0010 
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/ServiceRegistry/interface/Service.h"
0013 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
0014 
0015 #include <CLHEP/Random/RandGaussQ.h>
0016 #include <CLHEP/Units/SystemOfUnits.h>
0017 #include "G4ParticleDefinition.hh"
0018 #include "G4ParticleGun.hh"
0019 #include "globals.hh"  // Global Constants and typedefs
0020 
0021 LaserBeamsTEC1::LaserBeamsTEC1() : theParticleGun(nullptr), theDRand48Engine(nullptr) {
0022   G4int nPhotonsGun = 1;
0023   G4int nPhotonsBeam = 1;
0024   G4double Energy = 1.15 * CLHEP::eV;
0025   // call constructor with options
0026   LaserBeamsTEC1(nPhotonsGun, nPhotonsBeam, Energy);
0027 }
0028 
0029 LaserBeamsTEC1::LaserBeamsTEC1(G4int nPhotonsInGun, G4int nPhotonsInBeam, G4double PhotonEnergy)
0030     : thenParticleInGun(0), thenParticle(0), thePhotonEnergy(0) {
0031   /* *********************************************************************** */
0032   /*  initialize and configure the particle gun                              */
0033   /* *********************************************************************** */
0034 
0035   // the Photon energy
0036   thePhotonEnergy = PhotonEnergy;
0037 
0038   // number of particles in the Laser beam
0039   thenParticleInGun = nPhotonsInGun;
0040 
0041   // number of particles in one beam. ATTENTION: each beam contains
0042   // nParticleInGun with the same startpoint and direction. nParticle gives the
0043   // number of particles in the beam with a different startpoint. They are used
0044   // to simulate the gaussian beamprofile of the Laser Beams.
0045   thenParticle = nPhotonsInBeam;
0046 
0047   // create the particle gun
0048   theParticleGun = new G4ParticleGun(thenParticleInGun);
0049 
0050   // default kinematics
0051   G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
0052   G4ParticleDefinition *theOpticalPhoton = theParticleTable->FindParticle("opticalphoton");
0053 
0054   theParticleGun->SetParticleDefinition(theOpticalPhoton);
0055   theParticleGun->SetParticleTime(0.0 * CLHEP::ns);
0056   theParticleGun->SetParticlePosition(G4ThreeVector(-500.0 * CLHEP::cm, 0.0 * CLHEP::cm, 0.0 * CLHEP::cm));
0057   theParticleGun->SetParticleMomentumDirection(G4ThreeVector(5.0, 3.0, 0.0));
0058   theParticleGun->SetParticleEnergy(10.0 * CLHEP::keV);
0059   setOptPhotonPolar(90.0);
0060 
0061   // initialize the random number engine
0062   theDRand48Engine = new CLHEP::DRand48Engine();
0063 }
0064 
0065 LaserBeamsTEC1::~LaserBeamsTEC1() {
0066   if (theParticleGun != nullptr) {
0067     delete theParticleGun;
0068   }
0069   if (theDRand48Engine != nullptr) {
0070     delete theDRand48Engine;
0071   }
0072 }
0073 
0074 void LaserBeamsTEC1::GeneratePrimaries(G4Event *myEvent) {
0075   // this function is called at the beginning of an Event in
0076   // LaserAlignment::upDate(const BeginOfEvent * myEvent)
0077 
0078   // use the random number generator service of the framework
0079   edm::Service<edm::RandomNumberGenerator> rng;
0080   unsigned int seed = rng->mySeed();
0081 
0082   // set the seed
0083   theDRand48Engine->setSeed(seed);
0084 
0085   // number of LaserRings and Laserdiodes
0086   const G4int nLaserRings = 2;
0087   const G4int nLaserBeams = 8;
0088 
0089   // z position of the sixth Tracker Endcap Disc, where the Laserdiodes are
0090   // positioned
0091   using CLHEP::mm;
0092   G4double LaserPositionZ = 2057.5 * mm;
0093 
0094   // Radius of the inner and outer Laser ring
0095   G4double LaserRingRadius[nLaserRings] = {564.0 * mm, 840.0 * mm};
0096 
0097   // phi position of the first Laserdiode
0098   G4double LaserPhi0 = 0.392699082;
0099 
0100   // width of the LaserBeams
0101   G4double LaserRingSigmaX[nLaserRings] = {0.5 * mm, 0.5 * mm};
0102   G4double LaserRingSigmaY[nLaserRings] = {0.5 * mm, 0.5 * mm};
0103 
0104   // get the definition of the optical photon
0105   G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
0106   G4ParticleDefinition *theOpticalPhoton = theParticleTable->FindParticle("opticalphoton");
0107 
0108   // loop over the LaserRings
0109   for (int theRing = 0; theRing < nLaserRings; theRing++) {
0110     // loop over the LaserBeams
0111     for (int theBeam = 0; theBeam < nLaserBeams; theBeam++) {
0112       // code for forward and backward beam
0113       // calculate the right phi for the current beam
0114       G4double LaserPositionPhi = LaserPhi0 + G4double(theBeam * G4double(G4double(2 * M_PI) / nLaserBeams));
0115 
0116       // calculate x and y position for the current beam
0117       G4double LaserPositionX = cos(LaserPositionPhi) * LaserRingRadius[theRing];
0118       G4double LaserPositionY = sin(LaserPositionPhi) * LaserRingRadius[theRing];
0119 
0120       // loop over all the particles in one beam
0121       for (int theParticle = 0; theParticle < thenParticle; theParticle++) {
0122         // get randomnumbers  and calculate the position
0123         CLHEP::RandGaussQ aGaussObjX(*theDRand48Engine, LaserPositionX, LaserRingSigmaX[theRing]);
0124         CLHEP::RandGaussQ aGaussObjY(*theDRand48Engine, LaserPositionY, LaserRingSigmaY[theRing]);
0125 
0126         G4double theXPosition = aGaussObjX.fire();
0127         G4double theYPosition = aGaussObjY.fire();
0128         G4double theZPosition = LaserPositionZ;
0129 
0130         // set the properties of the newly created particle
0131         theParticleGun->SetParticleDefinition(theOpticalPhoton);
0132         theParticleGun->SetParticleTime(0.0 * CLHEP::ns);
0133         theParticleGun->SetParticlePosition(G4ThreeVector(theXPosition, theYPosition, theZPosition));
0134         theParticleGun->SetParticleEnergy(thePhotonEnergy);
0135 
0136         // loop over both directions of the beam
0137         for (int theDirection = 0; theDirection < 2; theDirection++) {
0138           // shoot in both beam directions ...
0139           if (theDirection == 0)  // shoot in forward direction (+z)
0140           {
0141             theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, 1.0));
0142             // set the polarization
0143             setOptPhotonPolar(90.0);
0144             // generate the particle
0145             theParticleGun->GeneratePrimaryVertex(myEvent);
0146           }
0147 
0148           if (theDirection == 1)  // shoot in backward direction (-z)
0149           {
0150             theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, -1.0));
0151             // set the polarization
0152             setOptPhotonPolar(90.0);
0153             // generate the particle
0154             theParticleGun->GeneratePrimaryVertex(myEvent);
0155           }
0156         }  // end loop over both beam directions
0157       }    // end loop over particles in beam
0158     }      // end loop over beams
0159   }        // end loop over rings
0160 }
0161 
0162 void LaserBeamsTEC1::setOptPhotonPolar(G4double Angle) {
0163   /* *********************************************************************** */
0164   /*   to get optical processes working properly, you have to make sure      *
0165    *   that the photon polarisation is defined.                              */
0166   /* *********************************************************************** */
0167 
0168   // first check if we have an optical photon
0169   if (theParticleGun->GetParticleDefinition()->GetParticleName() != "opticalphoton") {
0170     edm::LogWarning("SimLaserAlignment:LaserBeamsTEC1")
0171         << "<LaserBeamsTEC1::SetOptPhotonPolar()>: WARNING! The ParticleGun is "
0172            "not an optical photon";
0173     return;
0174   }
0175 
0176   //   G4cout << "  AC1CMS: The ParticleGun is an " <<
0177   //   theParticleGun->GetParticleDefinition()->GetParticleName();
0178   G4ThreeVector normal(1.0, 0.0, 0.0);
0179   G4ThreeVector kphoton = theParticleGun->GetParticleMomentumDirection();
0180   G4ThreeVector product = normal.cross(kphoton);
0181   G4double modul2 = product * product;
0182 
0183   G4ThreeVector e_perpendicular(0.0, 0.0, 1.0);
0184 
0185   if (modul2 > 0.0) {
0186     e_perpendicular = (1.0 / sqrt(modul2)) * product;
0187   }
0188 
0189   G4ThreeVector e_parallel = e_perpendicular.cross(kphoton);
0190 
0191   G4ThreeVector polar = cos(Angle) * e_parallel + sin(Angle) * e_perpendicular;
0192 
0193   //   G4cout << ", the polarization = " << polar << G4endl;
0194   theParticleGun->SetParticlePolarization(polar);
0195 }