LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
CosmicsGen_module.cc
Go to the documentation of this file.
1 
10 // ROOT includes
11 #include "TH1.h"
12 #include "TH2.h"
13 
14 // Framework includes
20 #include "art_root_io/TFileService.h"
21 #include "fhiclcpp/ParameterSet.h"
22 
23 // art extensions
25 
26 // larsoft includes
33 
34 namespace evgen {
35 
37  class CosmicsGen : public art::EDProducer {
38  public:
39  explicit CosmicsGen(fhicl::ParameterSet const& pset);
40 
41  private:
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;
52  TH1F* fPhotonCosQ;
53  TH1F* fPhotonEnergy;
56  TH1F* fPhotonsInTPC;
58 
64  TH1F* fElectronCosQ;
68  TH1F* fElectronsInTPC;
70 
72  TH2F* fMuonAngles;
73  TH2F* fMuonAnglesLo;
74  TH2F* fMuonAnglesMi;
75  TH2F* fMuonAnglesHi;
76  TH1F* fMuonCosQ;
77  TH1F* fMuonEnergy;
79  TH1F* fMuonsInCStat;
80  TH1F* fMuonsInTPC;
82  CLHEP::HepRandomEngine& fEngine;
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"
97  pset,
98  "Seed"))
99  , fCRYHelp{pset, fEngine, art::ServiceHandle<geo::Geometry const>{} -> GetWorldVolumeName()}
100  {
101  produces<std::vector<simb::MCTruth>>();
102  produces<sumdata::RunData, art::InRun>();
103  }
104 
105  //____________________________________________________________________________
107  {
109 
110  fPhotonAngles =
111  tfs->make<TH2F>("fPhotonAngles", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
113  tfs->make<TH2F>("fPhotonAnglesLo", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
115  tfs->make<TH2F>("fPhotonAnglesMi", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
117  tfs->make<TH2F>("fPhotonAnglesHi", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
118  fPhotonCosQ = tfs->make<TH1F>("fPhotonCosQ", ";cos#theta;tracks", 50, -1.0, 1.0);
119  fPhotonEnergy = tfs->make<TH1F>("fPhotonEnergy", ";E (GeV)", 5000, 0.0, 1000.0);
121  tfs->make<TH1F>("fPhotonsPerSample", ";Number Photons;Samples", 100, 0, 1000);
123  tfs->make<TH1F>("fPhotonsInCryostat", ";Number Photons;Samples", 100, 0, 1000);
124  fPhotonsInTPC = tfs->make<TH1F>("fPhotonsInTPC", ";Number Photons;Samples", 100, 0, 1000);
125 
127  tfs->make<TH2F>("fElectronAngles", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
129  tfs->make<TH2F>("fElectronAnglesLo", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
131  tfs->make<TH2F>("fElectronAnglesMi", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
133  tfs->make<TH2F>("fElectronAnglesHi", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
134  fElectronCosQ = tfs->make<TH1F>("fElectronCosQ", ";cos#theta;tracks", 50, -1.0, 1.0);
135  fElectronEnergy = tfs->make<TH1F>("fElectronEnergy", ";E (GeV)", 5000, 0.0, 1000.0);
137  tfs->make<TH1F>("fElectronsPerSample", ";Number Electrons;Samples", 100, 0, 1000);
139  tfs->make<TH1F>("fElectronsInCryotat", ";Number Electrons;Samples", 100, 0, 1000);
140  fElectronsInTPC = tfs->make<TH1F>("fElectronsInTPC", ";Number Electrons;Samples", 100, 0, 1000);
141 
142  fMuonAngles =
143  tfs->make<TH2F>("fMuonAngles", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
144  fMuonAnglesLo =
145  tfs->make<TH2F>("fMuonAnglesLo", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
146  fMuonAnglesMi =
147  tfs->make<TH2F>("fMuonAnglesMi", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
148  fMuonAnglesHi =
149  tfs->make<TH2F>("fMuonAnglesHi", ";#phi;cos#theta", 36, -180.0, 180.0, 50, -1.0, 1.0);
150  fMuonCosQ = tfs->make<TH1F>("fMuonCosQ", ";cos#theta;tracks", 50, -1.0, 1.0);
151  fMuonEnergy = tfs->make<TH1F>("fMuonEnergy", ";E (GeV)", 5000, 0.0, 1000.0);
152  fMuonsPerSample = tfs->make<TH1F>("fMuonsPerSample", ";Number Muons;Samples", 100, 0, 1000);
153  fMuonsInCStat = tfs->make<TH1F>("fMuonsInCryostat", ";Number Muons;Samples", 100, 0, 1000);
154  fMuonsInTPC = tfs->make<TH1F>("fMuonsInTPC", ";Number Muons;Samples", 100, 0, 1000);
155  }
156 
157  //____________________________________________________________________________
159  {
161  run.put(std::make_unique<sumdata::RunData>(geo->DetectorName()), art::fullRun());
162  }
163 
164  //____________________________________________________________________________
166  {
167  std::unique_ptr<std::vector<simb::MCTruth>> truthcol(new std::vector<simb::MCTruth>);
168 
169  // fill some histograms about this event
171 
172  int nCrossCryostat = 0;
173 
174  simb::MCTruth truth;
175 
176  while (nCrossCryostat < 1) {
177 
178  simb::MCTruth pretruth;
180  fCRYHelp.Sample(pretruth, geom->SurfaceY(), geom->DetLength(), 0);
181 
182  int numPhotons = 0;
183  int numElectrons = 0;
184  int numMuons = 0;
185  int allPhotons = 0;
186  int allElectrons = 0;
187  int allMuons = 0;
188 
189  // loop over particles in the truth object
190  for (int i = 0; i < pretruth.NParticles(); ++i) {
191  simb::MCParticle particle = pretruth.GetParticle(i);
192  const TLorentzVector& v4 = particle.Position();
193  const TLorentzVector& p4 = particle.Momentum();
194  double x0[3] = {v4.X(), v4.Y(), v4.Z()};
195  double dx[3] = {p4.Px(), p4.Py(), p4.Pz()};
196 
197  if (std::abs(particle.PdgCode()) == 13)
198  ++allMuons;
199  else if (std::abs(particle.PdgCode()) == 22)
200  ++allPhotons;
201  else if (std::abs(particle.PdgCode()) == 11)
202  ++allElectrons;
203 
204  TH1F* hCosQ = nullptr;
205  TH2F* hAngles = nullptr;
206  TH2F* hAnglesLo = nullptr;
207  TH2F* hAnglesMi = nullptr;
208  TH2F* hAnglesHi = nullptr;
209  TH1F* hEnergy = nullptr;
210  if (std::abs(particle.PdgCode()) == 13) {
211  hCosQ = fMuonCosQ;
212  hAngles = fMuonAngles;
213  hAnglesLo = fMuonAnglesLo;
214  hAnglesMi = fMuonAnglesMi;
215  hAnglesHi = fMuonAnglesHi;
216  hEnergy = fMuonEnergy;
217  }
218  else if (std::abs(particle.PdgCode()) == 22) {
219  hCosQ = fPhotonCosQ;
220  hAngles = fPhotonAngles;
221  hAnglesLo = fPhotonAnglesLo;
222  hAnglesMi = fPhotonAnglesMi;
223  hAnglesHi = fPhotonAnglesHi;
224  hEnergy = fPhotonEnergy;
225  }
226  else if (std::abs(particle.PdgCode()) == 11) {
227  hCosQ = fElectronCosQ;
228  hAngles = fElectronAngles;
229  hAnglesLo = fElectronAnglesLo;
230  hAnglesMi = fElectronAnglesMi;
231  hAnglesHi = fElectronAnglesHi;
232  hEnergy = fElectronEnergy;
233  }
234 
235  // now check if the particle goes through any cryostat in the detector
236  // if so, add it to the truth object.
237  for (auto const& cryostat : geom->Iterate<geo::CryostatGeo>()) {
238 
239  double bounds[6] = {0.};
240  cryostat.Boundaries(bounds);
241 
242  //add a buffer box around the cryostat bounds to increase the acceptance
243  //(geometrically) at the CRY level to make up for particles we will loose
244  //due to multiple scattering effects that pitch in during GEANT4 tracking
245  //By default, the buffer box has zero size
246  for (unsigned int cb = 0; cb < 6; cb++)
247  bounds[cb] = bounds[cb] + fbuffbox[cb];
248 
249  //calculate the intersection point with each cryostat surface
250  bool intersects_cryo = false;
251  for (int bnd = 0; bnd != 6; ++bnd) {
252  if (bnd < 2) {
253  double p2[3] = {bounds[bnd],
254  x0[1] + (dx[1] / dx[0]) * (bounds[bnd] - x0[0]),
255  x0[2] + (dx[2] / dx[0]) * (bounds[bnd] - x0[0])};
256  if (p2[1] >= bounds[2] && p2[1] <= bounds[3] && p2[2] >= bounds[4] &&
257  p2[2] <= bounds[5]) {
258  intersects_cryo = true;
259  break;
260  }
261  }
262  else if (bnd >= 2 && bnd < 4) {
263  double p2[3] = {x0[0] + (dx[0] / dx[1]) * (bounds[bnd] - x0[1]),
264  bounds[bnd],
265  x0[2] + (dx[2] / dx[1]) * (bounds[bnd] - x0[1])};
266  if (p2[0] >= bounds[0] && p2[0] <= bounds[1] && p2[2] >= bounds[4] &&
267  p2[2] <= bounds[5]) {
268  intersects_cryo = true;
269  break;
270  }
271  }
272  else if (bnd >= 4) {
273  double p2[3] = {x0[0] + (dx[0] / dx[2]) * (bounds[bnd] - x0[2]),
274  x0[1] + (dx[1] / dx[2]) * (bounds[bnd] - x0[2]),
275  bounds[bnd]};
276  if (p2[0] >= bounds[0] && p2[0] <= bounds[1] && p2[1] >= bounds[2] &&
277  p2[1] <= bounds[3]) {
278  intersects_cryo = true;
279  break;
280  }
281  }
282  }
283 
284  if (intersects_cryo) {
285  truth.Add(particle);
286 
287  if (std::abs(particle.PdgCode()) == 13)
288  ++numMuons;
289  else if (std::abs(particle.PdgCode()) == 22)
290  ++numPhotons;
291  else if (std::abs(particle.PdgCode()) == 11)
292  ++numElectrons;
293 
294  if (hCosQ != 0) {
295  double cosq = -p4.Py() / p4.P();
296  double phi = std::atan2(p4.Pz(), p4.Px());
297  phi *= 180 / M_PI;
298  hCosQ->Fill(cosq);
299  hAngles->Fill(phi, cosq);
300  if (p4.E() < 1.0)
301  hAnglesLo->Fill(phi, cosq);
302  else if (p4.E() < 10.0)
303  hAnglesMi->Fill(phi, cosq);
304  else
305  hAnglesHi->Fill(phi, cosq);
306  hEnergy->Fill(p4.E());
307  } //end if there is a cos(theta) histogram
308  break; //leave loop over cryostats to avoid adding particle multiple times
309  } // end if particle goes into a cryostat
310  } // end loop over cryostats in the detector
311 
312  } // loop on particles
313 
314  nCrossCryostat = truth.NParticles();
315 
316  fPhotonsPerSample->Fill(allPhotons);
317  fElectronsPerSample->Fill(allElectrons);
318  fMuonsPerSample->Fill(allMuons);
319 
320  fPhotonsInCStat->Fill(numPhotons);
321  fElectronsInCStat->Fill(numElectrons);
322  fMuonsInCStat->Fill(numMuons);
323  }
324 
325  truthcol->push_back(truth);
326  evt.put(std::move(truthcol));
327  } // end produce
328 
329 } // end namespace
330 
TH2F * fElectronAnglesHi
Electron rate vs angle, high momenta.
details::range_type< T > Iterate() const
Initializes the specified ID with the ID of the first cryostat.
Definition: GeometryCore.h:541
base_engine_t & createEngine(seed_t seed)
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:220
TH2F * fPhotonAngles
Photon rate vs angle.
int PdgCode() const
Definition: MCParticle.h:213
TH2F * fMuonAnglesHi
Muon rate vs angle, high momenta.
constexpr auto fullRun()
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.cc:6
TH2F * fMuonAnglesLo
Muon rate vs angle, low momenta.
TH1F * fPhotonEnergy
Photon energy (GeV)
constexpr auto abs(T v)
Returns the absolute value of the argument.
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Run.h:121
int NParticles() const
Definition: MCTruth.h:75
Geometry information for a single cryostat.
Definition: CryostatGeo.h:43
TH2F * fElectronAnglesMi
Electron rate vs angle, middle momenta.
Particle class.
void produce(art::Event &evt) override
TH1F * fMuonCosQ
Muon rate vs cos(Q)
Length_t DetLength(TPCID const &tpcid=tpc_zero) const
Returns the length of the active volume of the specified TPC.
Length_t SurfaceY() const
The position of the detector respect to earth surface.
Definition: GeometryCore.h:251
TH1F * fElectronsPerSample
number of electrons in the sampled time window
Definition: Run.h:37
CosmicsGen(fhicl::ParameterSet const &pset)
TH1F * fPhotonsPerSample
number of photons in the sampled time window
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
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:65
Interface to the CRY cosmic ray generator.
TH1F * fPhotonCosQ
Photon rate vs cos(Q)
CLHEP::HepRandomEngine & fEngine
art-managed random-number engine
void beginRun(art::Run &run) override
evgb::CRYHelper fCRYHelp
CRY generator object.
An art service to assist in the distribution of guaranteed unique seeds to all engines within an art ...
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.
void Add(simb::MCParticle const &part)
Definition: MCTruth.h:80
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:221
Definition: MVAAlg.h:12
TH2F * fPhotonAnglesHi
Photon rate vs angle, high momenta.
A module to check the results from the Monte Carlo generator.
std::string const & DetectorName() const
Returns a string with the name of the detector, as configured.
Definition: GeometryCore.h:203
TCEvent evt
Definition: DataStructs.cxx:8
Event generator information.
Definition: MCTruth.h:32
TH1F * fMuonsPerSample
number of muons in the sampled time window
Namespace collecting geometry-related classes utilities.
TH1F * fMuonEnergy
Muon energy (GeV)
Event Generation using GENIE, cosmics or single particles.
art framework interface to geometry description
Cosmic rays.
Definition: MCTruth.h:24