CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LaserBeamsBarrel.cc
Go to the documentation of this file.
1 
10 
14 
15 #include "CLHEP/Random/RandGaussQ.h"
16 #include "globals.hh" // Global Constants and typedefs
17 #include "G4ParticleDefinition.hh"
18 #include "G4ParticleGun.hh"
19 #include "G4SystemOfUnits.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  LaserBeamsBarrel(nPhotonsGun, nPhotonsBeam, Energy);
30 }
31 
32 LaserBeamsBarrel::LaserBeamsBarrel(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 }
71 
73 {
74  if ( theParticleGun != 0 ) { delete theParticleGun; }
75  if ( theDRand48Engine != 0 ) { delete theDRand48Engine; }
76 }
77 
79 {
80  // this function is called at the beginning of an Event in LaserAlignment::upDate(const BeginOfEvent * myEvent)
81 
82  // use the random number generator service of the framework
84  unsigned int seed = rng->mySeed();
85 
86  // set the seed
87  theDRand48Engine->setSeed(seed);
88 
89  // number of LaserBeams
90  const G4int nLaserBeams = 8;
91 
92  // z position of the Laserdiodes (value from design drawings)
93  G4double LaserPositionZ = 1137.0 * mm;
94 
95  // Radius of the Laser ring
96  G4double LaserRingRadius = 564.0 * mm;
97 
98  // phi positions of the Laserdiodes (from CMS Note 2001/053 or from http://abbaneo.home.cern.ch/abbaneo/cms/layout)
99  G4double LaserPhi[nLaserBeams] = { G4double(7.0/112.0) * G4double(2.0 * M_PI),
100  G4double(23.0/112.0) * G4double(2.0 * M_PI),
101  G4double(33.0/112.0) * G4double(2.0 * M_PI),
102  G4double(49.0/112.0) * G4double(2.0 * M_PI),
103  G4double(65.0/112.0) * G4double(2.0 * M_PI),
104  G4double(77.0/112.0) * G4double(2.0 * M_PI),
105  G4double(93.0/112.0) * G4double(2.0 * M_PI),
106  G4double(103.0/112.0) * G4double(2.0 * M_PI) };
107 
108  // width of the LaserBeams
109  G4double LaserBeamSigmaX = 0.5 * mm;
110  G4double LaserBeamSigmaY = 0.5 * mm;
111 
112  // get the definition of the optical photon
113  G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
114  G4ParticleDefinition * theOpticalPhoton = theParticleTable->FindParticle("opticalphoton");
115 
116  // loop over the LaserBeams
117  for (int theBeam = 0; theBeam < nLaserBeams; theBeam++)
118  {
119  // code for forward and backward beam
120  // calculate x and y position of the current laser diode
121  G4double LaserPositionX = cos(LaserPhi[theBeam]) * LaserRingRadius;
122  G4double LaserPositionY = sin(LaserPhi[theBeam]) * LaserRingRadius;
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, LaserBeamSigmaX );
129  CLHEP::RandGaussQ aGaussObjY( *theDRand48Engine, LaserPositionY, LaserBeamSigmaY );
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  else if (theDirection == 1) // shoot in backward direction (-z)
154  {
155  theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, -1.0));
156  // set the polarization
157  setOptPhotonPolar(90.0);
158  // generate the particle
159  theParticleGun->GeneratePrimaryVertex(myEvent);
160  }
161  } // end looop over both beam directions
162  } // end looop over particles in beam
163  } // end loop over beams
164 }
165 
167 {
168  /* *********************************************************************** */
169  /* to get optical processes working properly, you have to make sure *
170  * that the photon polarisation is defined. */
171  /* *********************************************************************** */
172 
173  // first check if we have an optical photon
174  if ( theParticleGun->GetParticleDefinition()->GetParticleName() != "opticalphoton" )
175  {
176  edm::LogWarning("SimLaserAlignment:LaserBeamsBarrel") << "<LaserBeamsBarrel::setOptPhotonPolar()>: WARNING! The ParticleGun is not an optical photon";
177  return;
178  }
179 
180 // G4cout << " AC1CMS: The ParticleGun is an " << theParticleGun->GetParticleDefinition()->GetParticleName();
181  G4ThreeVector normal(1.0, 0.0, 0.0);
182  G4ThreeVector kphoton = theParticleGun->GetParticleMomentumDirection();
183  G4ThreeVector product = normal.cross(kphoton);
184  G4double modul2 = product * product;
185 
186  G4ThreeVector e_perpendicular(0.0, 0.0, 1.0);
187 
188  if ( modul2 > 0.0 ) { e_perpendicular = (1.0 / sqrt(modul2)) * product; }
189 
190  G4ThreeVector e_parallel = e_perpendicular.cross(kphoton);
191 
192  G4ThreeVector polar = cos(Angle) * e_parallel + sin(Angle) * e_perpendicular;
193 
194 // G4cout << ", the polarization = " << polar << G4endl;
195  theParticleGun->SetParticlePolarization(polar);
196 }
void GeneratePrimaries(G4Event *myEvent)
shoot optical photons into the detector at the beginning of an event
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
G4double thePhotonEnergy
LaserBeamsBarrel()
default constructor
virtual std::uint32_t mySeed() const =0
G4ParticleGun * theParticleGun
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
~LaserBeamsBarrel()
destructor
#define M_PI
void setOptPhotonPolar(G4double Angle)
set the polarisation of the photons
Definition: Angle.h:17
CLHEP::DRand48Engine * theDRand48Engine