NUANCEGen_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 //
4 // NUANCE neutrino event generator
5 //
6 // brebel@fnal.gov
7 // saima@ksu.edu
8 //
9 ////////////////////////////////////////////////////////////////////////
10 
11 #include <cstdlib>
12 #include <string>
13 #include <iostream>
14 #include <iomanip>
15 #include <sstream>
16 #include <memory>
17 #include <unistd.h>
18 #include <stdio.h>
19 #include <fstream>
20 #include <memory> // std::unique_ptr<>
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"
36 #include "art_root_io/TFileService.h"
37 #include "art_root_io/TFileDirectory.h"
40 
41 // LArSoft includes
46 
47 // class TH1F;
48 // class TH2F;
49 
50 // namespace simb { class MCTruth; }
51 
52 namespace evgen {
53  /// A module to check the results from the Monte Carlo generator
54  class NUANCEGen : public art::EDProducer {
55  public:
56  explicit NUANCEGen(fhicl::ParameterSet const &pset);
57  virtual ~NUANCEGen();
58 
59  void produce(art::Event& evt);
60  void beginJob();
61  void beginRun(art::Run& run);
62  void reconfigure(fhicl::ParameterSet const& p);
63  void endJob();
64 
65  private:
66 
67  std::string ParticleStatus(int StatusCode);
68  std::string ReactionChannel(int ccnc,int mode);
69 
71 
74  std::unique_ptr<std::ifstream> fEventFile;
75  TStopwatch fStopwatch; ///keep track of how long it takes to run the job
76 
77 
79 
80  TH1F* fGenerated[6]; ///< Spectra as generated
81 
82  TH1F* fVertexX; ///< vertex location of generated events in x
83  TH1F* fVertexY; ///< vertex location of generated events in y
84  TH1F* fVertexZ; ///< vertex location of generated events in z
85 
86  TH2F* fVertexXY; ///< vertex location in xy
87  TH2F* fVertexXZ; ///< vertex location in xz
88  TH2F* fVertexYZ; ///< vertex location in yz
89 
90  TH1F* fDCosX; ///< direction cosine in x
91  TH1F* fDCosY; ///< direction cosine in y
92  TH1F* fDCosZ; ///< direction cosine in z
93 
94  TH1F* fMuMomentum; ///< momentum of outgoing muons
95  TH1F* fMuDCosX; ///< direction cosine of outgoing mu in x
96  TH1F* fMuDCosY; ///< direction cosine of outgoing mu in y
97  TH1F* fMuDCosZ; ///< direction cosine of outgoing mu in z
98 
99  TH1F* fEMomentum; ///< momentum of outgoing electrons
100  TH1F* fEDCosX; ///< direction cosine of outgoing e in x
101  TH1F* fEDCosY; ///< direction cosine of outgoing e in y
102  TH1F* fEDCosZ; ///< direction cosine of outgoing e in z
103 
104  TH1F* fCCMode; ///< CC interaction mode
105  TH1F* fNCMode; ///< CC interaction mode
106 
107  TH1F* fECons; ///< histogram to determine if energy is conserved in the event
108 
109  };
110 }
111 
112 namespace evgen{
113 
114  //____________________________________________________________________________
116  : EDProducer{pset}
117  , fNuanceFile (pset.get< std::string >("NuanceFile"))
118  , fBeamVerticalAngle(pset.get< double >("BeamVerticalAngle"))
119  {
120  fStopwatch.Start();
121 
122  produces< std::vector<simb::MCTruth> >();
123  produces< sumdata::RunData, art::InRun >();
124 
125  mf::LogInfo("NUANCEGen") << "Reading events from '" << fNuanceFile << "'";
126  fEventFile = std::make_unique<std::ifstream>(fNuanceFile.c_str());
127  if (!(fEventFile->good())) {
129  << "Failed to open input file '" << fNuanceFile << "'.\n";
130  }
131  }
132 
133  //____________________________________________________________________________
135  {
136  fStopwatch.Stop();
137  }
138 
139 //___________________________________________________________________________
141  {
142  // Get access to the TFile service.
144 
145  fGenerated[0] = tfs->make<TH1F>("fGenerated_necc","", 100, 0.0, 20.0);
146  fGenerated[1] = tfs->make<TH1F>("fGenerated_nebcc","", 100, 0.0, 20.0);
147  fGenerated[2] = tfs->make<TH1F>("fGenerated_nmcc","", 100, 0.0, 20.0);
148  fGenerated[3] = tfs->make<TH1F>("fGenerated_nmbcc","", 100, 0.0, 20.0);
149  fGenerated[4] = tfs->make<TH1F>("fGenerated_nnc","", 100, 0.0, 20.0);
150  fGenerated[5] = tfs->make<TH1F>("fGenerated_nbnc","", 100, 0.0, 20.0);
151 
152  fDCosX = tfs->make<TH1F>("fDCosX", ";dx/ds", 200, -1., 1.);
153  fDCosY = tfs->make<TH1F>("fDCosY", ";dy/ds", 200, -1., 1.);
154  fDCosZ = tfs->make<TH1F>("fDCosZ", ";dz/ds", 200, -1., 1.);
155 
156  fMuMomentum = tfs->make<TH1F>("fMuMomentum", ";p_{#mu} (GeV/c)", 500, 0., 50.);
157  fMuDCosX = tfs->make<TH1F>("fMuDCosX", ";dx/ds;", 200, -1., 1.);
158  fMuDCosY = tfs->make<TH1F>("fMuDCosY", ";dy/ds;", 200, -1., 1.);
159  fMuDCosZ = tfs->make<TH1F>("fMuDCosZ", ";dz/ds;", 200, -1., 1.);
160 
161  fEMomentum = tfs->make<TH1F>("fEMomentum", ";p_{e} (GeV/c)", 500, 0., 50.);
162  fEDCosX = tfs->make<TH1F>("fEDCosX", ";dx/ds;", 200, -1., 1.);
163  fEDCosY = tfs->make<TH1F>("fEDCosY", ";dy/ds;", 200, -1., 1.);
164  fEDCosZ = tfs->make<TH1F>("fEDCosZ", ";dz/ds;", 200, -1., 1.);
165 
166  fCCMode = tfs->make<TH1F>("fCCMode", ";CC Interaction Mode;", 5, 0., 5.);
167  fCCMode->GetXaxis()->SetBinLabel(1, "QE");
168  fCCMode->GetXaxis()->SetBinLabel(2, "Res");
169  fCCMode->GetXaxis()->SetBinLabel(3, "DIS");
170  fCCMode->GetXaxis()->SetBinLabel(4, "Coh");
171  fCCMode->GetXaxis()->SetBinLabel(5, "kInverseMuDecay");
172  fCCMode->GetXaxis()->CenterLabels();
173 
174  fNCMode = tfs->make<TH1F>("fNCMode", ";NC Interaction Mode;", 5, 0., 5.);
175  fNCMode->GetXaxis()->SetBinLabel(1, "QE");
176  fNCMode->GetXaxis()->SetBinLabel(2, "Res");
177  fNCMode->GetXaxis()->SetBinLabel(3, "DIS");
178  fNCMode->GetXaxis()->SetBinLabel(4, "Coh");
179  fNCMode->GetXaxis()->SetBinLabel(5, "kNuElectronElastic");
180  fNCMode->GetXaxis()->CenterLabels();
181 
182  //fDeltaE = tfs->make<TH1F>("fDeltaE", ";#Delta E_{#nu} (GeV);", 200, -1., 1.);
183  fECons = tfs->make<TH1F>("fECons", ";#Delta E(#nu,lepton);", 500, -5., 5.);
184 
186  double x = 2.1*geo->DetHalfWidth();
187  double y = 2.1*geo->DetHalfHeight();
188  double z = 2.*geo->DetLength();
189  int xdiv = TMath::Nint(2*x/5.);
190  int ydiv = TMath::Nint(2*y/5.);
191  int zdiv = TMath::Nint(2*z/5.);
192 
193  fVertexX = tfs->make<TH1F>("fVertexX", ";x (cm)", xdiv, -x, x);
194  fVertexY = tfs->make<TH1F>("fVertexY", ";y (cm)", ydiv, -y, y);
195  fVertexZ = tfs->make<TH1F>("fVertexZ", ";z (cm)", zdiv, -0.2*z, z);
196 
197  fVertexXY = tfs->make<TH2F>("fVertexXY", ";x (cm);y (cm)", xdiv, -x, x, ydiv, -y, y);
198  fVertexXZ = tfs->make<TH2F>("fVertexXZ", ";z (cm);x (cm)", zdiv, -0.2*z, z, xdiv, -x, x);
199  fVertexYZ = tfs->make<TH2F>("fVertexYZ", ";z (cm);y (cm)", zdiv, -0.2*z, z, ydiv, -y, y);
200 
201  }
202 
203  //____________________________________________________________________________
205  {
207  run.put(std::make_unique<sumdata::RunData>(geo->DetectorName()));
208  }
209 
210  //____________________________________________________________________________
212  {
213  fEventFile->close();
214  }
215  //____________________________________________________________________________
217  {
218 
219  std::cout << std::endl;
220  std::cout<<"------------------------------------------------------------------------------"<<std::endl;
221 // std::cout << "run : " << evt.Header().Run() << std::endl;
222 // std::cout << "subrun : " << evt.Header().Subrun() << std::endl;
223 // std::cout << "event : " << evt.Header().Event() << std::endl;
224  std::cout << "event : " << evt.id().event() << std::endl;
225 
226  double X = 0.; // vertex position from Nuance
227  double Y = 0.; // vertex position from Nuance
228  double Z = 0.; // vertex position from Nuance
229  double PDGCODE = -9999.;
230  double CHANNEL = -9999.;
231  int channel = -9999;
232  double energy = 0.; // in MeV from Nuance
233  double cosx = 0.;
234  double cosy = 0.;
235  double cosz = 0.;
236  std::string name, k, dollar;
237  int partnumber = 0;
238 
239  int trackid = -1; // set track id to -i as these are all primary particles and have id <= 0
240  std::string primary("primary");
241  int FirstMother = -1;
242  double Mass = -9999;
243  int Status = -9999;
244 
245 
246  int ccnc = -9999;
247  int mode = -9999;
248  int targetnucleusPdg = -9999;
249  int hitquarkPdg = -9999;
250 
251  double P; // momentum of MCParticle in GeV/c
252  TLorentzVector Neutrino;
253  TLorentzVector Lepton;
254  TLorentzVector Target;
255  TLorentzVector q;
256  TLorentzVector Hadron4mom;
257 
258  double InvariantMass = -9999;
259  double x = -9999;
260  double y = -9999;
261  double Q2 = -9999;
262 
263  int Tpdg = 0; // for target
264  double Tmass = 0;
265  int Tstatus = 11;
266  double Tcosx = 0., Tcosy = 0., Tcosz = 0., Tenergy = 0.;
267  TLorentzVector Tpos;
268  double M = 0;
269 
270  // rotated direction cosines
271  double CosX = 0;
272  double CosY = 0;
273  double CosZ = 0;
274 
275 
276  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
277  simb::MCTruth truth;
278 
279  if(!fEventFile->good())
280  std::cout << "NuanceFile: Problem reading nuance file" << std::endl;
281 
282  while(getline(*fEventFile,k)){
283  std::istringstream in;
284  in.clear();
285  in.str(k);
286 
287  in>>dollar>>name>>PDGCODE>>energy>>cosx>>cosy>>cosz>>partnumber;
288  //std::cout<<std::setprecision(9)<<dollar<<" "<<name<<" "<<PDGCODE<<" "<<energy<<" "<<cosx<<" "<<cosy<<" "<<cosz<<" "<<partnumber<<std::endl;
289 
290  //get the nuance channel number
291  if(name == "nuance"){
292  CHANNEL = PDGCODE;
293  channel = int (CHANNEL);
294  channel = abs(channel)+ simb::kNuanceOffset;
295  //std::cout << "channel = " << channel << std::endl;
296 
297  //set the interaction type; CC or NC
298  if ( abs(channel) == simb::kCCQE
299  || ( abs(channel) >= simb::kResCCNuProtonPiPlus && abs(channel) <= simb::kResCCNuNeutronPiPlus )
301  || ( abs(channel) >= simb::kResCCNuDeltaPlusPiPlus && abs(channel) <= simb::kResCCNuDelta2PlusPiMinus )
303  || ( abs(channel) >= simb::kResCCNuProtonRhoPlus && abs(channel) <= simb::kResCCNuNeutronRhoPlus )
304  || ( abs(channel) >= simb::kResCCNuBarNeutronRhoMinus && abs(channel) <= simb::kResCCNuBarNeutronRho0 )
305  || ( abs(channel) >= simb::kResCCNuSigmaPlusKaonPlus && abs(channel) <= simb::kResCCNuSigmaPlusKaon0 )
306  || ( abs(channel) >= simb::kResCCNuBarSigmaMinusKaon0 && abs(channel) <= simb::kResCCNuBarSigma0Kaon0 )
307  || abs(channel) == simb::kResCCNuProtonEta
308  || abs(channel) == simb::kResCCNuBarNeutronEta
309  || abs(channel) == simb::kResCCNuKaonPlusLambda0
310  || abs(channel) == simb::kResCCNuBarKaon0Lambda0
311  || ( abs(channel) >= simb::kResCCNuProtonPiPlusPiMinus && abs(channel) <= simb::kResCCNuProtonPi0Pi0 )
313  || abs(channel) == simb::kCCDIS || abs(channel) == simb::kCCQEHyperon
314  || abs(channel) == simb::kCCCOH || abs(channel) == simb::kInverseMuDecay )
315  ccnc = simb::kCC;
316 
317  else if ( abs(channel) != simb::kUnUsed1 && abs(channel) != simb::kUnUsed2 )
318  ccnc = simb::kNC;
319 
320  //set the interaction mode; QE, Res, DIS, Coh, kNuElectronElastic, kInverseMuDecay
321  if ( abs(channel) == simb::kCCQE || abs(channel) == simb::kNCQE || abs(channel) == simb::kCCQEHyperon )
322  mode = simb::kQE;
323  else if ( abs(channel) >= simb::kResCCNuProtonPiPlus && abs(channel) <= simb::kResCCNuBarProtonPi0Pi0 )
324  mode = simb::kRes;
325  else if ( abs(channel) == simb::kCCDIS || abs(channel) == simb::kNCDIS )
326  mode = simb::kDIS;
327  else if ( abs(channel) == simb::kNCCOH || abs(channel) == simb::kCCCOH )
328  mode = simb::kCoh;
329  else if ( abs(channel) == simb::kNuElectronElastic )
330  mode = 4;
331  else if ( abs(channel) == simb::kInverseMuDecay )
332  mode = 5;
333 
334  } // if(name == "nuance")
335 
336 
337  //get the nuance vertex position:
338  if(name == "vertex"){
339  X = PDGCODE;
340  Y = energy;
341  Z = cosx;
342  //std::cout << "vertex from nuance = " << X << " " << Y << " " << Z << std::endl;
343  }
344 
345  double PI = 3.14159265;
346 
347  CosX = cosx;
348  CosY = (TMath::Cos(fBeamVerticalAngle*PI/180)*cosy)+(TMath::Sin(fBeamVerticalAngle*PI/180)*cosz);
349  CosZ = (-TMath::Sin(fBeamVerticalAngle*PI/180)*cosy)+(TMath::Cos(fBeamVerticalAngle*PI/180)*cosz);
350 
351  //get the target info
352  if(name == "track" && (PDGCODE == 2212 || PDGCODE == 2112 || PDGCODE == 18040 || PDGCODE == 11) && partnumber == -1){
353  Tpdg = int (PDGCODE);
354 
355  if ( PDGCODE == 18040 ){
356  Tmass = 37.5593438; // GeV
357  }
358 
359  else {
360  const TDatabasePDG* databasePDG = TDatabasePDG::Instance();
361  const TParticlePDG* definition = databasePDG->GetParticle(Tpdg);
362  Tmass = definition->Mass();
363  }
364 
365  Tenergy = energy;
366  Tcosx = CosX;
367  Tcosy = CosY;
368  Tcosz = CosZ;
369  }
370 
371  //get mcparticles other than target
372  if(name == "track" && (
373  partnumber == 0 ||
374  (partnumber == -1 && (PDGCODE != 2212 && PDGCODE != 2112 && PDGCODE != 18040 && PDGCODE != 11))
375  )
376  ){
377 
378  int pdgcode = int (PDGCODE);
379  //std::cout << "pdgcode = " << pdgcode << std::endl;
380 
381  if ( PDGCODE == 18040 )
382  Mass = 37.5593438; // GeV
383 
384  else {
385  const TDatabasePDG* databasePDG = TDatabasePDG::Instance();
386  const TParticlePDG* definition = databasePDG->GetParticle(pdgcode);
387  Mass = definition->Mass(); // GeV
388  }
389 
390  if(partnumber == -1)
391  Status = 0;
392 
393  if(partnumber == 0)
394  Status = 1;
395 
396  simb::MCParticle mcpart(trackid,
397  pdgcode,
398  primary,
399  FirstMother,
400  Mass,
401  Status
402  );
403 
404  P = std::sqrt(pow(energy/1000,2.) - pow(Mass,2.)); // GeV/c
405  //std::cout << "Momentum = " << P << std::endl;
406 
408 
409  double X0 = X + geo->DetHalfWidth();
410  double Y0 = Y;
411  double Z0 = Z + 0.5*geo->DetLength();
412 
413  TLorentzVector pos(X0, Y0, Z0, 0);
414  Tpos = pos; // for target
415 
416  TLorentzVector mom(CosX*P, CosY*P, CosZ*P, energy/1000);
417 
418  mcpart.AddTrajectoryPoint(pos,mom);
419  truth.Add(mcpart);
420 
421 
422  if(name == "track" && (abs(pdgcode) == 14 || abs(pdgcode) == 12) && partnumber == -1)
423  Neutrino.SetPxPyPzE(CosX*P, CosY*P, CosZ*P, energy/1000);
424 
425  if(name == "track" && (abs(pdgcode) == 13 || abs(pdgcode) == 11 || abs(pdgcode) == 14 || abs(pdgcode) == 12) && partnumber == 0)
426  Lepton.SetPxPyPzE(CosX*P, CosY*P, CosZ*P, energy/1000);
427 
428  }// loop over particles in an event. (excluding target)
429 
430  if(name == "end"){
431  break;
432  }
433 
434  } // end while loop
435 
436 
437  /////////////////////////////////
438 
439  //adding the target to mcparticle
440  simb::MCParticle mcpart(trackid,
441  Tpdg,
442  primary,
443  FirstMother,
444  Tmass,
445  Tstatus
446  );
447 
448  TLorentzVector Tmom;
449  Tmom.SetPxPyPzE(Tcosx, Tcosy, Tcosz, Tenergy/1000); // target momentum coordinates make a unit vector, may be modified later, does not effect the event simulation
450 
451  mcpart.AddTrajectoryPoint(Tpos,Tmom);
452  truth.Add(mcpart);
453 
454 // std::cout<<"Neutrino 4-Momentum = "
455 // <<Neutrino.Px()<<" "
456 // << Neutrino.Py()<<" "
457 // <<Neutrino.Pz()<<" "
458 // <<Neutrino.E()<<std::endl;
459 
460 // std::cout << "Target 4-Momentum = "
461 // << Target.Px() << " "
462 // << Target.Py() << " "
463 // << Target.Pz()<< " "
464 // << Target.E() << std::endl;
465 
466 // std::cout << "4-momentum of Lepton = "
467 // << Lepton.Px() << " "
468 // << Lepton.Py() << " "
469 // << Lepton.Pz() << " "
470 // << Lepton.E() << std::endl;
471 
472  M = Tmass;
473  q = Neutrino - Lepton;
474  Q2 = -(q*q);
475 
476  double v = q.E();
477  x = Q2/ (2*M*v);
478  y = v/ Neutrino.E();
479  double W2 = M*M + 2*M*v - Q2;
480  InvariantMass = TMath::Sqrt(TMath::Max(0.,W2));
481 
482  if(mode == simb::kCoh){
483  x = -1;
484  y = -1;
485  InvariantMass = -1;
486  }
487 
489  truth.SetNeutrino(ccnc, mode, channel,
490  targetnucleusPdg,
491  Tpdg,
492  hitquarkPdg,
493  //InvariantMass, x, y, Q2
494  InvariantMass, x, y, Q2
495  );
496 
497  std::cout << truth.GetNeutrino() << std::endl;
498 
499  truthcol->push_back(truth);
500  FillHistograms(truth);
501  evt.put(std::move(truthcol));
502 
503  return;
504  }
505 
506 // //......................................................................
508  {
509  int code = StatusCode;
511 
512  switch(code)
513  {
514  case 0:
515  ParticleStatusName = "kIStInitialState";
516  break;
517  case 1:
518  ParticleStatusName = "kIStFinalState";
519  break;
520  case 11:
521  ParticleStatusName = "kIStNucleonTarget";
522  break;
523  default:
524  ParticleStatusName = "Status Unknown";
525  }
526  return ParticleStatusName;
527  }
528 
529 
530 // //......................................................................
532  {
533  std::string ReactionChannelName=" ";
534 
535  if(ccnc==0)
536  ReactionChannelName = "kCC";
537  else if(ccnc==1)
538  ReactionChannelName = "kNC";
539  else std::cout<<"Current mode unknown!! "<<std::endl;
540 
541  if(mode==0)
542  ReactionChannelName += "_kQE";
543  else if(mode==1)
544  ReactionChannelName += "_kRes";
545  else if(mode==2)
546  ReactionChannelName += "_kDIS";
547  else if(mode==3)
548  ReactionChannelName += "_kCoh";
549  else if(mode==4)
550  ReactionChannelName += "_kNuElectronElastic";
551  else if(mode==5)
552  ReactionChannelName += "_kInverseMuDecay";
553  else std::cout<<"interaction mode unknown!! "<<std::endl;
554 
555  return ReactionChannelName;
556  }
557 
558 // //......................................................................
560  {
561  // Decide which histograms to put the spectrum in
562  int id = -1;
563  if (mc.GetNeutrino().CCNC()==simb::kCC) {
564  fCCMode->Fill(mc.GetNeutrino().Mode());
565  if (mc.GetNeutrino().Nu().PdgCode() == 12) id = 0;
566  else if (mc.GetNeutrino().Nu().PdgCode() == -12) id = 1;
567  else if (mc.GetNeutrino().Nu().PdgCode() == 14) id = 2;
568  else if (mc.GetNeutrino().Nu().PdgCode() == -14) id = 3;
569  else return;
570  }
571  else {
572  fNCMode->Fill(mc.GetNeutrino().Mode());
573  if (mc.GetNeutrino().Nu().PdgCode() > 0) id = 4;
574  else id = 5;
575  }
576  if (id==-1) abort();
577 
578  // Fill the specta histograms
579  fGenerated[id]->Fill(mc.GetNeutrino().Nu().E() );
580 
581  //< fill the vertex histograms from the neutrino - that is always
582  //< particle 0 in the list
583  simb::MCNeutrino mcnu = mc.GetNeutrino();
584  const simb::MCParticle nu = mcnu.Nu();
585 
586  fVertexX->Fill(nu.Vx());
587  fVertexY->Fill(nu.Vy());
588  fVertexZ->Fill(nu.Vz());
589 
590  fVertexXY->Fill(nu.Vx(), nu.Vy());
591  fVertexXZ->Fill(nu.Vz(), nu.Vx());
592  fVertexYZ->Fill(nu.Vz(), nu.Vy());
593 
594  double mom = nu.P();
595  if(std::abs(mom) > 0.){
596  fDCosX->Fill(nu.Px()/mom);
597  fDCosY->Fill(nu.Py()/mom);
598  fDCosZ->Fill(nu.Pz()/mom);
599  }
600 
601 
602 // MF_LOG_DEBUG("GENIEInteractionInformation")
603 // << std::endl
604 // << "REACTION: " << ReactionChannel(mc.GetNeutrino().CCNC(),mc.GetNeutrino().Mode())
605 // << std::endl
606 // << "-----------> Particles in the Stack = " << mc.NParticles() << std::endl
607 // << std::setiosflags(std::ios::left)
608 // << std::setw(20) << "PARTICLE"
609 // << std::setiosflags(std::ios::left)
610 // << std::setw(32) << "STATUS"
611 // << std::setw(18) << "E (GeV)"
612 // << std::setw(18) << "m (GeV/c2)"
613 // << std::setw(18) << "Ek (GeV)"
614 // << std::endl << std::endl;
615 
616 // const TDatabasePDG* databasePDG = TDatabasePDG::Instance();
617 
618 // // Loop over the particle stack for this event
619 // for(int i = 0; i < mc.NParticles(); ++i){
620 // simb::MCParticle part(mc.GetParticle(i));
621 // std::string name = databasePDG->GetParticle(part.PdgCode())->GetName();
622 // int code = part.StatusCode();
623 // std::string status = ParticleStatus(code);
624 // double mass = part.Mass();
625 // double energy = part.E();
626 // double Ek = (energy-mass); // Kinetic Energy (GeV)
627 // if(status=="kIStFinalStB4Interactions")
628 // MF_LOG_DEBUG("GENIEFinalState")
629 // << std::setiosflags(std::ios::left) << std::setw(20) << name
630 // << std::setiosflags(std::ios::left) << std::setw(32) <<status
631 // << std::setw(18)<< energy
632 // << std::setw(18)<< mass
633 // << std::setw(18)<< Ek <<std::endl;
634 // else
635 // MF_LOG_DEBUG("GENIEFinalState")
636 // << std::setiosflags(std::ios::left) << std::setw(20) << name
637 // << std::setiosflags(std::ios::left) << std::setw(32) << status
638 // << std::setw(18) << energy
639 // << std::setw(18) << mass <<std::endl;
640 
641  std::cout << "REACTION: " << ReactionChannel(mc.GetNeutrino().CCNC(),mc.GetNeutrino().Mode()) << std::endl;
642  std::cout << "-----------> Particles in the Stack = " << mc.NParticles() << std::endl;
643  std::cout << std::setiosflags(std::ios::left)
644  << std::setw(20) << "PARTICLE"
645  << std::setiosflags(std::ios::left)
646  << std::setw(32) << "STATUS"
647  << std::setw(18) << "E (GeV)"
648  << std::setw(18) << "m (GeV/c2)"
649  << std::setw(18) << "Ek (GeV)"
650  << std::endl << std::endl;
651 
652  const TDatabasePDG* databasePDG = TDatabasePDG::Instance();
653 
654  // Loop over the particle stack for this event
655  for(int i = 0; i < mc.NParticles(); ++i){
656  simb::MCParticle part(mc.GetParticle(i));
658  if (part.PdgCode() == 18040)
659  name = "Ar40 18040";
660  else if (part.PdgCode() != -99999 )
661  {
662  name = databasePDG->GetParticle(part.PdgCode())->GetName();
663  }
664 
665  int code = part.StatusCode();
666  std::string status = ParticleStatus(code);
667  double mass = part.Mass();
668  double energy = part.E();
669  double Ek = (energy-mass); // Kinetic Energy (GeV)
670 
671  std::cout << std::setiosflags(std::ios::left) << std::setw(20) << name
672  << std::setiosflags(std::ios::left) << std::setw(32) <<status
673  << std::setw(18)<< energy
674  << std::setw(18)<< mass
675  << std::setw(18)<< Ek <<std::endl;
676  }
677 
678 
679  if(mc.GetNeutrino().CCNC() == simb::kCC){
680 
681  ///look for the outgoing lepton in the particle stack
682  ///just interested in the first one
683  for(int i = 0; i < mc.NParticles(); ++i){
684  simb::MCParticle part(mc.GetParticle(i));
685  if(abs(part.PdgCode()) == 11){
686  fEMomentum->Fill(part.P());
687  fEDCosX->Fill(part.Px()/part.P());
688  fEDCosY->Fill(part.Py()/part.P());
689  fEDCosZ->Fill(part.Pz()/part.P());
690  fECons->Fill(nu.E() - part.E());
691  break;
692  }
693  else if(abs(part.PdgCode()) == 13){
694  fMuMomentum->Fill(part.P());
695  fMuDCosX->Fill(part.Px()/part.P());
696  fMuDCosY->Fill(part.Py()/part.P());
697  fMuDCosZ->Fill(part.Pz()/part.P());
698  fECons->Fill(nu.E() - part.E());
699  break;
700  }
701  }// end loop over particles
702  }//end if CC interaction
703 
704  return;
705  }
706 
707 }
708 
709 namespace evgen{
710 
712 
713 }
static QCString name
Definition: declinfo.cpp:673
double E(const int i=0) const
Definition: MCParticle.h:232
TH1F * fVertexX
vertex location of generated events in x
neutral current quasi-elastic
Definition: MCNeutrino.h:97
TH2F * fVertexYZ
vertex location in yz
int PdgCode() const
Definition: MCParticle.h:211
int CCNC() const
Definition: MCNeutrino.h:148
TH1F * fMuDCosZ
direction cosine of outgoing mu in z
const simb::MCNeutrino & GetNeutrino() const
Definition: MCTruth.h:70
std::string fNuanceFile
double Py(const int i=0) const
Definition: MCParticle.h:230
geo::Length_t DetHalfWidth(geo::TPCID const &tpcid) const
Returns the half width of the active volume of the specified TPC.
charged current deep inelastic scatter
Definition: MCNeutrino.h:134
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:256
NUANCEGen(fhicl::ParameterSet const &pset)
resonant charged current, nubar p -> l+ p pi-
Definition: MCNeutrino.h:107
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1055
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:75
neutrino electron elastic scatter
Definition: MCNeutrino.h:140
std::string string
Definition: nybbler.cc:12
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
TH2F * fVertexXY
vertex location in xy
TH1F * fDCosY
direction cosine in y
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:146
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
TH1F * fEMomentum
momentum of outgoing electrons
constexpr T pow(T x)
Definition: pow.h:72
double Px(const int i=0) const
Definition: MCParticle.h:229
TH1F * fECons
histogram to determine if energy is conserved in the event
double Mass(Resonance_t res)
resonance mass (GeV)
offset to account for adding in Nuance codes to this enum
Definition: MCNeutrino.h:95
TH1F * fGenerated[6]
Spectra as generated.
TH1F * fNCMode
CC interaction mode.
TH2F * fVertexXZ
vertex location in xz
int NParticles() const
Definition: MCTruth.h:68
void SetNeutrino(int CCNC, int mode, int interactionType, int target, int nucleon, int quark, double w, double x, double y, double qsqr)
Definition: MCTruth.cxx:30
charged current quasi-elastic
Definition: MCNeutrino.h:96
std::string ParticleStatus(int StatusCode)
void reconfigure(fhicl::ParameterSet const &p)
TH1F * fDCosZ
direction cosine in z
constexpr double PI
Definition: Run.h:21
TH1F * fCCMode
CC interaction mode.
TH1F * fDCosX
direction cosine in x
double y
TH1F * fEDCosX
direction cosine of outgoing e in x
T abs(T value)
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:56
resonant charged current, nu n -> l- n pi+
Definition: MCNeutrino.h:100
void FillHistograms(simb::MCTruth mc)
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
TH1F * fMuDCosY
direction cosine of outgoing mu in y
std::string ReactionChannel(int ccnc, int mode)
double P(const int i=0) const
Definition: MCParticle.h:233
double z
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
CodeOutputInterface * code
std::unique_ptr< std::ifstream > fEventFile
TH1F * fEDCosZ
direction cosine of outgoing e in z
void produce(art::Event &evt)
TH1F * fMuMomentum
momentum of outgoing muons
charged current deep inelastic scatter
Definition: MCNeutrino.h:133
resonant charged current, nu p -> l- p pi+
Definition: MCNeutrino.h:98
Q_EXPORT QTSManip setw(int w)
Definition: qtextstream.h:331
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:69
p
Definition: test.py:223
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
double Vx(const int i=0) const
Definition: MCParticle.h:220
constexpr auto const & left(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:96
void Add(simb::MCParticle const &part)
Definition: MCTruth.h:73
charged current coherent pion
Definition: MCNeutrino.h:139
inverse muon decay
Definition: MCNeutrino.h:141
A module to check the results from the Monte Carlo generator.
double Pz(const int i=0) const
Definition: MCParticle.h:231
resonant charged current, nubar n -> l+ n pi-
Definition: MCNeutrino.h:105
TH1F * fMuDCosX
direction cosine of outgoing mu in x
double Vz(const int i=0) const
Definition: MCParticle.h:222
EventNumber_t event() const
Definition: EventID.h:116
TH1F * fEDCosY
direction cosine of outgoing e in y
list x
Definition: train.py:276
TCEvent evt
Definition: DataStructs.cxx:5
Event generator information.
Definition: MCTruth.h:30
ProductID put(std::unique_ptr< PROD > &&edp, FullSemantic< Level::Run > const semantic)
Definition: DataViewImpl.h:730
Event generator information.
Definition: MCNeutrino.h:18
Namespace collecting geometry-related classes utilities.
int Mode() const
Definition: MCNeutrino.h:149
Event Generation using GENIE, cosmics or single particles.
std::string ParticleStatusName(int code)
Describes the status of a particle (simb::MCParticle::StatusCode()).
void beginRun(art::Run &run)
TH1F * fVertexY
vertex location of generated events in y
EventID id() const
Definition: Event.cc:37
double Vy(const int i=0) const
Definition: MCParticle.h:221
std::string fNUANCEModuleLabel
keep track of how long it takes to run the job
QTextStream & endl(QTextStream &s)
unsigned int run
TH1F * fVertexZ
vertex location of generated events in z
Beam neutrinos.
Definition: MCTruth.h:21