CMS 3D CMS Logo

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