Go to the documentation of this file.00001
00009 #include "Alignment/LaserAlignmentSimulation/interface/LaserBeamsTEC1.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"
00017 #include "G4ParticleDefinition.hh"
00018 #include "G4ParticleGun.hh"
00019
00020 LaserBeamsTEC1::LaserBeamsTEC1()
00021 {
00022 G4int nPhotonsGun = 1;
00023 G4int nPhotonsBeam = 1;
00024 G4double Energy = 1.15 * eV;
00025
00026 LaserBeamsTEC1(nPhotonsGun, nPhotonsBeam, Energy);
00027 }
00028
00029 LaserBeamsTEC1::LaserBeamsTEC1(G4int nPhotonsInGun, G4int nPhotonsInBeam, G4double PhotonEnergy) : thenParticleInGun(0),
00030 thenParticle(0),
00031 thePhotonEnergy(0),
00032 theParticleGun(),
00033 theDRand48Engine()
00034 {
00035
00036
00037
00038
00039
00040 thePhotonEnergy = PhotonEnergy;
00041
00042
00043 thenParticleInGun = nPhotonsInGun;
00044
00045
00046
00047
00048 thenParticle = nPhotonsInBeam;
00049
00050
00051 theParticleGun = new G4ParticleGun(thenParticleInGun);
00052
00053
00054 G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
00055 G4ParticleDefinition * theOpticalPhoton = theParticleTable->FindParticle("opticalphoton");
00056
00057 theParticleGun->SetParticleDefinition(theOpticalPhoton);
00058 theParticleGun->SetParticleTime(0.0 * ns);
00059 theParticleGun->SetParticlePosition(G4ThreeVector(-500.0 * cm, 0.0 * cm, 0.0 * cm));
00060 theParticleGun->SetParticleMomentumDirection(G4ThreeVector(5.0, 3.0, 0.0));
00061 theParticleGun->SetParticleEnergy(10.0 * keV);
00062 setOptPhotonPolar(90.0);
00063
00064
00065 theDRand48Engine = new CLHEP::DRand48Engine();
00066
00067 }
00068
00069 LaserBeamsTEC1::~LaserBeamsTEC1()
00070 {
00071 if ( theParticleGun != 0 ) { delete theParticleGun; }
00072 if ( theDRand48Engine != 0 ) { delete theDRand48Engine; }
00073 }
00074
00075 void LaserBeamsTEC1::GeneratePrimaries(G4Event* myEvent)
00076 {
00077
00078
00079
00080 edm::Service<edm::RandomNumberGenerator> rng;
00081 unsigned int seed = rng->mySeed();
00082
00083
00084 theDRand48Engine->setSeed(seed);
00085
00086
00087 const G4int nLaserRings = 2;
00088 const G4int nLaserBeams = 8;
00089
00090
00091 G4double LaserPositionZ = 2057.5 * mm;
00092
00093
00094 G4double LaserRingRadius[nLaserRings] = {564.0 * mm, 840.0 * mm};
00095
00096
00097 G4double LaserPhi0 = 0.392699082;
00098
00099
00100 G4double LaserRingSigmaX[nLaserRings] = {0.5 * mm, 0.5 * mm};
00101 G4double LaserRingSigmaY[nLaserRings] = {0.5 * mm, 0.5 * mm};
00102
00103
00104 G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
00105 G4ParticleDefinition * theOpticalPhoton= theParticleTable->FindParticle("opticalphoton");
00106
00107
00108 for (int theRing = 0; theRing < nLaserRings; theRing++)
00109 {
00110
00111 for (int theBeam = 0; theBeam < nLaserBeams; theBeam++)
00112 {
00113
00114
00115 G4double LaserPositionPhi = LaserPhi0 + G4double(theBeam * G4double(G4double(2 * M_PI) / nLaserBeams));
00116
00117
00118 G4double LaserPositionX = cos(LaserPositionPhi) * LaserRingRadius[theRing];
00119 G4double LaserPositionY = sin(LaserPositionPhi) * LaserRingRadius[theRing];
00120
00121
00122 for (int theParticle = 0; theParticle < thenParticle; theParticle++)
00123 {
00124
00125 CLHEP::RandGaussQ aGaussObjX( *theDRand48Engine, LaserPositionX, LaserRingSigmaX[theRing] );
00126 CLHEP::RandGaussQ aGaussObjY( *theDRand48Engine, LaserPositionY, LaserRingSigmaY[theRing] );
00127
00128 G4double theXPosition = aGaussObjX.fire();
00129 G4double theYPosition = aGaussObjY.fire();
00130 G4double theZPosition = LaserPositionZ;
00131
00132
00133 theParticleGun->SetParticleDefinition(theOpticalPhoton);
00134 theParticleGun->SetParticleTime(0.0 * ns);
00135 theParticleGun->SetParticlePosition(G4ThreeVector(theXPosition, theYPosition, theZPosition));
00136 theParticleGun->SetParticleEnergy(thePhotonEnergy);
00137
00138
00139 for (int theDirection = 0; theDirection < 2; theDirection++)
00140 {
00141
00142 if (theDirection == 0)
00143 {
00144 theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, 1.0));
00145
00146 setOptPhotonPolar(90.0);
00147
00148 theParticleGun->GeneratePrimaryVertex(myEvent);
00149 }
00150
00151 if (theDirection == 1)
00152 {
00153 theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, -1.0));
00154
00155 setOptPhotonPolar(90.0);
00156
00157 theParticleGun->GeneratePrimaryVertex(myEvent);
00158 }
00159 }
00160 }
00161 }
00162 }
00163 }
00164
00165 void LaserBeamsTEC1::setOptPhotonPolar(G4double Angle)
00166 {
00167
00168
00169
00170
00171
00172
00173 if ( theParticleGun->GetParticleDefinition()->GetParticleName() != "opticalphoton" )
00174 {
00175 edm::LogWarning("SimLaserAlignment:LaserBeamsTEC1") << "<LaserBeamsTEC1::SetOptPhotonPolar()>: WARNING! The ParticleGun is not an optical photon";
00176 return;
00177 }
00178
00179
00180 G4ThreeVector normal(1.0, 0.0, 0.0);
00181 G4ThreeVector kphoton = theParticleGun->GetParticleMomentumDirection();
00182 G4ThreeVector product = normal.cross(kphoton);
00183 G4double modul2 = product * product;
00184
00185 G4ThreeVector e_perpendicular(0.0, 0.0, 1.0);
00186
00187 if ( modul2 > 0.0 ) { e_perpendicular = (1.0 / sqrt(modul2)) * product; }
00188
00189 G4ThreeVector e_parallel = e_perpendicular.cross(kphoton);
00190
00191 G4ThreeVector polar = cos(Angle) * e_parallel + sin(Angle) * e_perpendicular;
00192
00193
00194 theParticleGun->SetParticlePolarization(polar);
00195 }
00196