CMS 3D CMS Logo

LaserBeamsTEC2.cc
Go to the documentation of this file.
1 
10 
14 
15 #include "CLHEP/Random/RandGaussQ.h"
16 #include "CLHEP/Units/GlobalSystemOfUnits.h"
17 #include "G4ParticleDefinition.hh"
18 #include "G4ParticleGun.hh"
19 #include "globals.hh" // Global Constants and typedefs
20 
21 LaserBeamsTEC2::LaserBeamsTEC2() : theParticleGun(nullptr), theDRand48Engine(nullptr) {
22  G4int nPhotonsGun = 1;
23  G4int nPhotonsBeam = 1;
24  G4double Energy = 1.15 * eV;
25  // call constructor with options
26  LaserBeamsTEC2(nPhotonsGun, nPhotonsBeam, Energy);
27 }
28 
29 LaserBeamsTEC2::LaserBeamsTEC2(G4int nPhotonsInGun, G4int nPhotonsInBeam, G4double PhotonEnergy)
30  : thenParticleInGun(0), thenParticle(0), thePhotonEnergy(0), theParticleGun(), theDRand48Engine() {
31  /* *********************************************************************** */
32  /* initialize and configure the particle gun */
33  /* *********************************************************************** */
34 
35  // the Photon energy
36  thePhotonEnergy = PhotonEnergy;
37 
38  // number of particles in the Laser beam
39  thenParticleInGun = nPhotonsInGun;
40 
41  // number of particles in one beam. ATTENTION: each beam contains
42  // nParticleInGun with the same startpoint and direction. nParticle gives the
43  // number of particles in the beam with a different startpoint. They are used
44  // to simulate the gaussian beamprofile of the Laser Beams.
45  thenParticle = nPhotonsInBeam;
46 
47  // create the particle gun
48  theParticleGun = new G4ParticleGun(thenParticleInGun);
49 
50  // default kinematics
51  G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
52  G4ParticleDefinition *theOpticalPhoton = theParticleTable->FindParticle("opticalphoton");
53 
54  theParticleGun->SetParticleDefinition(theOpticalPhoton);
55  theParticleGun->SetParticleTime(0.0 * ns);
56  theParticleGun->SetParticlePosition(G4ThreeVector(-500.0 * cm, 0.0 * cm, 0.0 * cm));
57  theParticleGun->SetParticleMomentumDirection(G4ThreeVector(5.0, 3.0, 0.0));
58  theParticleGun->SetParticleEnergy(10.0 * keV);
59  setOptPhotonPolar(90.0);
60 
61  // initialize the random number engine
62  theDRand48Engine = new CLHEP::DRand48Engine();
63 }
64 
66  if (theParticleGun != nullptr) {
67  delete theParticleGun;
68  }
69  if (theDRand48Engine != nullptr) {
70  delete theDRand48Engine;
71  }
72 }
73 
75  // this function is called at the beginning of an Event in
76  // LaserAlignment::upDate(const BeginOfEvent * myEvent)
77 
78  // use the random number generator service of the framework
80  unsigned int seed = rng->mySeed();
81 
82  // set the seed
83  theDRand48Engine->setSeed(seed);
84 
85  // number of LaserRings and Laserdiodes
86  const G4int nLaserRings = 2;
87  const G4int nLaserBeams = 8;
88 
89  // z position of the sixth Tracker Endcap Disc, where the Laserdiodes are
90  // positioned
91  G4double LaserPositionZ = -2057.5 * mm;
92 
93  // Radius of the inner and outer Laser ring
94  G4double LaserRingRadius[nLaserRings] = {564.0 * mm, 840.0 * mm};
95 
96  // phi position of the first Laserdiode
97  G4double LaserPhi0 = 0.392699082;
98 
99  // width of the LaserBeams
100  G4double LaserRingSigmaX[nLaserRings] = {0.5 * mm, 0.5 * mm};
101  G4double LaserRingSigmaY[nLaserRings] = {0.5 * mm, 0.5 * mm};
102 
103  // get the definition of the optical photon
104  G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
105  G4ParticleDefinition *theOpticalPhoton = theParticleTable->FindParticle("opticalphoton");
106 
107  // loop over the LaserRings
108  for (int theRing = 0; theRing < nLaserRings; theRing++) {
109  // loop over the LaserBeams
110  for (int theBeam = 0; theBeam < nLaserBeams; theBeam++) {
111  // code for forward and backward beam
112  // calculate the right phi for the current beam
113  G4double LaserPositionPhi = LaserPhi0 + G4double(theBeam * G4double(G4double(2 * M_PI) / nLaserBeams));
114 
115  // calculate x and y position for the current beam
116  G4double LaserPositionX = cos(LaserPositionPhi) * LaserRingRadius[theRing];
117  G4double LaserPositionY = sin(LaserPositionPhi) * LaserRingRadius[theRing];
118 
119  // loop over all the particles in one beam
120  for (int theParticle = 0; theParticle < thenParticle; theParticle++) {
121  // get randomnumbers and calculate the position
122  CLHEP::RandGaussQ aGaussObjX(*theDRand48Engine, LaserPositionX, LaserRingSigmaX[theRing]);
123  CLHEP::RandGaussQ aGaussObjY(*theDRand48Engine, LaserPositionY, LaserRingSigmaY[theRing]);
124 
125  G4double theXPosition = aGaussObjX.fire();
126  G4double theYPosition = aGaussObjY.fire();
127  G4double theZPosition = LaserPositionZ;
128 
129  // set the properties of the newly created particle
130  theParticleGun->SetParticleDefinition(theOpticalPhoton);
131  theParticleGun->SetParticleTime(0.0 * ns);
132  theParticleGun->SetParticlePosition(G4ThreeVector(theXPosition, theYPosition, theZPosition));
133  theParticleGun->SetParticleEnergy(thePhotonEnergy);
134 
135  // loop over both directions of the beam
136  for (int theDirection = 0; theDirection < 2; theDirection++) {
137  // shoot in both beam directions ...
138  if (theDirection == 0) // shoot in forward direction (+z)
139  {
140  theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, 1.0));
141  // set the polarization
142  setOptPhotonPolar(90.0);
143  // generate the particle
144  theParticleGun->GeneratePrimaryVertex(myEvent);
145  }
146 
147  if (theDirection == 1) // shoot in backward direction (-z)
148  {
149  theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, -1.0));
150  // set the polarization
151  setOptPhotonPolar(90.0);
152  // generate the particle
153  theParticleGun->GeneratePrimaryVertex(myEvent);
154  }
155  } // end loop over both beam directions
156  } // end loop over particles in beam
157  } // end loop over beams
158  } // end loop over rings
159 }
160 
161 void LaserBeamsTEC2::setOptPhotonPolar(G4double Angle) {
162  /* *********************************************************************** */
163  /* to get optical processes working properly, you have to make sure *
164  * that the photon polarisation is defined. */
165  /* *********************************************************************** */
166 
167  // first check if we have an optical photon
168  if (theParticleGun->GetParticleDefinition()->GetParticleName() != "opticalphoton") {
169  edm::LogWarning("SimLaserAlignment:LaserBeamsTEC2")
170  << "<LaserBeamsTEC2::setOptPhotonPolar()>: WARNING! The ParticleGun is "
171  "not an optical photon";
172  return;
173  }
174 
175  // G4cout << " AC1CMS: The ParticleGun is an " <<
176  // theParticleGun->GetParticleDefinition()->GetParticleName();
177  G4ThreeVector normal(1.0, 0.0, 0.0);
178  G4ThreeVector kphoton = theParticleGun->GetParticleMomentumDirection();
179  G4ThreeVector product = normal.cross(kphoton);
180  G4double modul2 = product * product;
181 
182  G4ThreeVector e_perpendicular(0.0, 0.0, 1.0);
183 
184  if (modul2 > 0.0) {
185  e_perpendicular = (1.0 / sqrt(modul2)) * product;
186  }
187 
188  G4ThreeVector e_parallel = e_perpendicular.cross(kphoton);
189 
190  G4ThreeVector polar = cos(Angle) * e_parallel + sin(Angle) * e_perpendicular;
191 
192  // G4cout << ", the polarization = " << polar << G4endl;
193  theParticleGun->SetParticlePolarization(polar);
194 }
void setOptPhotonPolar(G4double Angle)
set the polarisation of the photons
G4ParticleGun * theParticleGun
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
LaserBeamsTEC2()
default constructor
G4int thenParticleInGun
virtual std::uint32_t mySeed() const =0
T sqrt(T t)
Definition: SSEVec.h:19
CLHEP::DRand48Engine * theDRand48Engine
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
~LaserBeamsTEC2() override
destructor
#define M_PI
G4double thePhotonEnergy
void GeneratePrimaries(G4Event *myEvent) override
shoot optical photons into the detector at the beginning of an event
Log< level::Warning, false > LogWarning