CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/Alignment/LaserAlignmentSimulation/src/LaserBeamsTEC2.cc

Go to the documentation of this file.
00001 
00009 #include "Alignment/LaserAlignmentSimulation/interface/LaserBeamsTEC2.h"
00010 
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 #include "FWCore/ServiceRegistry/interface/Service.h"
00013 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00014 
00015 #include "CLHEP/Random/RandGaussQ.h"
00016 #include "globals.hh"                        // Global Constants and typedefs
00017 #include "G4ParticleDefinition.hh"
00018 #include "G4ParticleGun.hh"
00019 
00020 LaserBeamsTEC2::LaserBeamsTEC2() :
00021   theParticleGun(0),
00022   theDRand48Engine(0)
00023 { 
00024   G4int nPhotonsGun = 1;
00025   G4int nPhotonsBeam = 1;
00026   G4double Energy = 1.15 * eV;
00027   // call constructor with options
00028   LaserBeamsTEC2(nPhotonsGun, nPhotonsBeam, Energy);
00029 }
00030 
00031 LaserBeamsTEC2::LaserBeamsTEC2(G4int nPhotonsInGun, G4int nPhotonsInBeam, G4double PhotonEnergy) : thenParticleInGun(0),
00032                                                                                                    thenParticle(0),
00033                                                                                                    thePhotonEnergy(0),
00034                                                                                                    theParticleGun(),
00035                                                                                                    theDRand48Engine()
00036 {
00037   /* *********************************************************************** */
00038   /*  initialize and configure the particle gun                              */
00039   /* *********************************************************************** */
00040 
00041   // the Photon energy
00042   thePhotonEnergy = PhotonEnergy;
00043 
00044   // number of particles in the Laser beam
00045   thenParticleInGun = nPhotonsInGun;
00046 
00047   // number of particles in one beam. ATTENTION: each beam contains nParticleInGun with the same
00048   // startpoint and direction. nParticle gives the number of particles in the beam with a different
00049   // startpoint. They are used to simulate the gaussian beamprofile of the Laser Beams.
00050   thenParticle = nPhotonsInBeam;
00051 
00052   // create the particle gun
00053   theParticleGun = new G4ParticleGun(thenParticleInGun);
00054 
00055   // default kinematics
00056   G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
00057   G4ParticleDefinition * theOpticalPhoton = theParticleTable->FindParticle("opticalphoton");
00058 
00059   theParticleGun->SetParticleDefinition(theOpticalPhoton);
00060   theParticleGun->SetParticleTime(0.0 * ns);
00061   theParticleGun->SetParticlePosition(G4ThreeVector(-500.0 * cm, 0.0 * cm, 0.0 * cm));
00062   theParticleGun->SetParticleMomentumDirection(G4ThreeVector(5.0, 3.0, 0.0));
00063   theParticleGun->SetParticleEnergy(10.0 * keV);
00064   setOptPhotonPolar(90.0);
00065 
00066   // initialize the random number engine
00067   theDRand48Engine = new CLHEP::DRand48Engine();
00068 }
00069 
00070 LaserBeamsTEC2::~LaserBeamsTEC2()
00071 {
00072   if ( theParticleGun != 0 )  { delete theParticleGun; }
00073   if ( theDRand48Engine != 0 ) { delete theDRand48Engine; }
00074 }
00075 
00076 void LaserBeamsTEC2::GeneratePrimaries(G4Event* myEvent)
00077 {
00078   // this function is called at the beginning of an Event in LaserAlignment::upDate(const BeginOfEvent * myEvent)
00079 
00080   // use the random number generator service of the framework
00081   edm::Service<edm::RandomNumberGenerator> rng;
00082   unsigned int seed = rng->mySeed();
00083 
00084   // set the seed
00085   theDRand48Engine->setSeed(seed);
00086 
00087 
00088   // number of LaserRings and Laserdiodes
00089   const G4int nLaserRings = 2;
00090   const G4int nLaserBeams = 8;
00091 
00092   // z position of the sixth Tracker Endcap Disc, where the Laserdiodes are positioned
00093   G4double LaserPositionZ = -2057.5 * mm;
00094 
00095   // Radius of the inner and outer Laser ring
00096   G4double LaserRingRadius[nLaserRings] = {564.0 * mm, 840.0 * mm};
00097 
00098   // phi position of the first Laserdiode
00099   G4double LaserPhi0 = 0.392699082;
00100 
00101   // width of the LaserBeams
00102   G4double LaserRingSigmaX[nLaserRings] = {0.5 * mm, 0.5 * mm};
00103   G4double LaserRingSigmaY[nLaserRings] = {0.5 * mm, 0.5 * mm};
00104 
00105   // get the definition of the optical photon
00106   G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
00107   G4ParticleDefinition * theOpticalPhoton= theParticleTable->FindParticle("opticalphoton");
00108 
00109   // loop over the LaserRings
00110   for (int theRing = 0; theRing < nLaserRings; theRing++)
00111     {
00112       // loop over the LaserBeams
00113       for (int theBeam = 0; theBeam < nLaserBeams; theBeam++)
00114         {
00115           // code for forward and backward beam
00116           // calculate the right phi for the current beam
00117           G4double LaserPositionPhi = LaserPhi0 + G4double(theBeam * G4double(G4double(2 * M_PI) / nLaserBeams));
00118 
00119           // calculate x and y position for the current beam
00120           G4double LaserPositionX = cos(LaserPositionPhi) * LaserRingRadius[theRing];
00121           G4double LaserPositionY = sin(LaserPositionPhi) * LaserRingRadius[theRing];
00122 
00123           // loop over all the particles in one beam
00124           for (int theParticle = 0; theParticle < thenParticle; theParticle++)
00125             {
00126               // get randomnumbers  and calculate the position
00127               CLHEP::RandGaussQ aGaussObjX( *theDRand48Engine, LaserPositionX, LaserRingSigmaX[theRing] );
00128               CLHEP::RandGaussQ aGaussObjY( *theDRand48Engine, LaserPositionY, LaserRingSigmaY[theRing] );
00129 
00130               G4double theXPosition = aGaussObjX.fire();
00131               G4double theYPosition = aGaussObjY.fire();
00132               G4double theZPosition = LaserPositionZ;
00133 
00134               // set the properties of the newly created particle
00135               theParticleGun->SetParticleDefinition(theOpticalPhoton);
00136               theParticleGun->SetParticleTime(0.0 * ns);
00137               theParticleGun->SetParticlePosition(G4ThreeVector(theXPosition, theYPosition, theZPosition));
00138               theParticleGun->SetParticleEnergy(thePhotonEnergy);
00139 
00140               // loop over both directions of the beam
00141               for (int theDirection = 0; theDirection < 2; theDirection++)
00142                 {
00143                   // shoot in both beam directions ...
00144                   if (theDirection == 0)  // shoot in forward direction (+z)
00145                     {
00146                       theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, 1.0));
00147                       // set the polarization
00148                       setOptPhotonPolar(90.0);
00149                       // generate the particle
00150                       theParticleGun->GeneratePrimaryVertex(myEvent);
00151                     }
00152 
00153                   if (theDirection == 1)  // shoot in backward direction (-z)
00154                     {
00155                       theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, -1.0));
00156                       // set the polarization
00157                       setOptPhotonPolar(90.0);
00158                       // generate the particle
00159                       theParticleGun->GeneratePrimaryVertex(myEvent);
00160                     }
00161                 } // end loop over both beam directions
00162             } // end loop over particles in beam
00163         } // end loop over beams
00164     } // end loop over rings
00165 }
00166 
00167 void LaserBeamsTEC2::setOptPhotonPolar(G4double Angle)
00168 {
00169   /* *********************************************************************** */
00170   /*   to get optical processes working properly, you have to make sure      *
00171    *   that the photon polarisation is defined.                              */
00172   /* *********************************************************************** */
00173 
00174   // first check if we have an optical photon
00175   if ( theParticleGun->GetParticleDefinition()->GetParticleName() != "opticalphoton" )
00176     { 
00177       edm::LogWarning("SimLaserAlignment:LaserBeamsTEC2") << "<LaserBeamsTEC2::setOptPhotonPolar()>: WARNING! The ParticleGun is not an optical photon";
00178       return;
00179     }
00180 
00181 //   G4cout << "  AC1CMS: The ParticleGun is an " << theParticleGun->GetParticleDefinition()->GetParticleName();
00182   G4ThreeVector normal(1.0, 0.0, 0.0);
00183   G4ThreeVector kphoton = theParticleGun->GetParticleMomentumDirection();
00184   G4ThreeVector product = normal.cross(kphoton);
00185   G4double modul2 = product * product;
00186 
00187   G4ThreeVector e_perpendicular(0.0, 0.0, 1.0);
00188   
00189   if ( modul2 > 0.0 ) { e_perpendicular = (1.0 / sqrt(modul2)) * product; }
00190   
00191   G4ThreeVector e_parallel = e_perpendicular.cross(kphoton);
00192 
00193   G4ThreeVector polar = cos(Angle) * e_parallel + sin(Angle) * e_perpendicular;
00194   
00195 //   G4cout << ", the polarization = " << polar << G4endl;
00196   theParticleGun->SetParticlePolarization(polar);
00197 }
00198