LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
GENIEGen_module.cc
Go to the documentation of this file.
1 //
3 //
4 // GENIE neutrino event generator
5 //
6 // brebel@fnal.gov
7 //
9 #ifndef EVGEN_GENIEGEN_H
10 #define EVGEN_GENIEGEN_H
11 
12 #include <cstdlib>
13 #include <string>
14 #include <iostream>
15 #include <iomanip>
16 #include <sstream>
17 #include <vector>
18 #include <map>
19 #include <memory>
20 #include <unistd.h>
21 
22 // ROOT includes
23 #include "TH1.h"
24 #include "TH2.h"
25 #include "TDatabasePDG.h"
26 #include "TSystem.h"
27 #include "TStopwatch.h"
28 
29 // Framework includes
33 #include "fhiclcpp/ParameterSet.h"
41 
42 // art extensions
44 
45 // LArSoft includes
46 #include "lardataalg/MCDumpers/MCDumpers.h" // sim::dump namespace
56 
58 namespace evgen {
90  class GENIEGen : public art::EDProducer {
91  public:
92  explicit GENIEGen(fhicl::ParameterSet const &pset);
93  virtual ~GENIEGen();
94 
95  void produce(art::Event& evt);
96  void beginJob();
97  void beginRun(art::Run& run);
98  void beginSubRun(art::SubRun& sr);
99  void endSubRun(art::SubRun& sr);
100 
101  private:
102 
103  std::string ParticleStatus(int StatusCode);
104  std::string ReactionChannel(int ccnc,int mode);
105 
106  void FillHistograms(simb::MCTruth mc);
107 
110  std::vector< double > fVtxPosHistRange;
111 
113  TStopwatch fStopwatch;
114 
118 
119  double fPrevTotPOT;
121 
122  TH1F* fGenerated[6];
123 
124  TH1F* fVertexX;
125  TH1F* fVertexY;
126  TH1F* fVertexZ;
127 
128  TH2F* fVertexXY;
129  TH2F* fVertexXZ;
130  TH2F* fVertexYZ;
131 
132  TH1F* fDCosX;
133  TH1F* fDCosY;
134  TH1F* fDCosZ;
135 
136  TH1F* fMuMomentum;
137  TH1F* fMuDCosX;
138  TH1F* fMuDCosY;
139  TH1F* fMuDCosZ;
140 
141  TH1F* fEMomentum;
142  TH1F* fEDCosX;
143  TH1F* fEDCosY;
144  TH1F* fEDCosZ;
145 
146  TH1F* fCCMode;
147  TH1F* fNCMode;
148 
149  TH1F* fDeltaE;
150  TH1F* fECons;
151 
152  };
153 }
154 
155 namespace evgen{
156 
157  //____________________________________________________________________________
159  : fGENIEHelp(0)
160  , fDefinedVtxHistRange (pset.get< bool >("DefinedVtxHistRange"))
161  , fVtxPosHistRange (pset.get< std::vector<double> >("VtxPosHistRange"))
162  , fPassEmptySpills (pset.get< bool >("PassEmptySpills"))
163  , fGlobalTimeOffset(pset.get< double >("GlobalTimeOffset",0))
164  , fRandomTimeOffset(pset.get< double >("RandomTimeOffset",1600.)) // BNB default value
165  , fBeamType(::sim::kBNB)
166  {
167  fStopwatch.Start();
168 
169  produces< std::vector<simb::MCTruth> >();
170  produces< std::vector<simb::MCFlux> >();
171  produces< std::vector<simb::GTruth> >();
172  produces< sumdata::RunData, art::InRun >();
173  produces< sumdata::POTSummary, art::InSubRun >();
174  produces< art::Assns<simb::MCTruth, simb::MCFlux> >();
175  produces< art::Assns<simb::MCTruth, simb::GTruth> >();
176  produces< std::vector<sim::BeamGateInfo> >();
177 
178  std::string beam_type_name = pset.get<std::string>("BeamName");
179 
180  if(beam_type_name == "numi")
181 
183 
184  else if(beam_type_name == "booster")
185 
187 
188  else
189 
191 
193 
194  signed int temp_seed; // the seed read by GENIEHelper is a signed integer...
195  fhicl::ParameterSet GENIEconfig(pset);
196  if (!GENIEconfig.get_if_present("RandomSeed", temp_seed)) { // TODO use has_key() when it becomes available
197  // no RandomSeed specified; check for the LArSoft-style "Seed" instead:
198  // obtain the random seed from a service,
199  // unless overridden in configuration with key "Seed"
200  unsigned int seed;
201  if (!GENIEconfig.get_if_present("Seed", seed))
202  seed = art::ServiceHandle<rndm::NuRandomService>()->getSeed();
203 
204  // The seed is not passed to RandomNumberGenerator,
205  // since GENIE uses a TRandom generator that is owned by the GENIEHelper.
206  // Instead, we explicitly configure the random seed for GENIEHelper:
207  GENIEconfig.put("RandomSeed", seed);
208  } // if no RandomSeed present
209 
210  fGENIEHelp = new evgb::GENIEHelper(GENIEconfig,
211  geo->ROOTGeoManager(),
212  geo->ROOTFile(),
213  geo->TotalMass(pset.get< std::string>("TopVolume").c_str()));
214 
215  }
216 
217  //____________________________________________________________________________
219  {
220  if(fGENIEHelp) delete fGENIEHelp;
221  fStopwatch.Stop();
222  mf::LogInfo("GENIEProductionTime") << "real time to produce file: " << fStopwatch.RealTime();
223  }
224 
225  //____________________________________________________________________________
228 
229  fPrevTotPOT = 0.;
230  fPrevTotGoodPOT = 0.;
231 
232  // Get access to the TFile service.
234 
235  fGenerated[0] = tfs->make<TH1F>("fGenerated_necc","", 100, 0.0, 20.0);
236  fGenerated[1] = tfs->make<TH1F>("fGenerated_nebcc","", 100, 0.0, 20.0);
237  fGenerated[2] = tfs->make<TH1F>("fGenerated_nmcc","", 100, 0.0, 20.0);
238  fGenerated[3] = tfs->make<TH1F>("fGenerated_nmbcc","", 100, 0.0, 20.0);
239  fGenerated[4] = tfs->make<TH1F>("fGenerated_nnc","", 100, 0.0, 20.0);
240  fGenerated[5] = tfs->make<TH1F>("fGenerated_nbnc","", 100, 0.0, 20.0);
241 
242  fDCosX = tfs->make<TH1F>("fDCosX", ";dx/ds", 200, -1., 1.);
243  fDCosY = tfs->make<TH1F>("fDCosY", ";dy/ds", 200, -1., 1.);
244  fDCosZ = tfs->make<TH1F>("fDCosZ", ";dz/ds", 200, -1., 1.);
245 
246  fMuMomentum = tfs->make<TH1F>("fMuMomentum", ";p_{#mu} (GeV/c)", 500, 0., 50.);
247  fMuDCosX = tfs->make<TH1F>("fMuDCosX", ";dx/ds;", 200, -1., 1.);
248  fMuDCosY = tfs->make<TH1F>("fMuDCosY", ";dy/ds;", 200, -1., 1.);
249  fMuDCosZ = tfs->make<TH1F>("fMuDCosZ", ";dz/ds;", 200, -1., 1.);
250 
251  fEMomentum = tfs->make<TH1F>("fEMomentum", ";p_{e} (GeV/c)", 500, 0., 50.);
252  fEDCosX = tfs->make<TH1F>("fEDCosX", ";dx/ds;", 200, -1., 1.);
253  fEDCosY = tfs->make<TH1F>("fEDCosY", ";dy/ds;", 200, -1., 1.);
254  fEDCosZ = tfs->make<TH1F>("fEDCosZ", ";dz/ds;", 200, -1., 1.);
255 
256  fCCMode = tfs->make<TH1F>("fCCMode", ";CC Interaction Mode;", 4, 0., 4.);
257  fCCMode->GetXaxis()->SetBinLabel(1, "QE");
258  fCCMode->GetXaxis()->SetBinLabel(2, "Res");
259  fCCMode->GetXaxis()->SetBinLabel(3, "DIS");
260  fCCMode->GetXaxis()->SetBinLabel(4, "Coh");
261  fCCMode->GetXaxis()->CenterLabels();
262 
263  fNCMode = tfs->make<TH1F>("fNCMode", ";NC Interaction Mode;", 4, 0., 4.);
264  fNCMode->GetXaxis()->SetBinLabel(1, "QE");
265  fNCMode->GetXaxis()->SetBinLabel(2, "Res");
266  fNCMode->GetXaxis()->SetBinLabel(3, "DIS");
267  fNCMode->GetXaxis()->SetBinLabel(4, "Coh");
268  fNCMode->GetXaxis()->CenterLabels();
269 
270  fDeltaE = tfs->make<TH1F>("fDeltaE", ";#Delta E_{#nu} (GeV);", 200, -1., 1.);
271  fECons = tfs->make<TH1F>("fECons", ";#Delta E(#nu,lepton);", 500, -5., 5.);
272 
274  double x = 2.1*geo->DetHalfWidth();
275  double y = 2.1*geo->DetHalfHeight();
276  double z = 2.*geo->DetLength();
277  int xdiv = TMath::Nint(2*x/5.);
278  int ydiv = TMath::Nint(2*y/5.);
279  int zdiv = TMath::Nint(2*z/5.);
280 
281  if (fDefinedVtxHistRange == false)
282  {
283  fVertexX = tfs->make<TH1F>("fVertexX", ";x (cm)", xdiv, -0.1*x, x);
284  fVertexY = tfs->make<TH1F>("fVertexY", ";y (cm)", ydiv, -y, y);
285  fVertexZ = tfs->make<TH1F>("fVertexZ", ";z (cm)", zdiv, -0.1*z, z);
286 
287  fVertexXY = tfs->make<TH2F>("fVertexXY", ";x (cm);y (cm)", xdiv, -0.1*x, x, ydiv, -y, y);
288  fVertexXZ = tfs->make<TH2F>("fVertexXZ", ";z (cm);x (cm)", zdiv, -0.2*z, z, xdiv, -0.1*x, x);
289  fVertexYZ = tfs->make<TH2F>("fVertexYZ", ";z (cm);y (cm)", zdiv, -0.2*z, z, ydiv, -y, y);
290  }
291  else
292  {
293  fVertexX = tfs->make<TH1F>("fVertexX", ";x (cm)", xdiv, fVtxPosHistRange[0], fVtxPosHistRange[1]);
294  fVertexY = tfs->make<TH1F>("fVertexY", ";y (cm)", ydiv, fVtxPosHistRange[2], fVtxPosHistRange[3]);
295  fVertexZ = tfs->make<TH1F>("fVertexZ", ";z (cm)", zdiv, fVtxPosHistRange[4], fVtxPosHistRange[5]);
296 
297  fVertexXY = tfs->make<TH2F>("fVertexXY", ";x (cm);y (cm)", xdiv, fVtxPosHistRange[0], fVtxPosHistRange[1], ydiv, fVtxPosHistRange[2], fVtxPosHistRange[3]);
298  fVertexXZ = tfs->make<TH2F>("fVertexXZ", ";z (cm);x (cm)", zdiv, fVtxPosHistRange[4], fVtxPosHistRange[5], xdiv, fVtxPosHistRange[0], fVtxPosHistRange[1]);
299  fVertexYZ = tfs->make<TH2F>("fVertexYZ", ";z (cm);y (cm)", zdiv, fVtxPosHistRange[4], fVtxPosHistRange[5], ydiv, fVtxPosHistRange[2], fVtxPosHistRange[3]);
300  }
301 
302  }
303 
304  //____________________________________________________________________________
306  {
307 
308  // grab the geometry object to see what geometry we are using
310  std::unique_ptr<sumdata::RunData> runcol(new sumdata::RunData(geo->DetectorName()));
311 
312  run.put(std::move(runcol));
313 
314  return;
315  }
316 
317  //____________________________________________________________________________
319  {
320 
323 
324  return;
325  }
326 
327  //____________________________________________________________________________
329  {
330 
331  std::unique_ptr<sumdata::POTSummary> p(new sumdata::POTSummary());
332 
335 
336  sr.put(std::move(p));
337 
338  return;
339  }
340 
341  //____________________________________________________________________________
343  {
344  std::unique_ptr< std::vector<simb::MCTruth> > truthcol (new std::vector<simb::MCTruth>);
345  std::unique_ptr< std::vector<simb::MCFlux> > fluxcol (new std::vector<simb::MCFlux >);
346  std::unique_ptr< std::vector<simb::GTruth> > gtruthcol (new std::vector<simb::GTruth >);
347  std::unique_ptr< art::Assns<simb::MCTruth, simb::MCFlux> > tfassn(new art::Assns<simb::MCTruth, simb::MCFlux>);
348  std::unique_ptr< art::Assns<simb::MCTruth, simb::GTruth> > tgtassn(new art::Assns<simb::MCTruth, simb::GTruth>);
349  std::unique_ptr< std::vector<sim::BeamGateInfo> > gateCollection(new std::vector<sim::BeamGateInfo>);
350 
351  while(truthcol->size() < 1){
352  while(!fGENIEHelp->Stop()){
353 
354  simb::MCTruth truth;
355  simb::MCFlux flux;
356  simb::GTruth gTruth;
357 
358  // GENIEHelper returns a false in the sample method if
359  // either no neutrino was generated, or the interaction
360  // occurred beyond the detector's z extent - ie something we
361  // would never see anyway.
362  if(fGENIEHelp->Sample(truth, flux, gTruth)){
363 
364  truthcol ->push_back(truth);
365  fluxcol ->push_back(flux);
366  gtruthcol->push_back(gTruth);
367  util::CreateAssn(*this, evt, *truthcol, *fluxcol, *tfassn, fluxcol->size()-1, fluxcol->size());
368  util::CreateAssn(*this, evt, *truthcol, *gtruthcol, *tgtassn, gtruthcol->size()-1, gtruthcol->size());
369 
370  FillHistograms(truth);
371 
372  // check that the process code is not unsupported by GENIE
373  // (see issue #18025 for reference);
374  // if it is, print all the information we can about this truth record
375  if (truth.NeutrinoSet() && (truth.GetNeutrino().InteractionType() == simb::kNuanceOffset)) {
376  mf::LogWarning log("GENIEmissingProcessMapping");
377  log << "Found an interaction that is not represented by the interaction type code in GENIE:"
378  "\nMCTruth record:"
379  "\n"
380  ;
381  sim::dump::DumpMCTruth(log, truth, 2U); // 2 trajectory points per line
382  log <<
383  "\nGENIE truth record:"
384  "\n"
385  ;
386  sim::dump::DumpGTruth(log, gTruth);
387  } // if
388 
389  }// end if genie was able to make an event
390 
391  }// end event generation loop
392 
393  // check to see if we are to pass empty spills
394  if(truthcol->size() < 1 && fPassEmptySpills){
395  LOG_DEBUG("GENIEGen") << "no events made for this spill but "
396  << "passing it on and ending the event anyway";
397  break;
398  }
399 
400  }// end loop while no interactions are made
401 
402  // Create a simulated "beam gate" for these neutrino events.
403  // We're creating a vector of these because, in a
404  // distant-but-possible future, we may be generating more than one
405  // beam gate within a simulated time window.
406  gateCollection->push_back(sim::BeamGateInfo( fGlobalTimeOffset, fRandomTimeOffset, fBeamType ));
407 
408  // put the collections in the event
409  evt.put(std::move(truthcol));
410  evt.put(std::move(fluxcol));
411  evt.put(std::move(gtruthcol));
412  evt.put(std::move(tfassn));
413  evt.put(std::move(tgtassn));
414  evt.put(std::move(gateCollection));
415 
416  return;
417  }
418 
419  //......................................................................
420  std::string GENIEGen::ParticleStatus(int StatusCode)
421  {
422  int code = StatusCode;
423  std::string ParticleStatusName;
424 
425  switch(code)
426  {
427  case -1:
428  ParticleStatusName = "kIStUndefined";
429  break;
430  case 0:
431  ParticleStatusName = "kIStInitialState";
432  break;
433  case 1:
434  ParticleStatusName = "kIStStableFinalState";
435  break;
436  case 2:
437  ParticleStatusName = "kIStIntermediateState";
438  break;
439  case 3:
440  ParticleStatusName = "kIStDecayedState";
441  break;
442  case 11:
443  ParticleStatusName = "kIStNucleonTarget";
444  break;
445  case 12:
446  ParticleStatusName = "kIStDISPreFragmHadronicState";
447  break;
448  case 13:
449  ParticleStatusName = "kIStPreDecayResonantState";
450  break;
451  case 14:
452  ParticleStatusName = "kIStHadronInTheNucleus";
453  break;
454  case 15:
455  ParticleStatusName = "kIStFinalStateNuclearRemnant";
456  break;
457  case 16:
458  ParticleStatusName = "kIStNucleonClusterTarget";
459  break;
460  default:
461  ParticleStatusName = "Status Unknown";
462  }
463  return ParticleStatusName;
464  }
465 
466  //......................................................................
467  std::string GENIEGen::ReactionChannel(int ccnc,int mode)
468  {
469  std::string ReactionChannelName=" ";
470 
471  if(ccnc==0)
472  ReactionChannelName = "kCC";
473  else if(ccnc==1)
474  ReactionChannelName = "kNC";
475  else std::cout<<"Current mode unknown!! "<<std::endl;
476 
477  if(mode==0)
478  ReactionChannelName += "_kQE";
479  else if(mode==1)
480  ReactionChannelName += "_kRes";
481  else if(mode==2)
482  ReactionChannelName += "_kDIS";
483  else if(mode==3)
484  ReactionChannelName += "_kCoh";
485  else std::cout<<"interaction mode unknown!! "<<std::endl;
486 
487  return ReactionChannelName;
488  }
489 
490  //......................................................................
492  {
493  // Decide which histograms to put the spectrum in
494  int id = -1;
495  if (mc.GetNeutrino().CCNC()==simb::kCC) {
496  fCCMode->Fill(mc.GetNeutrino().Mode());
497  if (mc.GetNeutrino().Nu().PdgCode() == 12) id = 0;
498  else if (mc.GetNeutrino().Nu().PdgCode() == -12) id = 1;
499  else if (mc.GetNeutrino().Nu().PdgCode() == 14) id = 2;
500  else if (mc.GetNeutrino().Nu().PdgCode() == -14) id = 3;
501  else return;
502  }
503  else {
504  fNCMode->Fill(mc.GetNeutrino().Mode());
505  if (mc.GetNeutrino().Nu().PdgCode() > 0) id = 4;
506  else id = 5;
507  }
508  if (id==-1) abort();
509 
510  // Fill the specta histograms
511  fGenerated[id]->Fill(mc.GetNeutrino().Nu().E() );
512 
515  simb::MCNeutrino mcnu = mc.GetNeutrino();
516  const simb::MCParticle nu = mcnu.Nu();
517 
518  fVertexX->Fill(nu.Vx());
519  fVertexY->Fill(nu.Vy());
520  fVertexZ->Fill(nu.Vz());
521 
522  fVertexXY->Fill(nu.Vx(), nu.Vy());
523  fVertexXZ->Fill(nu.Vz(), nu.Vx());
524  fVertexYZ->Fill(nu.Vz(), nu.Vy());
525 
526  double mom = nu.P();
527  if(std::abs(mom) > 0.){
528  fDCosX->Fill(nu.Px()/mom);
529  fDCosY->Fill(nu.Py()/mom);
530  fDCosZ->Fill(nu.Pz()/mom);
531  }
532 
533 
534  LOG_DEBUG("GENIEInteractionInformation")
535  << std::endl
536  << "REACTION: " << ReactionChannel(mc.GetNeutrino().CCNC(),mc.GetNeutrino().Mode())
537  << std::endl
538  << "-----------> Particles in the Stack = " << mc.NParticles() << std::endl
539  << std::setiosflags(std::ios::left)
540  << std::setw(20) << "PARTICLE"
541  << std::setiosflags(std::ios::left)
542  << std::setw(32) << "STATUS"
543  << std::setw(18) << "E (GeV)"
544  << std::setw(18) << "m (GeV/c2)"
545  << std::setw(18) << "Ek (GeV)"
546  << std::endl << std::endl;
547 
548  const TDatabasePDG* databasePDG = TDatabasePDG::Instance();
549 
550  // Loop over the particle stack for this event
551  for(int i = 0; i < mc.NParticles(); ++i){
553  std::string name = databasePDG->GetParticle(part.PdgCode())->GetName();
554  int code = part.StatusCode();
555  std::string status = ParticleStatus(code);
556  double mass = part.Mass();
557  double energy = part.E();
558  double Ek = (energy-mass); // Kinetic Energy (GeV)
559  if(status=="kIStStableFinalState"||status=="kIStHadronInTheNucleus")
560  LOG_DEBUG("GENIEFinalState")
561  << std::setiosflags(std::ios::left) << std::setw(20) << name
562  << std::setiosflags(std::ios::left) << std::setw(32) <<status
563  << std::setw(18)<< energy
564  << std::setw(18)<< mass
565  << std::setw(18)<< Ek <<std::endl;
566  else
567  LOG_DEBUG("GENIEFinalState")
568  << std::setiosflags(std::ios::left) << std::setw(20) << name
569  << std::setiosflags(std::ios::left) << std::setw(32) << status
570  << std::setw(18) << energy
571  << std::setw(18) << mass <<std::endl;
572  }
573 
574 
575  if(mc.GetNeutrino().CCNC() == simb::kCC){
576 
579  for(int i = 0; i < mc.NParticles(); ++i){
581  if(abs(part.PdgCode()) == 11){
582  fEMomentum->Fill(part.P());
583  fEDCosX->Fill(part.Px()/part.P());
584  fEDCosY->Fill(part.Py()/part.P());
585  fEDCosZ->Fill(part.Pz()/part.P());
586  fECons->Fill(nu.E() - part.E());
587  break;
588  }
589  else if(abs(part.PdgCode()) == 13){
590  fMuMomentum->Fill(part.P());
591  fMuDCosX->Fill(part.Px()/part.P());
592  fMuDCosY->Fill(part.Py()/part.P());
593  fMuDCosZ->Fill(part.Pz()/part.P());
594  fECons->Fill(nu.E() - part.E());
595  break;
596  }
597  }// end loop over particles
598  }//end if CC interaction
599 
600  return;
601  }
602 
603 }
604 
605 namespace evgen{
606 
608 
609 }
610 
611 #endif
TH1F * fVertexX
vertex location of generated events in x
std::string ParticleStatus(int StatusCode)
Float_t x
Definition: compare.C:6
double E(const int i=0) const
Definition: MCParticle.h:237
TH1F * fMuDCosZ
direction cosine of outgoing mu in z
TH2F * fVertexXZ
vertex location in xz
int PdgCode() const
Definition: MCParticle.h:216
int CCNC() const
Definition: MCNeutrino.h:152
const simb::MCNeutrino & GetNeutrino() const
Definition: MCTruth.h:74
double Py(const int i=0) const
Definition: MCParticle.h:235
geo::Length_t DetHalfWidth(geo::TPCID const &tpcid) const
Returns the half width of the active volume of the specified TPC.
void DumpGTruth(Stream &&out, simb::GTruth const &truth, std::string indent, std::string firstIndent)
Dumps the content of the GENIE truth in the output stream.
Definition: MCDumpers.h:379
TH1F * fDeltaE
difference in neutrino energy from MCTruth::Enu() vs TParticle
TH2F * fVertexXY
vertex location in xy
void beginSubRun(art::SubRun &sr)
std::string ReactionChannel(int ccnc, int mode)
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
void beginRun(art::Run &run)
Unknown view.
Definition: geo_types.h:83
TH1F * fMuDCosY
direction cosine of outgoing mu in y
Float_t y
Definition: compare.C:6
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:150
bool Sample(simb::MCTruth &truth, simb::MCFlux &flux, simb::GTruth &gtruth)
Double_t z
Definition: plot.C:279
double Px(const int i=0) const
Definition: MCParticle.h:234
STL namespace.
offset to account for adding in Nuance codes to this enum
Definition: MCNeutrino.h:99
void endSubRun(art::SubRun &sr)
int NParticles() const
Definition: MCTruth.h:72
art::ProductID put(std::unique_ptr< PROD > &&)
Definition: Run.h:148
TGeoManager * ROOTGeoManager() const
Access to the ROOT geometry description manager.
TH1F * fNCMode
CC interaction mode.
object containing MC flux information
std::string ROOTFile() const
Returns the full directory path to the geometry file source.
Definition: Run.h:30
::sim::BeamType_t fBeamType
The width of a simulated "beam gate".
NuMI.
Definition: BeamTypes.h:12
double TotalMass() const
Returns the total mass [kg] of the specified volume (default: world).
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
TH1F * fDCosZ
direction cosine in z
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
int InteractionType() const
Definition: MCNeutrino.h:154
void produce(art::Event &evt)
geo::Length_t DetHalfHeight(geo::TPCID const &tpcid) const
Returns the half height of the active volume of the specified TPC.
TH1F * fEDCosY
direction cosine of outgoing e in y
Utility functions to print MC truth information.
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
GENIEGen(fhicl::ParameterSet const &pset)
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
long seed
Definition: chem4.cc:68
int fPassEmptySpills
whether or not to kill evnets with no interactions
TH1F * fCCMode
CC interaction mode.
double fPrevTotGoodPOT
Total good POT from subruns previous to current subrun.
evgb::GENIEHelper * fGENIEHelp
GENIEHelper object.
double P(const int i=0) const
Definition: MCParticle.h:238
T get(std::string const &key) const
Definition: ParameterSet.h:231
TString part[npart]
Definition: Style.C:32
double fGlobalTimeOffset
keep track of how long it takes to run the job
TH1F * fEDCosX
direction cosine of outgoing e in x
double energy
Definition: plottest35.C:25
TH1F * fMuDCosX
direction cosine of outgoing mu in x
void DumpMCTruth(Stream &&out, simb::MCTruth const &truth, unsigned int pointsPerLine, std::string indent, std::string firstIndent)
Dumps the content of the specified MC truth in the output stream.
Definition: MCDumpers.h:346
TH1F * fGenerated[6]
Spectra as generated.
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
bool get_if_present(std::string const &key, T &value) const
Definition: ParameterSet.h:208
bool CreateAssn(PRODUCER const &prod, art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t indx=UINT_MAX)
Creates a single one-to-one association.
Wrapper for generating neutrino interactions with GENIE.
double fPrevTotPOT
The type of beam.
Monte Carlo Simulation.
An art service to assist in the distribution of guaranteed unique seeds to all engines within an art ...
TH1F * fDCosY
direction cosine in y
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:73
double Vx(const int i=0) const
Definition: MCParticle.h:225
constexpr auto const & left(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:104
TH1F * fEMomentum
momentum of outgoing electrons
T * make(ARGS...args) const
BNB.
Definition: BeamTypes.h:11
TH2F * fVertexYZ
vertex location in yz
Utility object to perform functions of association.
double TotalExposure() const
Definition: GENIEHelper.h:66
TH1F * fECons
histogram to determine if energy is conserved in the event
ProductID put(std::unique_ptr< PROD > &&)
double Pz(const int i=0) const
Definition: MCParticle.h:236
double Vz(const int i=0) const
Definition: MCParticle.h:227
#define LOG_DEBUG(id)
TH1F * fEDCosZ
direction cosine of outgoing e in z
TH1F * fDCosX
direction cosine in x
bool NeutrinoSet() const
Definition: MCTruth.h:75
TCEvent evt
Definition: DataStructs.cxx:5
Event generator information.
Definition: MCTruth.h:30
TStopwatch fStopwatch
Event generator information.
Definition: MCNeutrino.h:18
Namespace collecting geometry-related classes utilities.
TH1F * fMuMomentum
momentum of outgoing muons
int Mode() const
Definition: MCNeutrino.h:153
Event Generation using GENIE, cosmics or single particles.
std::string ParticleStatusName(int code)
Describes the status of a particle (simb::MCParticle::StatusCode()).
TH1F * fVertexY
vertex location of generated events in y
BeamType_t
Defines category of beams to be stored in sim::BeamGateInfo.
Definition: BeamTypes.h:9
void FillHistograms(simb::MCTruth mc)
std::vector< double > fVtxPosHistRange
use defined hist range; it is useful to have for asymmetric ranges like in DP FD. ...
TH1F * fVertexZ
vertex location of generated events in z
double Vy(const int i=0) const
Definition: MCParticle.h:226
void put(std::string const &key)
art framework interface to geometry description
A module to check the results from the Monte Carlo generator.
double fRandomTimeOffset
The start of a simulated "beam gate".