LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
NDKGen_module.cc
Go to the documentation of this file.
1 //
3 //
4 // NDK neutrino event generator
5 //
6 // echurch@fnal.gov
7 //
9 #include <cstdlib>
10 #include <string>
11 #include <iostream>
12 #include <iomanip>
13 #include <sstream>
14 #include <vector>
15 #include <map>
16 #include <memory>
17 #include <unistd.h>
18 #include <stdio.h>
19 #include <fstream>
20 
21 // ROOT includes
22 #include "TH1.h"
23 #include "TH2.h"
24 #include "TDatabasePDG.h"
25 #include "TSystem.h"
26 
27 #include "CLHEP/Random/RandFlat.h"
28 
29 // Framework includes
32 #include "fhiclcpp/ParameterSet.h"
39 
40 // art extensions
42 
43 
44 // LArSoft includes
51 
53 
54 #include "TStopwatch.h"
55 
56 class TH1F;
57 class TH2F;
58 
59 namespace simb { class MCTruth; }
60 
61 
62 namespace evgen {
64  class NDKGen : public art::EDProducer {
65  public:
66  explicit NDKGen(fhicl::ParameterSet const &pset);
67  virtual ~NDKGen();
68 
69  void produce(art::Event& evt);
70  void beginJob();
71  void beginRun(art::Run& run);
72  void reconfigure(fhicl::ParameterSet const& p);
73  void endJob();
74 
75  private:
76 
77  std::string ParticleStatus(int StatusCode);
78  std::string ReactionChannel(int ccnc,int mode);
79 
80  void FillHistograms(simb::MCTruth mc);
81 
82  std::string fNdkFile;
83  std::ifstream *fEventFile;
84  TStopwatch fStopwatch;
85 
86  std::string fNDKModuleLabel;
87 
88  TH1F* fGenerated[6];
89 
90  TH1F* fVertexX;
91  TH1F* fVertexY;
92  TH1F* fVertexZ;
93 
94  TH2F* fVertexXY;
95  TH2F* fVertexXZ;
96  TH2F* fVertexYZ;
97 
98  TH1F* fDCosX;
99  TH1F* fDCosY;
100  TH1F* fDCosZ;
101 
102  TH1F* fMuMomentum;
103  TH1F* fMuDCosX;
104  TH1F* fMuDCosY;
105  TH1F* fMuDCosZ;
106 
107  TH1F* fEMomentum;
108  TH1F* fEDCosX;
109  TH1F* fEDCosY;
110  TH1F* fEDCosZ;
111 
112  TH1F* fCCMode;
113  TH1F* fNCMode;
114 
115  // for c2: fDeltaE is no longer used
116  //TH1F* fDeltaE; ///< difference in neutrino energy from MCTruth::Enu() vs TParticle
117  TH1F* fECons;
118 
119  };
120 } // namespace
121 
122 namespace evgen{
123 
124  //____________________________________________________________________________
125  NDKGen::NDKGen(fhicl::ParameterSet const& pset)
126  {
127  this->reconfigure(pset);
128  fStopwatch.Start();
129 
130  produces< std::vector<simb::MCTruth> >();
131  produces< sumdata::RunData, art::InRun >();
132 
133  fEventFile = new std::ifstream(fNdkFile.c_str());
134  if(!fEventFile->good())
135  exit(0);
136  // create a default random engine; obtain the random seed from NuRandomService,
137  // unless overridden in configuration with key "Seed"
139  ->createEngine(*this, pset, "Seed");
140 
141  }
142 
143  //____________________________________________________________________________
144  NDKGen::~NDKGen()
145  {
146  fStopwatch.Stop();
147  }
148 
149  //____________________________________________________________________________
150  void NDKGen::reconfigure(fhicl::ParameterSet const& p)
151  {
152  fNdkFile =(p.get< std::string >("NdkFile"));
153  return;
154  }
155 //___________________________________________________________________________
156 
158 
159  // Get access to the TFile service.
161 
162  fGenerated[0] = tfs->make<TH1F>("fGenerated_necc","", 100, 0.0, 20.0);
163  fGenerated[1] = tfs->make<TH1F>("fGenerated_nebcc","", 100, 0.0, 20.0);
164  fGenerated[2] = tfs->make<TH1F>("fGenerated_nmcc","", 100, 0.0, 20.0);
165  fGenerated[3] = tfs->make<TH1F>("fGenerated_nmbcc","", 100, 0.0, 20.0);
166  fGenerated[4] = tfs->make<TH1F>("fGenerated_nnc","", 100, 0.0, 20.0);
167  fGenerated[5] = tfs->make<TH1F>("fGenerated_nbnc","", 100, 0.0, 20.0);
168 
169  fDCosX = tfs->make<TH1F>("fDCosX", ";dx/ds", 200, -1., 1.);
170  fDCosY = tfs->make<TH1F>("fDCosY", ";dy/ds", 200, -1., 1.);
171  fDCosZ = tfs->make<TH1F>("fDCosZ", ";dz/ds", 200, -1., 1.);
172 
173  fMuMomentum = tfs->make<TH1F>("fMuMomentum", ";p_{#mu} (GeV/c)", 500, 0., 50.);
174  fMuDCosX = tfs->make<TH1F>("fMuDCosX", ";dx/ds;", 200, -1., 1.);
175  fMuDCosY = tfs->make<TH1F>("fMuDCosY", ";dy/ds;", 200, -1., 1.);
176  fMuDCosZ = tfs->make<TH1F>("fMuDCosZ", ";dz/ds;", 200, -1., 1.);
177 
178  fEMomentum = tfs->make<TH1F>("fEMomentum", ";p_{e} (GeV/c)", 500, 0., 50.);
179  fEDCosX = tfs->make<TH1F>("fEDCosX", ";dx/ds;", 200, -1., 1.);
180  fEDCosY = tfs->make<TH1F>("fEDCosY", ";dy/ds;", 200, -1., 1.);
181  fEDCosZ = tfs->make<TH1F>("fEDCosZ", ";dz/ds;", 200, -1., 1.);
182 
183  fCCMode = tfs->make<TH1F>("fCCMode", ";CC Interaction Mode;", 5, 0., 5.);
184  fCCMode->GetXaxis()->SetBinLabel(1, "QE");
185  fCCMode->GetXaxis()->SetBinLabel(2, "Res");
186  fCCMode->GetXaxis()->SetBinLabel(3, "DIS");
187  fCCMode->GetXaxis()->SetBinLabel(4, "Coh");
188  fCCMode->GetXaxis()->SetBinLabel(5, "kInverseMuDecay");
189  fCCMode->GetXaxis()->CenterLabels();
190 
191  fNCMode = tfs->make<TH1F>("fNCMode", ";NC Interaction Mode;", 5, 0., 5.);
192  fNCMode->GetXaxis()->SetBinLabel(1, "QE");
193  fNCMode->GetXaxis()->SetBinLabel(2, "Res");
194  fNCMode->GetXaxis()->SetBinLabel(3, "DIS");
195  fNCMode->GetXaxis()->SetBinLabel(4, "Coh");
196  fNCMode->GetXaxis()->SetBinLabel(5, "kNuElectronElastic");
197  fNCMode->GetXaxis()->CenterLabels();
198 
199  //fDeltaE = tfs->make<TH1F>("fDeltaE", ";#Delta E_{#nu} (GeV);", 200, -1., 1.);
200  fECons = tfs->make<TH1F>("fECons", ";#Delta E(#nu,lepton);", 500, -5., 5.);
201 
203  double x = 2.1*geo->DetHalfWidth();
204  double y = 2.1*geo->DetHalfHeight();
205  double z = 2.*geo->DetLength();
206  int xdiv = TMath::Nint(2*x/5.);
207  int ydiv = TMath::Nint(2*y/5.);
208  int zdiv = TMath::Nint(2*z/5.);
209 
210  fVertexX = tfs->make<TH1F>("fVertexX", ";x (cm)", xdiv, -x, x);
211  fVertexY = tfs->make<TH1F>("fVertexY", ";y (cm)", ydiv, -y, y);
212  fVertexZ = tfs->make<TH1F>("fVertexZ", ";z (cm)", zdiv, -0.2*z, z);
213 
214  fVertexXY = tfs->make<TH2F>("fVertexXY", ";x (cm);y (cm)", xdiv, -x, x, ydiv, -y, y);
215  fVertexXZ = tfs->make<TH2F>("fVertexXZ", ";z (cm);x (cm)", zdiv, -0.2*z, z, xdiv, -x, x);
216  fVertexYZ = tfs->make<TH2F>("fVertexYZ", ";z (cm);y (cm)", zdiv, -0.2*z, z, ydiv, -y, y);
217 
218  }
219 
220  //____________________________________________________________________________
221  void NDKGen::beginRun(art::Run& run)
222  {
223 
224  // grab the geometry object to see what geometry we are using
226  std::unique_ptr<sumdata::RunData> runcol(new sumdata::RunData(geo->DetectorName()));
227 
228  run.put(std::move(runcol));
229 
230  return;
231  }
232 
233  //____________________________________________________________________________
234  void NDKGen::endJob()
235  {
236  fEventFile->close();
237  }
238  //____________________________________________________________________________
239  void NDKGen::produce(art::Event& evt)
240  {
241 
242  std::cout << std::endl;
243  std::cout<<"------------------------------------------------------------------------------"<<std::endl;
244  //std::cout << "run : " << evt.Header().Run() << std::endl;
245  //std::cout << "subrun : " << evt.Header().Subrun() << std::endl;
246  //std::cout << "event : " << evt.Header().Event() << std::endl;
247  std::cout << "event : " << evt.id().event() << std::endl;
248 
249  // TODO: fill more quantities out, as below.
250  /*
251  double X; // vertex position from Ndk
252  double Y; // vertex position from Ndk
253  double Z; // vertex position from Ndk
254  double PDGCODE = -9999.;
255  double CHANNEL = -9999.;
256  int channel = -9999;
257  double energy = 0.; // in MeV from Ndk
258  double cosx = 0.;
259  double cosy = 0.;
260  double cosz = 0.;
261 
262  int partnumber = 0;
263  */
264 
265  std::string name, k, dollar;
266 
267 
268  // event dump format on file output by the two commands ....
269  // gevgen_ndcy -g 1000180400 -m 8 -n 400 -o ndk
270  // gevdump -f ndk.1000.ghep.root > ndk.out
271  std::string Name;
272  int Idx, Ist, PDG, Mother1, Mother2, Daughter1 ,Daughter2;
273  double Px, Py, Pz, E, m ;
274  std::string p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14;
275 
276 
277  int trackid = -1; // set track id to -i as these are all primary particles and have id <= 0
278  std::string primary("primary");
279  int FirstMother = -1;
280  double Mass = -9999;
281  int Status = -9999;
282 
283  double P; // momentum of MCParticle IN GeV/c
284 
285  // TODO: Could perhaps imagine using these in NDk.
286  /*
287  int targetnucleusPdg = -9999;
288  int hitquarkPdg = -9999;
289  double Q2 = -9999;
290  */
291  TLorentzVector Neutrino;
292  TLorentzVector Lepton;
293  TLorentzVector Target;
294  TLorentzVector q;
295  TLorentzVector Hadron4mom;
296 
297  // TODO: Could perhaps imagine using these in NDk.
298  /*
299  int Tpdg = 0; // for target
300  double Tmass = 0;
301  int Tstatus = 11;
302  double Tcosx, Tcosy, Tcosz, Tenergy;
303  */
304 
305  TLorentzVector Tpos;
306 
307 
308  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
309  simb::MCTruth truth;
310 
313  CLHEP::HepRandomEngine &engine = rng->getEngine();
314  CLHEP::RandFlat flat(engine);
315 
316  double fvCut (5.0); // force vtx to be this far from any wall.
317 
318  // Find boundary of active volume
319  double minx = 1e9;
320  double maxx = -1e9;
321  double miny = 1e9;
322  double maxy = -1e9;
323  double minz = 1e9;
324  double maxz = -1e9;
325  for (size_t i = 0; i<geo->NTPC(); ++i)
326  {
327  double local[3] = {0.,0.,0.};
328  double world[3] = {0.,0.,0.};
329  const geo::TPCGeo &tpc = geo->TPC(i);
330  tpc.LocalToWorld(local,world);
331  if (minx>world[0]-geo->DetHalfWidth(i))
332  minx = world[0]-geo->DetHalfWidth(i);
333  if (maxx<world[0]+geo->DetHalfWidth(i))
334  maxx = world[0]+geo->DetHalfWidth(i);
335  if (miny>world[1]-geo->DetHalfHeight(i))
336  miny = world[1]-geo->DetHalfHeight(i);
337  if (maxy<world[1]+geo->DetHalfHeight(i))
338  maxy = world[1]+geo->DetHalfHeight(i);
339  if (minz>world[2]-geo->DetLength(i)/2.)
340  minz = world[2]-geo->DetLength(i)/2.;
341  if (maxz<world[2]+geo->DetLength(i)/2.)
342  maxz = world[2]+geo->DetLength(i)/2.;
343  }
344 
345  // Assign vertice position
346  double X0 = 0.0 + flat.fire( minx+fvCut , maxx-fvCut );
347  double Y0 = 0.0 + flat.fire( miny+fvCut , maxy-fvCut );
348  double Z0 = 0.0 + flat.fire( minz+fvCut , maxz-fvCut );
349 
350  std::cout << "NDKGen_module: X, Y, Z of vtx: " << X0 << ", "<< Y0 << ", "<< Z0 << std::endl;
351 
352  int GenieEvt = -999;
353 
354  if(!fEventFile->good())
355  std::cout << "NdkFile: Problem reading Ndk file" << std::endl;
356 
357  while(getline(*fEventFile,k)){
358 
359  if (k.find("** Event:")!= std::string::npos) {
360  std::istringstream in;
361  in.clear();
362  in.str(k);
363  std::string dummy;
364  in>> dummy>> dummy>> dummy >> dummy>> dummy>> dummy>> dummy >> dummy>> dummy>> dummy >> GenieEvt;
365  std::cout<<"Genie Evt "<< GenieEvt <<" art evt "<<evt.id().event()<<"\n";
366  }
367 
368  if (GenieEvt+1 != static_cast<int>(evt.id().event()))
369  continue;
370  else {
371 
372  if (!k.compare(0,25,"GENIE Interaction Summary")) // testing for new event.
373  break;
374  if (k.compare(0,1,"|") || k.compare(1,2," ")) continue; // uninteresting line if it doesn't start with "|" and if second and third characters aren't spaces.
375  if (k.find("Fin-Init") != std::string::npos) continue; // Meh.
376  if (k.find("Ar") != std::string::npos) continue; // Meh.
377  if (k.find("Cl") != std::string::npos) continue; // ignore chlorine nucleus in nnbar events
378  if (k.find("HadrBlob") != std::string::npos) continue; // Meh.
379  if (k.find("NucBindE") != std::string::npos) continue; // Meh. atmo
380  if (k.find("FLAGS") != std::string::npos) break; // Event end. gevgen_ndcy
381  if (k.find("Vertex") != std::string::npos) break; // Event end. atmo
382 
383  // if (!k.compare(26,1,"3") || !k.compare(26,1,"1")) ; // New event or stable particles.
384  if (!k.compare(26,1,"1")) // New event or stable particles.
385  {
386 
387  std::istringstream in;
388  in.clear();
389  in.str(k);
390 
391  in>>p1>> Idx >>p2>> Name >>p3>> Ist >>p4>> PDG >>p5>>Mother1 >> p6 >> Mother2 >>p7>> Daughter1 >>p8>> Daughter2 >>p9>>Px>>p10>>Py>>p11>>Pz>>p12>>E>>p13>> m>>p14;
392  //std::cout<<std::setprecision(9)<<dollar<<" "<<name<<" "<<PDGCODE<<" "<<energy<<" "<<cosx<<" "<<cosy<<" "<<cosz<<" "<<partnumber<<std::endl;
393  if (Ist!=1) continue;
394 
395  std::cout << "PDG = " << PDG << std::endl;
396 
397  const TDatabasePDG* databasePDG = TDatabasePDG::Instance();
398  const TParticlePDG* definition = databasePDG->GetParticle(PDG);
399  Mass = definition->Mass(); // GeV
400  if (E-Mass < 0.001) continue; // KE is too low.
401 
402  // if(partnumber == -1)
403  Status = 1;
404 
405  simb::MCParticle mcpart(trackid,
406  PDG,
407  primary,
408  FirstMother,
409  Mass,
410  Status
411  );
412 
413  P = std::sqrt(pow(E,2.) - pow(Mass,2.)); // GeV/c
414  std::cout << "Momentum = " << P << std::endl;
415 
416  TLorentzVector pos(X0, Y0, Z0, 0);
417 
418  Tpos = pos; // for target
419 
420  TLorentzVector mom(Px,Py,Pz, E);
421 
422  mcpart.AddTrajectoryPoint(pos,mom);
423  truth.Add(mcpart);
424 
425 
426  }// loop over particles in an event
427  truth.SetOrigin(simb::kUnknown);
428 
429  //if (!k.compare(1,1,"FLAGS")) // end of event
430  // break;
431 
432  }
433  } // end while loop
434 
436  std::cout << "NDKGen.cxx: Putting " << truth.NParticles() << " tracks on stack." << std::endl;
437  truthcol->push_back(truth);
438  //FillHistograms(truth);
439  evt.put(std::move(truthcol));
440 
441  return;
442  }
443 
444 // //......................................................................
445  std::string NDKGen::ParticleStatus(int StatusCode)
446  {
447  int code = StatusCode;
448  std::string ParticleStatusName;
449 
450  switch(code)
451  {
452  case 0:
453  ParticleStatusName = "kIStInitialState";
454  break;
455  case 1:
456  ParticleStatusName = "kIStFinalState";
457  break;
458  case 11:
459  ParticleStatusName = "kIStNucleonTarget";
460  break;
461  default:
462  ParticleStatusName = "Status Unknown";
463  }
464  return ParticleStatusName;
465  }
466 
467 
468 // //......................................................................
469  std::string NDKGen::ReactionChannel(int ccnc,int mode)
470  {
471  std::string ReactionChannelName=" ";
472 
473  if(ccnc==0)
474  ReactionChannelName = "kCC";
475  else if(ccnc==1)
476  ReactionChannelName = "kNC";
477  else std::cout<<"Current mode unknown!! "<<std::endl;
478 
479  if(mode==0)
480  ReactionChannelName += "_kQE";
481  else if(mode==1)
482  ReactionChannelName += "_kRes";
483  else if(mode==2)
484  ReactionChannelName += "_kDIS";
485  else if(mode==3)
486  ReactionChannelName += "_kCoh";
487  else if(mode==4)
488  ReactionChannelName += "_kNuElectronElastic";
489  else if(mode==5)
490  ReactionChannelName += "_kInverseMuDecay";
491  else std::cout<<"interaction mode unknown!! "<<std::endl;
492 
493  return ReactionChannelName;
494  }
495 
496 // //......................................................................
497  void NDKGen::FillHistograms(simb::MCTruth mc)
498  {
499  // Decide which histograms to put the spectrum in
500  int id = -1;
501  if (mc.GetNeutrino().CCNC()==simb::kCC) {
502  fCCMode->Fill(mc.GetNeutrino().Mode());
503  if (mc.GetNeutrino().Nu().PdgCode() == 12) id = 0;
504  else if (mc.GetNeutrino().Nu().PdgCode() == -12) id = 1;
505  else if (mc.GetNeutrino().Nu().PdgCode() == 14) id = 2;
506  else if (mc.GetNeutrino().Nu().PdgCode() == -14) id = 3;
507  else return;
508  }
509  else {
510  fNCMode->Fill(mc.GetNeutrino().Mode());
511  if (mc.GetNeutrino().Nu().PdgCode() > 0) id = 4;
512  else id = 5;
513  }
514  if (id==-1) abort();
515 
516  // Fill the specta histograms
517  fGenerated[id]->Fill(mc.GetNeutrino().Nu().E() );
518 
519  //< fill the vertex histograms from the neutrino - that is always
520  //< particle 0 in the list
521  simb::MCNeutrino mcnu = mc.GetNeutrino();
522  const simb::MCParticle nu = mcnu.Nu();
523 
524  fVertexX->Fill(nu.Vx());
525  fVertexY->Fill(nu.Vy());
526  fVertexZ->Fill(nu.Vz());
527 
528  fVertexXY->Fill(nu.Vx(), nu.Vy());
529  fVertexXZ->Fill(nu.Vz(), nu.Vx());
530  fVertexYZ->Fill(nu.Vz(), nu.Vy());
531 
532  double mom = nu.P();
533  if(std::abs(mom) > 0.){
534  fDCosX->Fill(nu.Px()/mom);
535  fDCosY->Fill(nu.Py()/mom);
536  fDCosZ->Fill(nu.Pz()/mom);
537  }
538 
539 
540 // LOG_DEBUG("GENIEInteractionInformation")
541 // << std::endl
542 // << "REACTION: " << ReactionChannel(mc.GetNeutrino().CCNC(),mc.GetNeutrino().Mode())
543 // << std::endl
544 // << "-----------> Particles in the Stack = " << mc.NParticles() << std::endl
545 // << std::setiosflags(std::ios::left)
546 // << std::setw(20) << "PARTICLE"
547 // << std::setiosflags(std::ios::left)
548 // << std::setw(32) << "STATUS"
549 // << std::setw(18) << "E (GeV)"
550 // << std::setw(18) << "m (GeV/c2)"
551 // << std::setw(18) << "Ek (GeV)"
552 // << std::endl << std::endl;
553 
554 // const TDatabasePDG* databasePDG = TDatabasePDG::Instance();
555 
556 // // Loop over the particle stack for this event
557 // for(int i = 0; i < mc.NParticles(); ++i){
558 // simb::MCParticle part(mc.GetParticle(i));
559 // std::string name = databasePDG->GetParticle(part.PdgCode())->GetName();
560 // int code = part.StatusCode();
561 // std::string status = ParticleStatus(code);
562 // double mass = part.Mass();
563 // double energy = part.E();
564 // double Ek = (energy-mass); // Kinetic Energy (GeV)
565 // if(status=="kIStFinalStB4Interactions")
566 // LOG_DEBUG("GENIEFinalState")
567 // << std::setiosflags(std::ios::left) << std::setw(20) << name
568 // << std::setiosflags(std::ios::left) << std::setw(32) <<status
569 // << std::setw(18)<< energy
570 // << std::setw(18)<< mass
571 // << std::setw(18)<< Ek <<std::endl;
572 // else
573 // LOG_DEBUG("GENIEFinalState")
574 // << std::setiosflags(std::ios::left) << std::setw(20) << name
575 // << std::setiosflags(std::ios::left) << std::setw(32) << status
576 // << std::setw(18) << energy
577 // << std::setw(18) << mass <<std::endl;
578 
579  std::cout << "REACTION: " << ReactionChannel(mc.GetNeutrino().CCNC(),mc.GetNeutrino().Mode()) << std::endl;
580  std::cout << "-----------> Particles in the Stack = " << mc.NParticles() << std::endl;
581  std::cout << std::setiosflags(std::ios::left)
582  << std::setw(20) << "PARTICLE"
583  << std::setiosflags(std::ios::left)
584  << std::setw(32) << "STATUS"
585  << std::setw(18) << "E (GeV)"
586  << std::setw(18) << "m (GeV/c2)"
587  << std::setw(18) << "Ek (GeV)"
588  << std::endl << std::endl;
589 
590  const TDatabasePDG* databasePDG = TDatabasePDG::Instance();
591 
592  // Loop over the particle stack for this event
593  for(int i = 0; i < mc.NParticles(); ++i){
595  std::string name;
596  if (part.PdgCode() == 18040)
597  name = "Ar40 18040";
598  else if (part.PdgCode() != -99999 )
599  {
600  name = databasePDG->GetParticle(part.PdgCode())->GetName();
601  }
602 
603  int code = part.StatusCode();
604  std::string status = ParticleStatus(code);
605  double mass = part.Mass();
606  double energy = part.E();
607  double Ek = (energy-mass); // Kinetic Energy (GeV)
608 
609  std::cout << std::setiosflags(std::ios::left) << std::setw(20) << name
610  << std::setiosflags(std::ios::left) << std::setw(32) <<status
611  << std::setw(18)<< energy
612  << std::setw(18)<< mass
613  << std::setw(18)<< Ek <<std::endl;
614  }
615 
616  if(mc.GetNeutrino().CCNC() == simb::kCC){
617 
620  for(int i = 0; i < mc.NParticles(); ++i){
622  if(abs(part.PdgCode()) == 11){
623  fEMomentum->Fill(part.P());
624  fEDCosX->Fill(part.Px()/part.P());
625  fEDCosY->Fill(part.Py()/part.P());
626  fEDCosZ->Fill(part.Pz()/part.P());
627  fECons->Fill(nu.E() - part.E());
628  break;
629  }
630  else if(abs(part.PdgCode()) == 13){
631  fMuMomentum->Fill(part.P());
632  fMuDCosX->Fill(part.Px()/part.P());
633  fMuDCosY->Fill(part.Py()/part.P());
634  fMuDCosZ->Fill(part.Pz()/part.P());
635  fECons->Fill(nu.E() - part.E());
636  break;
637  }
638  }// end loop over particles
639  }//end if CC interaction
640 
641  return;
642  }
643 
645 
646 } // namespace
Float_t x
Definition: compare.C:6
double E(const int i=0) const
Definition: MCParticle.h:237
TH1F * fMuMomentum
momentum of outgoing muons
TH2F * fVertexXY
vertex location in xy
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
TH1F * fMuDCosY
direction cosine of outgoing mu in y
TH1F * fMuDCosZ
direction cosine of outgoing mu in z
geo::Length_t DetHalfWidth(geo::TPCID const &tpcid) const
Returns the half width of the active volume of the specified TPC.
TH1F * fECons
histogram to determine if energy is conserved in the event
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:261
Float_t E
Definition: plot.C:23
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:78
Float_t y
Definition: compare.C:6
TH1F * fEDCosZ
direction cosine of outgoing e in z
std::string fNdkFile
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:150
TH1F * fEMomentum
momentum of outgoing electrons
Double_t z
Definition: plot.C:279
TH2F * fVertexXZ
vertex location in xz
Geometry information for a single TPC.
Definition: TPCGeo.h:37
double Px(const int i=0) const
Definition: MCParticle.h:234
TH1F * fVertexY
vertex location of generated events in y
int NParticles() const
Definition: MCTruth.h:72
art::ProductID put(std::unique_ptr< PROD > &&)
Definition: Run.h:148
Particle class.
void Add(simb::MCParticle &part)
Definition: MCTruth.h:77
TH1F * fCCMode
CC interaction mode.
TH1F * fEDCosY
direction cosine of outgoing e in y
Definition: Run.h:30
TStopwatch fStopwatch
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
base_engine_t & getEngine() const
A module to check the results from the Monte Carlo generator.
geo::Length_t DetHalfHeight(geo::TPCID const &tpcid) const
Returns the half height of the active volume of the specified TPC.
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
void beginJob()
Definition: Breakpoints.cc:14
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
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
TH1F * fMuDCosX
direction cosine of outgoing mu in x
std::ifstream * fEventFile
double energy
Definition: plottest35.C:25
TH2F * fVertexYZ
vertex location in yz
TH1F * fVertexZ
vertex location of generated events in z
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
TH1F * fDCosZ
direction cosine in z
std::string fNDKModuleLabel
keep track of how long it takes to run the job
An art service to assist in the distribution of guaranteed unique seeds to all engines within an art ...
Framework includes.
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:73
TH1F * fEDCosX
direction cosine of outgoing e in x
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
T * make(ARGS...args) const
ifstream in
Definition: comparison.C:7
double Pz(const int i=0) const
Definition: MCParticle.h:236
double Vz(const int i=0) const
Definition: MCParticle.h:227
EventNumber_t event() const
Definition: EventID.h:117
TPCGeo const & TPC(unsigned int const tpc=0, unsigned int const cstat=0) const
Returns the specified TPC.
TCEvent evt
Definition: DataStructs.cxx:5
Event generator information.
Definition: MCTruth.h:30
TH1F * fNCMode
CC interaction mode.
TH1F * fDCosY
direction cosine in y
Event generator information.
Definition: MCNeutrino.h:18
Namespace collecting geometry-related classes utilities.
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 * fVertexX
vertex location of generated events in x
void LocalToWorld(const double *tpc, double *world) const
Transform point from local TPC frame to world frame.
Definition: TPCGeo.h:490
TH1F * fDCosX
direction cosine in x
EventID id() const
Definition: Event.h:56
double Vy(const int i=0) const
Definition: MCParticle.h:226
art framework interface to geometry description