ATLAS Offline Software
Loading...
Searching...
No Matches
CosmicGenerator.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// -------------------------------------------------------------
6// File: CosmicGenerator/CosmicGenerator.cxx
7// Description:
8
9// The output will be stored in the transient event store so it can be
10// passed to the simulation.
11//
12// AuthorList:
13// W. Seligman: Initial Code 08-Nov-2002,
14// based on work by M. Shapiro and I. Hinchliffe
15//
16
17// Modification for increasing efficiency of muon hitting the detector:
18// H. Ma. March 17, 2006
19// Property: ExzCut:
20// if true, the method exzCut(...) will be called to apply a
21// energy dependent position cut on the surface.
22// This rejects low energy muons at large distance.
23// Property: RMax
24// Used by exzCut to reject non-projective muons, which are
25// too far out on the surface
26
27
28// Modifications to accomodate Pixel EndCap C Cosmic Test needs
29// Marian Zdrazil June 7, 2006 mzdrazil@lbl.gov
30//
31// Modifications to accomodate replacement of Pixel EndCap C by a Pixel EndCap A
32// Marian Zdrazil November 24, 2006 mzdrazil@lbl.gov
33//
34// Description:
35// ------------
36// It is easier and actually more useful to leave the EndCap A
37// in the vertical position (the way it is positioned in the ATLAS detector)
38// instead of rotating it clockwise by 90deg which corresponds to the
39// placement during the Pixel EndCap A cosmic test in SR1 in November 2006.
40// This is why we will generate cosmic muons coming from the positive Z-axis
41// direction better than rotating the whole setup in PixelGeoModel.
42
43// Modifications July 3rd 2007, Rob McPherson
44// - Fix mu+/mu- bug (always present in Athena versions)
45// - Fix sign of Py (since tag CosmicGenerator-00-00-21, muons only upward-going)
46
47// Optimize selection of events passed to Geant4 for full simulation:
48// - cut on energy based on pathlength in rock
49// - reweighting of generated cosmic rays
50// - geometrical cut in plane of pixel detector
51//
52// Juerg Beringer November 2007 JBeringer@lgl.gov
53// Robert Cahn November 2007 RNCahn@lbl.gov
54
55
59
60#include "CLHEP/Vector/ThreeVector.h"
61#include "CLHEP/Geometry/Normal3D.h"
62#include "CLHEP/Units/PhysicalConstants.h"
63#include "CLHEP/Random/RandFlat.h"
65
66#include <limits>
67#include <cmath>
68#include <vector>
69#include <string>
70#include <fstream>
71
72
73CLHEP::HepRandomEngine* CosmicGenerator::COSMIC_RANDOM_ENGINE{};
74extern "C" float cosmicrndm_(int* /*dummy*/)
75{
76 return CLHEP::RandFlat::shoot(CosmicGenerator::COSMIC_RANDOM_ENGINE);
77}
78
79//--------------------------------------------------------------------------
80CosmicGenerator::CosmicGenerator(const std::string& name,
81 ISvcLocator* pSvcLocator)
82 : GenModule(name,pSvcLocator)
83//--------------------------------------------------------------------------
84{
85}
86
87//---------------------------------------------------------------------------
89//---------------------------------------------------------------------------
90
91 // Initialize event count.
92
93 m_events = 0;
94
95 m_accepted=0;
96 m_rejected=0;
97
98 if(m_infile.value()=="NONE") {
101
103 gun->SetCosCut(m_ctcut);
105 float flux_withCT = gun->InitializeGenerator();
106
107 ATH_MSG_INFO( "Initialisation cosmic gun done." );
108 ATH_MSG_INFO( "Accepted diff flux after E and cos(theta) cuts = " << flux_withCT << " /cm^2/s" );
109 if (! m_doReweighting) {
110 // The following is only correct w/o reweighting
111 ATH_MSG_INFO( "Accepted total flux after E and cos(theta) cuts = " <<
112 flux_withCT*(m_xhig-m_xlow)/m_mm*(m_zhig-m_zlow)/m_mm << " /s" );
113 }
114
115 }
116 else {
117 ATH_MSG_INFO( "Cosmics are read from file " << m_infile );
118 m_ffile.open(m_infile.value().c_str());
119 if(!m_ffile) {
120 ATH_MSG_FATAL( "Could not open input file - stop! " );
121 return StatusCode::FAILURE;
122 }
123 m_readfile = true;
124 }
125
126 m_center=CLHEP::Hep3Vector(m_IPx, m_IPy, m_IPz);
127
128 return StatusCode::SUCCESS;
129
130}
131
132CLHEP::HepLorentzVector CosmicGenerator::generateVertex(void) {
133
134 // Get the pointer to the engine of the stream named SINGLE. If the
135 // stream does not exist is created automaticaly
136 CLHEP::HepRandomEngine* engine = COSMIC_RANDOM_ENGINE;
137
138 // Generate a random number according to the distribution.
139
140 float x_val = CLHEP::RandFlat::shoot(engine, m_xlow, m_xhig);
141 float z_val = CLHEP::RandFlat::shoot(engine, m_zlow, m_zhig);
142
143 // Generate a random number for time offset
144 float t_val = m_tmin; // Assign defined value
145 if(m_tmin < m_tmax){
146 t_val = CLHEP::RandFlat::shoot(engine, m_tmin, m_tmax);
147 }
148 else if(m_tmin.value() == m_tmax.value()){
149 t_val = m_tmin;
150 }
151 else ATH_MSG_FATAL("You specified m_tmin = " << m_tmin << " and m_tmax " << m_tmax);
152 CLHEP::HepLorentzVector p(x_val,m_yval,z_val, t_val*CLHEP::c_light);
153
154 return p;
155
156}
157
158CLHEP::HepLorentzVector CosmicGenerator::generateVertexReweighted(void) {
159
160 // Get the pointer to the engine of the stream named SINGLE. If the
161 // stream does not exist is created automaticaly
162 CLHEP::HepRandomEngine* engine = COSMIC_RANDOM_ENGINE;
163
164 // Generate non-uniform distribution of vertices to reflect azimuthal
165 // angle subtended by the sphere of radiusm m_radius
166 // Inside m_radius, the density of vertices is proportional to 2 pi r dr
167 // Outside m_radius, the density is proportional to 2r arcsin (m_radius/r)
168 // We approximate the latter by its maximum: m_radius * pi
169 // We generate vertices out to m_rvertmax.
170 // Integrating the approximated distribution gives
171 // pi r**2 for r < m_radius and pi m_radius r for r> m_radius
172 // So with ran in (0,1) we take r=max_r * ran for ran>m_radius/max_r
173 // and r= sqrt(m_radius*max_r*ran) for ran<m_radius/max_r
174 // for r>m_radius we use acceptance/rejection by comparing
175 // m_radius * pi * new_ran with 2r arcsin (m_radius/r)
176 int accept=0;
177 float max_r = m_rvertmax;
178 float r_val = 0.;
179 while(accept==0){
180 float ran_one = CLHEP::RandFlat::shoot(engine,0.,1.);
181 if(ran_one>(m_radius/max_r)){
182 r_val = ran_one*max_r;
183 float ran_two = CLHEP::RandFlat::shoot(engine,0.,1.);
184 if(m_radius*M_PI*ran_two<2*r_val*asin(m_radius/r_val)){
185 accept=1;
186 }
187 }
188 else
189 {
190 r_val = std::sqrt(m_radius*max_r*ran_one);
191 accept=1;
192 }
193 }
194 float ran_three= CLHEP::RandFlat::shoot(engine, 0.,2*M_PI);
195 float x_val = r_val*cos(ran_three);
196 float z_val = r_val*sin(ran_three);
197
198 // Generate a random number for time offset
199 float t_val = m_tmin; // Assign defined value
200 if(m_tmin < m_tmax){
201 t_val = CLHEP::RandFlat::shoot(engine, m_tmin, m_tmax);
202 }
203 else if(m_tmin.value() == m_tmax.value()){
204 t_val = m_tmin;
205 }
206 else ATH_MSG_FATAL( " You specified m_tmin = " << m_tmin << " and m_tmax " << m_tmax );
207 CLHEP::HepLorentzVector p(x_val,m_yval,z_val, t_val*CLHEP::c_light);
208
209 return p;
210}
211
212
213//---------------------------------------------------------------------------
215//---------------------------------------------------------------------------
216
217 ++m_events;
218 ATH_MSG_DEBUG( "Event #" << m_events);
219
220 assert(COSMIC_RANDOM_ENGINE != 0);
221 //Re-seed the random number stream
222 long seeds[7];
223 const EventContext& ctx = Gaudi::Hive::currentContext();
224 ATHRNG::calculateSeedsMC21(seeds, "COSMICS", ctx.eventID().event_number(), m_dsid, m_randomSeed);
225 COSMIC_RANDOM_ENGINE->setSeeds(seeds, 0); // NOT THREAD-SAFE
226 CLHEP::HepRandomEngine* engine = COSMIC_RANDOM_ENGINE;
227
228 // clear up the vectors
229 m_fourPos.clear();
230 m_fourMom.clear();
231 m_polarization.clear();
232 m_pdgCode.clear();
233
234
235 if(m_readfile)
236 {
237 if(!m_ffile.eof())
238 {
240 m_ffile >> evt;
241
242 ATH_MSG_VERBOSE( evt );
243
244 HepMC::Polarization thePolarization(0.0,0.0);
245 m_polarization.push_back(thePolarization);
246
247 //
248 // units are already converted to MeV's and mm.
249 //
250 m_fourPos.push_back(evt.Vertex());
251 m_fourMom.push_back(evt.Momentum());
252 m_pdgCode.push_back(evt.pdgID());
253
254 }
255 else
256 {
257 ATH_MSG_FATAL( "End of file reached - stop " );
258 exit(1);
259 return StatusCode::FAILURE;
260 }
261 }
262 else
263 {
264
265 bool accepted=false;
266 CLHEP::HepLorentzVector pp;
268 CLHEP::HepLorentzVector vert;
269 CLHEP::Hep3Vector vert3;
270 double theta1;
271 double phi1;
272 double mag1;
273
274 while(!accepted){
275
276 if (m_doReweighting && m_cavOpt) {
277 // The code here doesn't make sense without the sphere cut in the
278 // cavern optimization that is selected by OptimizeForCavern=True
280 vert3 = CLHEP::Hep3Vector(vert.x(),vert.y(),vert.z());
281
282 double vert_radius=sqrt(vert3.x()*vert3.x() + vert3.z()*vert3.z());
283
284 pp = gun->GenerateEvent();
285
286 theta1=pp.theta();
287 phi1=pp.phi();
288 mag1=pp.rho();
289
290 if (vert_radius>m_radius) {
291 phi1=atan2(vert.z(),vert.x())+M_PI;
292 float delta_phi=std::asin(m_radius/vert_radius);
293 phi1=phi1+CLHEP::RandFlat::shoot(engine, -delta_phi, delta_phi);
294 }
295 pp.setX(mag1*sin(theta1)*std::cos(phi1));
296 pp.setY(mag1*sin(theta1)*std::sin(phi1));
297
298 } else {
299 vert = generateVertex();
300 vert3 = CLHEP::Hep3Vector(vert.x(),vert.y(),vert.z());
301
302 pp = gun->GenerateEvent();
303
304 theta1=pp.theta();
305 phi1=pp.phi();
306 mag1=pp.rho();
307 }
308
309 CLHEP::Hep3Vector pp_corr(mag1*sin(theta1)*std::cos(phi1),
310 -mag1*std::cos(theta1),
311 mag1*std::sin(theta1)*std::sin(phi1));
312 CLHEP::Hep3Vector direction(pp_corr.x(),pp_corr.y(), pp_corr.z());
313
314 // if optimization activated, check for the direction of the generated muon
315 if(m_cavOpt) {
316
317 CLHEP::Hep3Vector center_dir=m_center-vert3;
318 double beta=direction.angle(center_dir);
319 double alpha=std::asin(m_radius/center_dir.r());
320
321 if(std::abs(beta)<alpha) {
322
323 if(m_exzCut) {
324 // Old optimization code - is it still useful?
325 CLHEP::HepLorentzVector pp2(pp_corr.x(),pp_corr.y(), pp_corr.z(), pp.e());
326 if( exzCut(vert3,pp2) ) {
327 accepted=true;
328 }
329 } else {
330
331 accepted = true;
332
333 ATH_MSG_DEBUG( "x0 = " << vert3.x()
334 << ", y0 = " << vert3.y()
335 << ", z0 = " << vert3.z()
336 << ", theta = " << pp.theta()
337 << ", phi = " << pp.phi()
338 << ", energy = " << pp.e()*m_GeV );
339
340 if (m_doPathlengthCut) {
341 double path = pathLengthInRock(vert3.x(),vert3.y(),vert3.z(),pp.theta(),pp.phi());
342 double energyLoss = 2.33e-3 * 244. * path; //FIXME Hardcoded values!
343 ATH_MSG_DEBUG( "Energy loss is " << energyLoss
344 << " --> " << (energyLoss>pp.e()*m_GeV ? "REJECTED" : "ACCEPTED") << " by pathlength cut");
345 if (energyLoss-m_energyCutThreshold > pp.e()*m_GeV) accepted = false;
346 }
347
349 bool aimedAtPixels = pointsAtPixels(vert3.x(),vert3.y(),vert3.z(),pp.theta(),pp.phi());
350 ATH_MSG_DEBUG( (aimedAtPixels ? "AIMED AT PIXELS" : "NOT AIMED AT PIXELS") );
351 if (!aimedAtPixels) accepted = false;
352 }
353
354 // FOR DEBUGGING ONLY
355 if (accepted) {
356 ATH_MSG_VERBOSE("The following event has been accepted for simulation:");
357 ATH_MSG_VERBOSE( "x0 = " << vert3.x() << ", y0 = " << vert3.y() << ", z0 = " << vert3.z()
358 << ", theta = " << pp.theta() << ", phi = " << pp.phi() << ", energy = " << pp.e()*m_GeV );
359
360 if (m_doPathlengthCut) {
361 double path = pathLengthInRock(vert3.x(),vert3.y(),vert3.z(),pp.theta(),pp.phi());
362 double energyLoss = 2.33e-3 * 244. * path;
363 ATH_MSG_VERBOSE( "Energy loss is " << energyLoss
364 << " --> " << (energyLoss>pp.e()*m_GeV ? "REJECTED" : "ACCEPTED") << " by pathlength cut" );
365 }
366
368 bool aimedAtPixels = pointsAtPixels(vert3.x(),vert3.y(),vert3.z(),pp.theta(),pp.phi());
369 ATH_MSG_VERBOSE( (aimedAtPixels ? "AIMED AT PIXELS" : "NOT AIMED AT PIXELS") );
370 }
371
372 }
373
374 }
375 }
376
377 if(accepted) {
378 m_accepted++;
379 } else {
380 ATH_MSG_VERBOSE("Rejected muon due to cavern optimization request!");
381 m_rejected++;
382 }
383 }
384 else if(m_srOneOpt == 1) {
385 CLHEP::Hep3Vector srOneVec(direction.x(), direction.y(), direction.z());
386 if(mag1 < 0) // Check if momentum vector is flipped.
387 srOneVec *= -1;
388
389 if( (srOneVec.phi() >= -2.25) && (srOneVec.phi() <= -1.7) &&
390 (srOneVec.theta() >= 0.85) && (srOneVec.theta() <= 2.25) ) { //FIXME Hardcoded values!
391 accepted = true;
392 m_accepted++;
393 ATH_MSG_DEBUG("Muon accepted by SR1 SCT/TRT optimization!");
394 } else {
395 ATH_MSG_DEBUG("Rejected muon due to SR1 SCT/TRT optimization request!");
396 m_rejected++;
397 }
398 }
399 else if(m_srOneOpt == 2) {
400 CLHEP::Hep3Vector srOneVec(direction.x(), direction.y(), direction.z());
401 if(mag1 < 0) // Check if momentum vector is flipped.
402 srOneVec *= -1;
403
404 if( (srOneVec.phi() >= -1.68) && (srOneVec.phi() <= -1.08) &&
405 (srOneVec.theta() >= 0.29) && (srOneVec.theta() <= 0.72) ) { //FIXME Hardcoded values!
406 accepted = true;
407 m_accepted++;
408 ATH_MSG_DEBUG("Muon accepted by SR1 SCT/TRT EndCapC optimization!");
409 } else {
410 ATH_MSG_DEBUG("Rejected muon due to SR1 SCT/TRT EndcapC optimization request!");
411 m_rejected++;
412 }
413 }
414
415 else if(m_srOnePixECOpt) {
416 CLHEP::Hep3Vector srOneVec(direction.x(), direction.y(), direction.z());
417 if(mag1 < 0) // Check if momentum vector is flipped.
418 srOneVec *= -1;
419
420 if( (srOneVec.phi() >= m_phimin) && (srOneVec.phi() <= m_phimax) &&
421 (srOneVec.theta() >= m_thetamin) && (srOneVec.theta() <= m_thetamax) ) {
422 accepted = true;
423 m_accepted++;
424 ATH_MSG_DEBUG("Muon accepted by SR1 Pixel EndCap optimization!");
425 } else {
426 ATH_MSG_DEBUG("Rejected muon due to SR1 Pixel EndCap optimization request!");
427 m_rejected++;
428 }
429 }
430
431 else if (m_muonECOpt) {
432 double coor_x, coor_y, coor_z;
433 coor_z = m_zpos; // defined in jobOpt.
434 coor_x = direction.x()*(coor_z - vert.z())/direction.z() +vert.x();
435 coor_y = direction.y()*(coor_z - vert.z())/direction.z() +vert.y();
436 if( ((coor_x)*(coor_x) + (coor_y)*(coor_y)) <= m_radius*m_radius ) {
437 accepted = true;
438 m_accepted++;
439 } else {
440 coor_z = -m_zpos;
441 coor_x = direction.x()*(coor_z - vert.z())/direction.z() +vert.x();
442 coor_y = direction.y()*(coor_z - vert.z())/direction.z() +vert.y();
443 if( ((coor_x)*(coor_x) + (coor_y)*(coor_y)) <= m_radius*m_radius ) {
444 accepted = true;
445 m_accepted++;
446 } else {
447 ATH_MSG_DEBUG("Rejected muon due to Muon EndCap optimization request!");
448 m_rejected++;
449 }
450 }
451 }
452
453 else accepted=true; // if no opt required accept the first muon
454 }
455
456 pp.setX(pp.x()*m_GeV);
457 pp.setY(pp.y()*m_GeV);
458 pp.setZ(pp.z()*m_GeV);
459 pp.setT(pp.t()*m_GeV);
460
461 // Get the mass of the particle to be generated
462 int charge = gun->GetMuonCharge();
463 // m_pdgCode.push_back(charge*13);
464 m_pdgCode.push_back(charge*-13);
465
466 const HepPDT::ParticleData* particle = particleData(std::abs(m_pdgCode.back()));
467 if (particle==nullptr){
468 ATH_MSG_FATAL( "Particle with PDG ID=" << std::abs(m_pdgCode.back()) << " returned a nullptr" );
469 return StatusCode::FAILURE;
470 }
471
472 double mass = particle->mass().value();
473
474 // Compute the kinematic values. First, the vertex 4-vector:
475 double x = vert.x();
476 double y = vert.y();
477 double z = vert.z();
478 double t = vert.t();
479
480 // Do we need to swap Y- and Z-axis for the PixelEndCap A Cosmic test ?
481 // if not...do nothing...if so, invert position of y- and z- coordinate
482 //
483 // but not only that...change also the direction of the incoming cosmic muon(s),
484 // they must go towards the pixel endcap A, i.e. y -> -y
485 //
486 if(!m_swapYZAxis)
487 m_fourPos.push_back(CLHEP::HepLorentzVector(x,y,z,t));
488 else
489 m_fourPos.push_back(CLHEP::HepLorentzVector(x,z,y,t));
490
491 // Set the polarization. Realistically, this is going to be zero
492 // for most studies.
493 HepMC::Polarization thePolarization(0.0);
494 m_polarization.push_back(thePolarization);
495
496
497 // The method of calculating e, theta, and phi depends on the user's
498 // commands. Let the KinematicManager handle it.
499 double e = pp.e();
500 double theta = pp.theta();
501 double phi = pp.phi();
502
503 // At this point, we have e, theta, and phi. Put them together to
504 // get the four-momentum.
505
506 double p2 = e*e - mass*mass;
507 if ( p2 < 0 )
508 {
509 ATH_MSG_ERROR( "Event #" << m_events
510 << " E=" << e << ", mass=" << mass
511 << " -- you have generated a tachyon! Increase energy or change particle ID." );
512 return StatusCode::FAILURE;
513 }
514
515 double p = std::sqrt(p2);
516 double px = p*std::sin(theta)*std::cos(phi);
517 double pz = p*std::sin(theta)*std::sin(phi);
518 double py = -p*std::cos(theta);
519
520 // Do we need to swap Y- and Z-axis for the PixelEndCap C Cosmic test ?
521 // if not...do nothing...if so, invert position of y- and z- coordinate
522 //
523 // well and don't forget about the direction of the incoming cosmic muon(s) either
524 // that means: y -> -y
525 //
526 if(!m_swapYZAxis) {
527 // Line below corrupted py sign and forces muons to be upwards, not downwards.
528 // m_fourMom.push_back(CLHEP::HepLorentzVector(px,-py,pz,pp.e()));
529 m_fourMom.push_back(CLHEP::HepLorentzVector(px,py,pz,pp.e()));
530 }
531 else
532 m_fourMom.push_back(CLHEP::HepLorentzVector(px,pz,-py,pp.e()));
533
535 " (x,y,z,t) = ("
536 << m_fourPos.back().x() << ","
537 << m_fourPos.back().y() << ","
538 << m_fourPos.back().z() << ","
539 << m_fourPos.back().t() << "), (Px,Py,Pz,E) = ("
540 << m_fourMom.back().px() << ","
541 << m_fourMom.back().py() << ","
542 << m_fourMom.back().pz() << ","
543 << m_fourMom.back().e() << ")" );
545 " (theta,phi) = (" << theta << "," << phi << "), "
546 << "polarization(theta,phi) = ("<< m_polarization.back().theta() << ","<< m_polarization.back().phi() << ")" );
547 }
548 return StatusCode::SUCCESS;
549
550}
551
552
553//---------------------------------------------------------------------------
555//---------------------------------------------------------------------------
556 // Get the KinematicManager.
557
558 if(m_cavOpt){
559
560 ATH_MSG_INFO("********************************************");
561 ATH_MSG_INFO("** you have run CosmicGenerator with some ");
562 ATH_MSG_INFO("** optimizations for cavern simulation");
563 ATH_MSG_INFO("** "<<m_accepted<<" muons were accepted");
564 ATH_MSG_INFO("** "<<m_rejected<<" muons were rejected");
565 ATH_MSG_INFO("********************************************");
566
567 }
568 if(m_srOneOpt == 1){
569
570 ATH_MSG_INFO("**********************************************");
571 ATH_MSG_INFO("** you have run CosmicGenerator with some ");
572 ATH_MSG_INFO("** optimizations for SR1 SCT/TRT simulation");
573 ATH_MSG_INFO("** "<<m_accepted<<" muons were accepted");
574 ATH_MSG_INFO("** "<<m_rejected<<" muons were rejected");
575 ATH_MSG_INFO("**********************************************");
576
577 }
578 if(m_srOneOpt == 2){
579
580 ATH_MSG_INFO("**********************************************");
581 ATH_MSG_INFO("** you have run CosmicGenerator with some ");
582 ATH_MSG_INFO("** optimizations for SR1 SCT/TRT EndcapC simulation");
583 ATH_MSG_INFO("** "<<m_accepted<<" muons were accepted");
584 ATH_MSG_INFO("** "<<m_rejected<<" muons were rejected");
585 ATH_MSG_INFO("**********************************************");
586
587 }
588
589 if(m_srOnePixECOpt){
590
591 ATH_MSG_INFO("***************************************************");
592 ATH_MSG_INFO("** you have run CosmicGenerator with some ");
593 ATH_MSG_INFO("** optimizations for SR1 Pixel EndCap simulation");
594 ATH_MSG_INFO("** "<<m_accepted<<" muons were accepted");
595 ATH_MSG_INFO("** "<<m_rejected<<" muons were rejected");
596 ATH_MSG_INFO("***************************************************");
597
598 if(m_swapYZAxis){
599 ATH_MSG_INFO("***************************************************");
600 ATH_MSG_INFO(" You have swapped Y- and Z-axis, i.e. muons are ");
601 ATH_MSG_INFO(" not coming from the top any more !!! ");
602 ATH_MSG_INFO("***************************************************");
603 }
604
605 }
606
607 if(m_muonECOpt) {
608
609 ATH_MSG_INFO("***************************************************");
610 ATH_MSG_INFO("** you have run CosmicGenerator with some " );
611 ATH_MSG_INFO("** filters for cosmic muon simulation" );
612 ATH_MSG_INFO("** "<<m_accepted<<" muons were accepted" );
613 ATH_MSG_INFO("** "<<m_rejected<<" muons were rejected" );
614 ATH_MSG_INFO("***************************************************");
615 }
616
617 return StatusCode::SUCCESS;
618}
619
620
621//---------------------------------------------------------------------------
622StatusCode CosmicGenerator::fillEvt(HepMC::GenEvent* event) {
623//---------------------------------------------------------------------------
624
625
626 // loop over generated vertices
627 if(m_fourMom.size()==m_fourPos.size()&&m_fourMom.size()==m_polarization.size()){
628
629 for(std::size_t v=0;v<m_fourMom.size();++v){
630
631 // Note: The vertex and particle are owned by the event, so the
632 // event is responsible for those pointers.
633
634 // Create the particle, and specify its polarization.
635
636 HepMC::GenParticlePtr particle = HepMC::newGenParticlePtr( HepMC::FourVector(m_fourMom[v].px(),m_fourMom[v].py(),m_fourMom[v].pz(),m_fourMom[v].e()), m_pdgCode[v], 1);
638
639 // Create the vertex, and add the particle to the vertex.
640 HepMC::GenVertexPtr vertex = HepMC::newGenVertexPtr(HepMC::FourVector(m_fourPos[v].x(),m_fourPos[v].y(),m_fourPos[v].z(),m_fourPos[v].t()));
641 vertex->add_particle_out( std::move(particle) );
642
643 // Add the vertex to the event.
644 event->add_vertex( std::move(vertex) );
645
646 }
647
648 event->set_event_number(m_events); // Set the event number
649 if (event->weights().empty()){
650 event->weights().push_back(1.0);
651 }
652 return StatusCode::SUCCESS;
653 } else {
654 ATH_MSG_ERROR("Wrong different number of vertexes/momenta/polaritazions!");
655 return StatusCode::FAILURE;
656 }
657
658}
659
660
661// Energy dependent position cut on the surface.
662bool CosmicGenerator::exzCut(const CLHEP::Hep3Vector& pos,const CLHEP::HepLorentzVector& p)
663{
664// p is in GeV...
665
666 double r =0;
667 bool cut = false;
668 if(pos.z()<0){
669 r = std::sqrt((std::pow(pos.x(),2)+std::pow(pos.z()+28000,2))) ; //FIXME Hardcoded values!
670 double e = 0.45238*r+5000 ; //FIXME Hardcoded values!
671 cut = p.e()*m_GeV>e;
672 }
673 else
674 {
675 r = std::sqrt((std::pow(pos.x(),2)+std::pow(pos.z()-20000,2))) ; //FIXME Hardcoded values!
676 if(r<15000) { //FIXME Hardcoded values!
677 cut = true;
678 } else
679 {
680 double e = 0.461538*(r-15000)+10000 ; //FIXME Hardcoded values!
681 cut = p.e()*m_GeV>e;
682 ATH_MSG_VERBOSE("z>0 r , e, p.e = "<<r <<" " <<e <<" " <<p.e()*m_GeV);
683 }
684 }
685
686 cut = cut && r < m_rmax ;
687
688 return cut;
689}
690
691
692// Estimate pathlength in rock towards the pixel detector, taking into account both
693// the PX14 and PX16 shafts. The shaft positions are currently hard-coded.
694double CosmicGenerator::pathLengthInRock(double xgen, double ygen, double zgen, double theta, double phi) {
695 // y is vertical direction, z along beam, major shaft has z>0
696
697 // Definition of shafts and cavern
698 double p14_x = 1700; //FIXME Hardcoded values!
699 double p14_z = 13500; //FIXME Hardcoded values!
700 double p14_radius = 9000.; //FIXME Hardcoded values!
701 double p16_x = 1700; //FIXME Hardcoded values!
702 double p16_z = -20000; //FIXME Hardcoded values!
703 double p16_radius = 6300.; //FIXME Hardcoded values!
704 double y1 = 26400; // ! mm cavern height above IP //FIXME Hardcoded values!
705
706 // direction of trajectory
707 // x=x0 - t sinth cosphi; y=y0 + t costh; z=z0 - t sinth sinphi
708 double cosphi = std::cos(phi);
709 double sinphi = std::sin(phi);
710 double costh = std::cos(theta);
711 double sinth = std::sin(theta);
712
713 double y0 = m_ysurface;
714 double t = (ygen-y0)/costh;
715 double x0 = xgen + t*sinth*cosphi; // x position at y=0
716 double z0 = zgen + t*sinth*sinphi; // z position at y=0
717
718 // full path length ignoring shaft
719 double full_distance = (y0-y1)/costh;
720
721 // does trajectory intersect p14 cylinder?
722 double z_mid14 = (x0-p14_x)*sinphi-(z0-p14_z)*cosphi;
723 double min_dist14 = std::abs(z_mid14); //minimum distance of line from center
724 double shaft_distance14 = 0.;
725 if (min_dist14<p14_radius) {
726
727 // z values at intersections
728 double z_plus14 = -cosphi*z_mid14+sinphi*std::sqrt(std::pow(p14_radius,2.)-std::pow(z_mid14,2.)) + p14_z;
729 double z_minus14 = -cosphi*z_mid14-sinphi*std::sqrt(std::pow(p14_radius,2.)-std::pow(z_mid14,2.)) + p14_z;
730
731 // y values at intersections
732 double y_plus14 = y0-costh*(z_plus14-z0)/sinth/sinphi;
733 double y_minus14 = y0-costh*(z_minus14-z0)/sinth/sinphi;
734 double y_great14 = y_plus14>y_minus14 ? y_plus14 : y_minus14;
735 double y_less14 = y_plus14>y_minus14 ? y_minus14 : y_plus14;
736
737 // top intersection must occur above bottom of shaft
738 if ( (y_great14>y1) && (y_less14<y0) ) {
739 double y_top14 = y_great14<y0 ? y_great14 : y0;
740 double y_bottom14 = y_less14>y1 ? y_less14 : y1;
741 shaft_distance14 = (y_top14-y_bottom14)/costh;
742 }
743 }
744
745 // does trajectory intersect p16 cylinder?
746 double z_mid16 = (x0-p16_x)*sinphi-(z0-p16_z)*cosphi;
747 double min_dist16 = std::abs(z_mid16);
748 double shaft_distance16 = 0.;
749 if (min_dist16<p16_radius) {
750
751 // z values at intersections
752 double z_plus16 = -cosphi*z_mid16+sinphi*std::sqrt(std::pow(p16_radius,2.)-std::pow(z_mid16,2.)) + p16_z;
753 double z_minus16 = -cosphi*z_mid16-sinphi*std::sqrt(std::pow(p16_radius,2.)-std::pow(z_mid16,2.)) + p16_z;
754
755 // determine y values at intersections
756 double y_plus16 = y0-costh*(z_plus16-z0)/sinth/sinphi;
757 double y_minus16 = y0-costh*(z_minus16-z0)/sinth/sinphi;
758 double y_great16 = y_plus16>y_minus16 ? y_plus16 : y_minus16;
759 double y_less16 = y_plus16>y_minus16 ? y_minus16 : y_plus16;
760
761 // top intersection must occur above bottom of shaft
762 if ( (y_great16>y1) && (y_less16<y0) ) {
763 double y_top16 = y_great16<y0 ? y_great16 : y0;
764 double y_bottom16 = y_less16>y1 ? y_less16 : y1;
765 shaft_distance16 = (y_top16-y_bottom16)/costh;
766 }
767 }
768
769 double rock_distance = full_distance - shaft_distance14-shaft_distance16;
770 return rock_distance;
771}
772
773
774// Check if trajectory points towards a horizontal rectangle centered at the pixel detector
775bool CosmicGenerator::pointsAtPixels(double xgen, double ygen, double zgen, double theta, double phi) {
776 // y is vertical direction, z along beam, major shaft has z<0
777 bool does = false;
778
779 // direction of trajectory
780 // x=xgen+ t sinth cosphi; y=ygen+t costh; z=zgen+t sinth sinphi
781 double cosphi = std::cos(phi);
782 double sinphi = std::sin(phi);
783 double costh = std::cos(theta);
784 double sinth = std::sin(theta);
785 double t = ygen/costh; //for parameterized trajectory
786 double x_pos = xgen + t*sinth*cosphi; //x position at y=0
787 double z_pos = zgen + t*sinth*sinphi; //z position at y=0
788
789 ATH_MSG_VERBOSE("x_pos = " << x_pos << ", z_pos = " << z_pos);
790
791 if((std::abs(x_pos)<m_pixelplanemaxx)&&(std::abs(z_pos)<m_pixelplanemaxz)){
792 does=true;
793 }
794
795 return does;
796}
#define M_PI
Scalar phi() const
phi method
Scalar theta() const
theta method
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
double charge(const T &p)
Definition AtlasPID.h:997
float cosmicrndm_(int *)
#define y
#define x
#define z
BooleanProperty m_doAimedAtPixelsCut
FloatProperty m_thetamax
FloatProperty m_zhig
std::vector< CLHEP::HepLorentzVector > m_fourPos
FloatProperty m_IPz
FloatProperty m_xhig
FloatProperty m_ctcut
BooleanProperty m_doPathlengthCut
BooleanProperty m_cavOpt
IntegerProperty m_printMod
BooleanProperty m_doReweighting
FloatProperty m_phimin
CLHEP::HepLorentzVector generateVertexReweighted(void)
FloatProperty m_phimax
virtual StatusCode genInitialize()
For initializing the generator, if required.
FloatProperty m_tmin
FloatProperty m_IPy
DoubleProperty m_ysurface
bool exzCut(const CLHEP::Hep3Vector &pos, const CLHEP::HepLorentzVector &p)
DoubleProperty m_pixelplanemaxx
FloatProperty m_emin
std::ifstream m_ffile
static constexpr float m_mm
FloatProperty m_zpos
IntegerProperty m_srOneOpt
bool pointsAtPixels(double xgen, double ygen, double zgen, double theta, double phi)
IntegerProperty m_printEvent
FloatProperty m_tmax
virtual StatusCode fillEvt(HepMC::GenEvent *evt)
For filling the HepMC event object.
StringProperty m_infile
BooleanProperty m_srOnePixECOpt
DoubleProperty m_energyCutThreshold
static CLHEP::HepRandomEngine * COSMIC_RANDOM_ENGINE
Static pointer to random number generator for use by.
CLHEP::HepLorentzVector generateVertex(void)
BooleanProperty m_exzCut
FloatProperty m_zlow
FloatProperty m_IPx
std::vector< HepMC::Polarization > m_polarization
DoubleProperty m_rvertmax
std::vector< CLHEP::HepLorentzVector > m_fourMom
CLHEP::Hep3Vector m_center
DoubleProperty m_pixelplanemaxz
FloatProperty m_emax
FloatProperty m_thetamin
virtual StatusCode callGenerator()
For calling the generator on each iteration of the event loop.
FloatProperty m_yval
CosmicGenerator(const std::string &name, ISvcLocator *pSvcLocator)
FloatProperty m_radius
BooleanProperty m_muonECOpt
static constexpr float m_GeV
BooleanProperty m_swapYZAxis
std::vector< int > m_pdgCode
FloatProperty m_rmax
virtual StatusCode genFinalize()
For finalising the generator, if required.
FloatProperty m_xlow
double pathLengthInRock(double xgen, double ygen, double zgen, double theta, double phi)
IntegerProperty m_dsid
void PrintLevel(int printevt, int printmod)
Definition CosmicGun.cxx:94
void SetEnergyRange(float emin, float emax)
int GetMuonCharge(void)
float InitializeGenerator()
Definition CosmicGun.cxx:85
CLHEP::HepLorentzVector GenerateEvent(void)
void SetCosCut(float ctcut)
static CosmicGun * GetCosmicGun(void)
Definition CosmicGun.cxx:53
const HepPDT::ParticleData * particleData(int pid) const
Access an element in the particle data table.
Definition GenBase.h:126
GenModule(const std::string &name, ISvcLocator *pSvcLocator)
Constructor.
Definition GenModule.cxx:14
CLHEP::HepRandomEngine * getRandomEngineDuringInitialize(const std::string &streamName, unsigned long int randomSeedOffset, unsigned int conditionsRun=1, unsigned int lbn=1) const
Definition GenModule.cxx:53
IntegerProperty m_randomSeed
Seed for random number engine.
Definition GenModule.h:84
int r
Definition globals.cxx:22
void calculateSeedsMC21(long *seeds, const std::string &algName, uint64_t ev, uint64_t run, uint64_t offset=0)
Set the random seed using a string (e.g.
void set_polarization(T &a, const Polarization &b)
HepMC::GenVertex * GenVertexPtr
Definition GenVertex.h:59
GenVertexPtr newGenVertexPtr(const HepMC::FourVector &pos=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), const int i=0)
Definition GenVertex.h:64
GenParticlePtr newGenParticlePtr(const HepMC::FourVector &mom=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), int pid=0, int status=0)
Definition GenParticle.h:39
GenParticle * GenParticlePtr
Definition GenParticle.h:37