CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
LaserBeamsBarrel.cc
Go to the documentation of this file.
1 
10 
14 
15 #include "CLHEP/Random/RandGaussQ.h"
16 #include "G4ParticleDefinition.hh"
17 #include "G4ParticleGun.hh"
18 #include "G4SystemOfUnits.hh"
19 #include "globals.hh" // Global Constants and typedefs
20 
21 LaserBeamsBarrel::LaserBeamsBarrel() : theParticleGun(nullptr), theDRand48Engine(nullptr) {
22  G4int nPhotonsGun = 1;
23  G4int nPhotonsBeam = 1;
24  G4double Energy = 1.15 * eV;
25  // call constructor with options
26  LaserBeamsBarrel(nPhotonsGun, nPhotonsBeam, Energy);
27 }
28 
29 LaserBeamsBarrel::LaserBeamsBarrel(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 LaserBeams
86  const G4int nLaserBeams = 8;
87 
88  // z position of the Laserdiodes (value from design drawings)
89  G4double LaserPositionZ = 1137.0 * mm;
90 
91  // Radius of the Laser ring
92  G4double LaserRingRadius = 564.0 * mm;
93 
94  // phi positions of the Laserdiodes (from CMS Note 2001/053 or from
95  // http://abbaneo.home.cern.ch/abbaneo/cms/layout)
96  G4double LaserPhi[nLaserBeams] = {G4double(7.0 / 112.0) * G4double(2.0 * M_PI),
97  G4double(23.0 / 112.0) * G4double(2.0 * M_PI),
98  G4double(33.0 / 112.0) * G4double(2.0 * M_PI),
99  G4double(49.0 / 112.0) * G4double(2.0 * M_PI),
100  G4double(65.0 / 112.0) * G4double(2.0 * M_PI),
101  G4double(77.0 / 112.0) * G4double(2.0 * M_PI),
102  G4double(93.0 / 112.0) * G4double(2.0 * M_PI),
103  G4double(103.0 / 112.0) * G4double(2.0 * M_PI)};
104 
105  // width of the LaserBeams
106  G4double LaserBeamSigmaX = 0.5 * mm;
107  G4double LaserBeamSigmaY = 0.5 * mm;
108 
109  // get the definition of the optical photon
110  G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
111  G4ParticleDefinition *theOpticalPhoton = theParticleTable->FindParticle("opticalphoton");
112 
113  // loop over the LaserBeams
114  for (int theBeam = 0; theBeam < nLaserBeams; theBeam++) {
115  // code for forward and backward beam
116  // calculate x and y position of the current laser diode
117  G4double LaserPositionX = cos(LaserPhi[theBeam]) * LaserRingRadius;
118  G4double LaserPositionY = sin(LaserPhi[theBeam]) * LaserRingRadius;
119 
120  // loop over all the particles in one beam
121  for (int theParticle = 0; theParticle < thenParticle; theParticle++) {
122  // get randomnumbers and calculate the position
123  CLHEP::RandGaussQ aGaussObjX(*theDRand48Engine, LaserPositionX, LaserBeamSigmaX);
124  CLHEP::RandGaussQ aGaussObjY(*theDRand48Engine, LaserPositionY, LaserBeamSigmaY);
125 
126  G4double theXPosition = aGaussObjX.fire();
127  G4double theYPosition = aGaussObjY.fire();
128  G4double theZPosition = LaserPositionZ;
129 
130  // set the properties of the newly created particle
131  theParticleGun->SetParticleDefinition(theOpticalPhoton);
132  theParticleGun->SetParticleTime(0.0 * ns);
133  theParticleGun->SetParticlePosition(G4ThreeVector(theXPosition, theYPosition, theZPosition));
134  theParticleGun->SetParticleEnergy(thePhotonEnergy);
135 
136  // loop over both directions of the beam
137  for (int theDirection = 0; theDirection < 2; theDirection++) {
138  // shoot in both beam directions ...
139  if (theDirection == 0) // shoot in forward direction (+z)
140  {
141  theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, 1.0));
142  // set the polarization
143  setOptPhotonPolar(90.0);
144  // generate the particle
145  theParticleGun->GeneratePrimaryVertex(myEvent);
146  } else if (theDirection == 1) // shoot in backward direction (-z)
147  {
148  theParticleGun->SetParticleMomentumDirection(G4ThreeVector(0.0, 0.0, -1.0));
149  // set the polarization
150  setOptPhotonPolar(90.0);
151  // generate the particle
152  theParticleGun->GeneratePrimaryVertex(myEvent);
153  }
154  } // end looop over both beam directions
155  } // end looop over particles in beam
156  } // end loop over beams
157 }
158 
160  /* *********************************************************************** */
161  /* to get optical processes working properly, you have to make sure *
162  * that the photon polarisation is defined. */
163  /* *********************************************************************** */
164 
165  // first check if we have an optical photon
166  if (theParticleGun->GetParticleDefinition()->GetParticleName() != "opticalphoton") {
167  edm::LogWarning("SimLaserAlignment:LaserBeamsBarrel")
168  << "<LaserBeamsBarrel::setOptPhotonPolar()>: WARNING! The ParticleGun "
169  "is not an optical photon";
170  return;
171  }
172 
173  // G4cout << " AC1CMS: The ParticleGun is an " <<
174  // theParticleGun->GetParticleDefinition()->GetParticleName();
175  G4ThreeVector normal(1.0, 0.0, 0.0);
176  G4ThreeVector kphoton = theParticleGun->GetParticleMomentumDirection();
177  G4ThreeVector product = normal.cross(kphoton);
178  G4double modul2 = product * product;
179 
180  G4ThreeVector e_perpendicular(0.0, 0.0, 1.0);
181 
182  if (modul2 > 0.0) {
183  e_perpendicular = (1.0 / sqrt(modul2)) * product;
184  }
185 
186  G4ThreeVector e_parallel = e_perpendicular.cross(kphoton);
187 
188  G4ThreeVector polar = cos(Angle) * e_parallel + sin(Angle) * e_perpendicular;
189 
190  // G4cout << ", the polarization = " << polar << G4endl;
191  theParticleGun->SetParticlePolarization(polar);
192 }
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:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
~LaserBeamsBarrel() override
destructor
#define M_PI
void setOptPhotonPolar(G4double Angle)
set the polarisation of the photons
CLHEP::DRand48Engine * theDRand48Engine
Log< level::Warning, false > LogWarning
void GeneratePrimaries(G4Event *myEvent) override
shoot optical photons into the detector at the beginning of an event