CosmicsGen_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 /// \file CosmicsGen_module.cc
3 /// \brief Generator for cosmic-rays
4 ///
5 /// Module to produce cosmic ray MC using CRY
6 ///
7 /// \author brebel@fnal.gov
8 ////////////////////////////////////////////////////////////////////////
9 
10 // ROOT includes
11 #include "TH1.h"
12 #include "TH2.h"
13 
14 // Framework includes
18 #include "fhiclcpp/ParameterSet.h"
20 #include "art_root_io/TFileService.h"
21 
22 // art extensions
23 #include "nurandom/RandomUtils/NuRandomService.h"
24 
25 // larsoft includes
28 #include "nugen/EventGeneratorBase/evgenbase.h"
32 
33 namespace evgen {
34 
35  /// A module to check the results from the Monte Carlo generator
36  class CosmicsGen : public art::EDProducer {
37  public:
38  explicit CosmicsGen(fhicl::ParameterSet const& pset);
39 
40  private:
41 
42  void produce(art::Event& evt) override;
43  void beginJob() override;
44  void beginRun(art::Run& run) override;
45 
46  std::vector<double> fbuffbox;
47 
48  TH2F* fPhotonAngles; ///< Photon rate vs angle
49  TH2F* fPhotonAnglesLo; ///< Photon rate vs angle, low momenta
50  TH2F* fPhotonAnglesMi; ///< Photon rate vs angle, middle momenta
51  TH2F* fPhotonAnglesHi; ///< Photon rate vs angle, high momenta
52  TH1F* fPhotonCosQ; ///< Photon rate vs cos(Q)
53  TH1F* fPhotonEnergy; ///< Photon energy (GeV)
54  TH1F* fPhotonsPerSample; ///< number of photons in the sampled time window
55  TH1F* fPhotonsInCStat; ///< number of photons in the cryostat during
56  ///< the sampled time window
57  TH1F* fPhotonsInTPC; ///< number of photons in the tpc during
58  ///< the sampled time window
59 
60  TH2F* fElectronAngles; ///< Electron rate vs angle
61  TH2F* fElectronAnglesLo; ///< Electron rate vs angle, low momenta
62  TH2F* fElectronAnglesMi; ///< Electron rate vs angle, middle momenta
63  TH2F* fElectronAnglesHi; ///< Electron rate vs angle, high momenta
64  TH1F* fElectronCosQ; ///< Electron rate vs cos(Q)
65  TH1F* fElectronEnergy; ///< Electron energy (GeV)
66  TH1F* fElectronsPerSample; ///< number of electrons in the sampled time window
67  TH1F* fElectronsInCStat; ///< number of electrons in the cryostat during
68  ///< the sampled time window
69  TH1F* fElectronsInTPC; ///< number of electrons in the tpc during
70  ///< the sampled time window
71 
72  TH2F* fMuonAngles; ///< Muon rate vs angle
73  TH2F* fMuonAnglesLo; ///< Muon rate vs angle, low momenta
74  TH2F* fMuonAnglesMi; ///< Muon rate vs angle, middle momenta
75  TH2F* fMuonAnglesHi; ///< Muon rate vs angle, high momenta
76  TH1F* fMuonCosQ; ///< Muon rate vs cos(Q)
77  TH1F* fMuonEnergy; ///< Muon energy (GeV)
78  TH1F* fMuonsPerSample; ///< number of muons in the sampled time window
79  TH1F* fMuonsInCStat; ///< number of muons in the cryostat during
80  ///< the sampled time window
81  TH1F* fMuonsInTPC; ///< number of muons in the tpc during
82  ///< the sampled time window
83  CLHEP::HepRandomEngine& fEngine; ///< art-managed random-number engine
84  evgb::CRYHelper fCRYHelp; ///< CRY generator object
85  };
86 }
87 
88 namespace evgen{
89 
90  //____________________________________________________________________________
92  : art::EDProducer{pset}
93  , fbuffbox{pset.get<std::vector<double>>("BufferBox",{0.0, 0.0, 0.0, 0.0, 0.0, 0.0})}
94  // create a default random engine; obtain the random seed from NuRandomService,
95  // unless overridden in configuration with key "Seed"
96  , fEngine(art::ServiceHandle<rndm::NuRandomService>()->createEngine(*this, pset, "Seed"))
97  , fCRYHelp{pset,
98  fEngine,
100  {
101  produces< std::vector<simb::MCTruth> >();
102  produces< sumdata::RunData, art::InRun >();
103  }
104 
105  //____________________________________________________________________________
107  {
109 
110  fPhotonAngles = tfs->make<TH2F>("fPhotonAngles", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
111  fPhotonAnglesLo = tfs->make<TH2F>("fPhotonAnglesLo", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
112  fPhotonAnglesMi = tfs->make<TH2F>("fPhotonAnglesMi", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
113  fPhotonAnglesHi = tfs->make<TH2F>("fPhotonAnglesHi", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
114  fPhotonCosQ = tfs->make<TH1F>("fPhotonCosQ", ";cos#theta;tracks", 50,-1.0,1.0);
115  fPhotonEnergy = tfs->make<TH1F>("fPhotonEnergy", ";E (GeV)", 5000,0.0,1000.0);
116  fPhotonsPerSample = tfs->make<TH1F>("fPhotonsPerSample", ";Number Photons;Samples", 100, 0, 1000);
117  fPhotonsInCStat = tfs->make<TH1F>("fPhotonsInCryostat", ";Number Photons;Samples", 100, 0, 1000);
118  fPhotonsInTPC = tfs->make<TH1F>("fPhotonsInTPC", ";Number Photons;Samples", 100, 0, 1000);
119 
120  fElectronAngles = tfs->make<TH2F>("fElectronAngles", ";#phi;cos#theta",36,-180.0,180.0,50,-1.0,1.0);
121  fElectronAnglesLo = tfs->make<TH2F>("fElectronAnglesLo", ";#phi;cos#theta",36,-180.0,180.0,50,-1.0,1.0);
122  fElectronAnglesMi = tfs->make<TH2F>("fElectronAnglesMi", ";#phi;cos#theta",36,-180.0,180.0,50,-1.0,1.0);
123  fElectronAnglesHi = tfs->make<TH2F>("fElectronAnglesHi", ";#phi;cos#theta",36,-180.0,180.0,50,-1.0,1.0);
124  fElectronCosQ = tfs->make<TH1F>("fElectronCosQ", ";cos#theta;tracks",50,-1.0,1.0);
125  fElectronEnergy = tfs->make<TH1F>("fElectronEnergy", ";E (GeV)", 5000,0.0,1000.0);
126  fElectronsPerSample = tfs->make<TH1F>("fElectronsPerSample", ";Number Electrons;Samples", 100, 0, 1000);
127  fElectronsInCStat = tfs->make<TH1F>("fElectronsInCryotat", ";Number Electrons;Samples", 100, 0, 1000);
128  fElectronsInTPC = tfs->make<TH1F>("fElectronsInTPC", ";Number Electrons;Samples", 100, 0, 1000);
129 
130  fMuonAngles = tfs->make<TH2F>("fMuonAngles", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
131  fMuonAnglesLo = tfs->make<TH2F>("fMuonAnglesLo", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
132  fMuonAnglesMi = tfs->make<TH2F>("fMuonAnglesMi", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
133  fMuonAnglesHi = tfs->make<TH2F>("fMuonAnglesHi", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
134  fMuonCosQ = tfs->make<TH1F>("fMuonCosQ", ";cos#theta;tracks",50,-1.0,1.0);
135  fMuonEnergy = tfs->make<TH1F>("fMuonEnergy", ";E (GeV)", 5000,0.0,1000.0);
136  fMuonsPerSample = tfs->make<TH1F>("fMuonsPerSample", ";Number Muons;Samples", 100, 0, 1000);
137  fMuonsInCStat = tfs->make<TH1F>("fMuonsInCryostat", ";Number Muons;Samples", 100, 0, 1000);
138  fMuonsInTPC = tfs->make<TH1F>("fMuonsInTPC", ";Number Muons;Samples", 100, 0, 1000);
139  }
140 
141  //____________________________________________________________________________
143  {
145  run.put(std::make_unique<sumdata::RunData>(geo->DetectorName()));
146  }
147 
148  //____________________________________________________________________________
150  {
151  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
152 
153  // fill some histograms about this event
155 
156  int nCrossCryostat = 0;
157 
158  simb::MCTruth truth;
159 
160  while(nCrossCryostat < 1){
161 
162  simb::MCTruth pretruth;
164  fCRYHelp.Sample(pretruth,
165  geom->SurfaceY(),
166  geom->DetLength(),
167  0);
168 
169  int numPhotons = 0;
170  int numElectrons = 0;
171  int numMuons = 0;
172  int allPhotons = 0;
173  int allElectrons = 0;
174  int allMuons = 0;
175 
176  // loop over particles in the truth object
177  for(int i = 0; i < pretruth.NParticles(); ++i){
178  simb::MCParticle particle = pretruth.GetParticle(i);
179  const TLorentzVector& v4 = particle.Position();
180  const TLorentzVector& p4 = particle.Momentum();
181  double x0[3] = {v4.X(), v4.Y(), v4.Z() };
182  double dx[3] = {p4.Px(), p4.Py(), p4.Pz()};
183 
184  if (std::abs(particle.PdgCode())==13) ++allMuons;
185  else if (std::abs(particle.PdgCode())==22) ++allPhotons;
186  else if (std::abs(particle.PdgCode())==11) ++allElectrons;
187 
188  TH1F* hCosQ = 0;
189  TH2F* hAngles = 0;
190  TH2F* hAnglesLo = 0;
191  TH2F* hAnglesMi = 0;
192  TH2F* hAnglesHi = 0;
193  TH1F* hEnergy = 0;
194  if (std::abs(particle.PdgCode())==13) {
195  hCosQ = fMuonCosQ;
196  hAngles = fMuonAngles;
197  hAnglesLo = fMuonAnglesLo;
198  hAnglesMi = fMuonAnglesMi;
199  hAnglesHi = fMuonAnglesHi;
200  hEnergy = fMuonEnergy;
201  }
202  else if (std::abs(particle.PdgCode())==22) {
203  hCosQ = fPhotonCosQ;
204  hAngles = fPhotonAngles;
205  hAnglesLo = fPhotonAnglesLo;
206  hAnglesMi = fPhotonAnglesMi;
207  hAnglesHi = fPhotonAnglesHi;
208  hEnergy = fPhotonEnergy;
209  }
210  else if (std::abs(particle.PdgCode())==11) {
211  hCosQ = fElectronCosQ;
212  hAngles = fElectronAngles;
213  hAnglesLo = fElectronAnglesLo;
214  hAnglesMi = fElectronAnglesMi;
215  hAnglesHi = fElectronAnglesHi;
216  hEnergy = fElectronEnergy;
217  }
218 
219  // now check if the particle goes through any cryostat in the detector
220  // if so, add it to the truth object.
221  for(unsigned int c = 0; c < geom->Ncryostats(); ++c){
222 
223  double bounds[6] = {0.};
224  geom->CryostatBoundaries(bounds, c);
225 
226  //add a buffer box around the cryostat bounds to increase the acceptance
227  //(geometrically) at the CRY level to make up for particles we will loose
228  //due to multiple scattering effects that pitch in during GEANT4 tracking
229  //By default, the buffer box has zero size
230  for (unsigned int cb=0; cb<6; cb++)
231  bounds[cb] = bounds[cb]+fbuffbox[cb];
232 
233  //calculate the intersection point with each cryostat surface
234  bool intersects_cryo = false;
235  for (int bnd=0; bnd!=6; ++bnd) {
236  if (bnd<2) {
237  double p2[3] = {bounds[bnd], x0[1] + (dx[1]/dx[0])*(bounds[bnd] - x0[0]), x0[2] + (dx[2]/dx[0])*(bounds[bnd] - x0[0])};
238  if ( p2[1] >= bounds[2] && p2[1] <= bounds[3] &&
239  p2[2] >= bounds[4] && p2[2] <= bounds[5] ) {
240  intersects_cryo = true;
241  break;
242  }
243  }
244  else if (bnd>=2 && bnd<4) {
245  double p2[3] = {x0[0] + (dx[0]/dx[1])*(bounds[bnd] - x0[1]), bounds[bnd], x0[2] + (dx[2]/dx[1])*(bounds[bnd] - x0[1])};
246  if ( p2[0] >= bounds[0] && p2[0] <= bounds[1] &&
247  p2[2] >= bounds[4] && p2[2] <= bounds[5] ) {
248  intersects_cryo = true;
249  break;
250  }
251  }
252  else if (bnd>=4) {
253  double p2[3] = {x0[0] + (dx[0]/dx[2])*(bounds[bnd] - x0[2]), x0[1] + (dx[1]/dx[2])*(bounds[bnd] - x0[2]), bounds[bnd]};
254  if ( p2[0] >= bounds[0] && p2[0] <= bounds[1] &&
255  p2[1] >= bounds[2] && p2[1] <= bounds[3] ) {
256  intersects_cryo = true;
257  break;
258  }
259  }
260  }
261 
262 
263  if (intersects_cryo) {
264  truth.Add(particle);
265 
266  if (std::abs(particle.PdgCode())==13) ++numMuons;
267  else if (std::abs(particle.PdgCode())==22) ++numPhotons;
268  else if (std::abs(particle.PdgCode())==11) ++numElectrons;
269 
270  //The following code no longer works now that we require intersection with the cryostat boundary
271  //For example, the particle could intersect this cryostat but miss its TPC, but intersect a TPC
272  //in another cryostat
273  /*try{
274  unsigned int tpc = 0;
275  unsigned int cstat = 0;
276  geom->PositionToTPC(x2, tpc, cstat);
277  if (std::abs(particle.PdgCode())==13) ++tpcMuons;
278  else if (std::abs(particle.PdgCode())==22) ++tpcPhotons;
279  else if (std::abs(particle.PdgCode())==11) ++tpcElectrons;
280  }
281  catch(cet::exception &e){
282  MF_LOG_DEBUG("CosmicsGen") << "current particle does not go through any tpc";
283  }*///
284 
285  if (hCosQ!=0) {
286  double cosq = -p4.Py()/p4.P();
287  double phi = std::atan2(p4.Pz(),p4.Px());
288  phi *= 180/M_PI;
289  hCosQ->Fill(cosq);
290  hAngles->Fill(phi,cosq);
291  if (p4.E()<1.0) hAnglesLo->Fill(phi,cosq);
292  else if (p4.E()<10.0) hAnglesMi->Fill(phi,cosq);
293  else hAnglesHi->Fill(phi,cosq);
294  hEnergy->Fill(p4.E());
295  }//end if there is a cos(theta) histogram
296  break; //leave loop over cryostats to avoid adding particle multiple times
297  }// end if particle goes into a cryostat
298  }// end loop over cryostats in the detector
299 
300  }// loop on particles
301 
302  nCrossCryostat = truth.NParticles();
303 
304  fPhotonsPerSample ->Fill(allPhotons);
305  fElectronsPerSample->Fill(allElectrons);
306  fMuonsPerSample ->Fill(allMuons);
307 
308  fPhotonsInCStat ->Fill(numPhotons);
309  fElectronsInCStat->Fill(numElectrons);
310  fMuonsInCStat ->Fill(numMuons);
311 
312  /*fPhotonsInTPC ->Fill(tpcPhotons);
313  fElectronsInTPC->Fill(tpcElectrons);
314  fMuonsInTPC ->Fill(tpcMuons);*/
315  }
316 
317  truthcol->push_back(truth);
318  evt.put(std::move(truthcol));
319  }// end produce
320 
321 }// end namespace
322 
323 
324 namespace evgen{
325 
327 
328 }
TH2F * fElectronAnglesHi
Electron rate vs angle, high momenta.
void beginJob() override
Interface to the CRY cosmic-ray generator.
Definition: CRYHelper.h:26
const TLorentzVector & Position(const int i=0) const
Definition: MCParticle.h:219
TH2F * fPhotonAngles
Photon rate vs angle.
int PdgCode() const
Definition: MCParticle.h:212
TH2F * fMuonAnglesHi
Muon rate vs angle, high momenta.
TH1F * fElectronEnergy
Electron energy (GeV)
TH2F * fElectronAnglesLo
Electron rate vs angle, low momenta.
TH2F * fElectronAngles
Electron rate vs angle.
TH1F * fElectronCosQ
Electron rate vs cos(Q)
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:82
TH2F * fMuonAngles
Muon rate vs angle.
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
TH2F * fMuonAnglesLo
Muon rate vs angle, low momenta.
TH1F * fPhotonEnergy
Photon energy (GeV)
const std::string GetWorldVolumeName() const
Return the name of the world volume (needed by Geant4 simulation)
int NParticles() const
Definition: MCTruth.h:75
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
TH2F * fElectronAnglesMi
Electron rate vs angle, middle momenta.
Particle class.
void produce(art::Event &evt) override
TH1F * fMuonCosQ
Muon rate vs cos(Q)
TH1F * fElectronsPerSample
number of electrons in the sampled time window
void CryostatBoundaries(double *boundaries, geo::CryostatID const &cid) const
Returns the boundaries of the specified cryostat.
art framework interface to geometry description
Definition: Run.h:17
CosmicsGen(fhicl::ParameterSet const &pset)
TH1F * fPhotonsPerSample
number of photons in the sampled time window
T abs(T value)
double Sample(simb::MCTruth &mctruth, double const &surfaceY, double const &detectorLength, double *w, double rantime=0)
Definition: CRYHelper.cxx:97
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
Interface to the CRY cosmic ray generator.
def move(depos, offset)
Definition: depos.py:107
TH1F * fPhotonCosQ
Photon rate vs cos(Q)
CLHEP::HepRandomEngine & fEngine
art-managed random-number engine
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
void beginRun(art::Run &run) override
geo::Length_t SurfaceY() const
The position of the detector respect to earth surface.
evgb::CRYHelper fCRYHelp
CRY generator object.
BoundingBox bounds(int x, int y, int w, int h)
Definition: main.cpp:37
TH2F * fMuonAnglesMi
Muon rate vs angle, middle momenta.
std::vector< double > fbuffbox
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:76
TH2F * fPhotonAnglesMi
Photon rate vs angle, middle momenta.
TH2F * fPhotonAnglesLo
Photon rate vs angle, low momenta.
#define M_PI
Definition: includeROOT.h:54
void Add(simb::MCParticle const &part)
Definition: MCTruth.h:80
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:220
TH2F * fPhotonAnglesHi
Photon rate vs angle, high momenta.
A module to check the results from the Monte Carlo generator.
TCEvent evt
Definition: DataStructs.cxx:7
Event generator information.
Definition: MCTruth.h:32
TH1F * fMuonsPerSample
number of muons in the sampled time window
LArSoft geometry interface.
Definition: ChannelGeo.h:16
TH1F * fMuonEnergy
Muon energy (GeV)
Event Generation using GENIE, cosmics or single particles.
Cosmic rays.
Definition: MCTruth.h:24