LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
LArG4Ana_module.cc
Go to the documentation of this file.
1 #ifndef LARG4_LARG4ANA_H
8 #define LARG4_LARG4ANA_H
9 
13 
14 // Framework includes
16 #include "fhiclcpp/ParameterSet.h"
24 #include "cetlib_except/exception.h"
25 
26 // LArSoft Includes
32 
33 // ROOT includes
34 #include "TH1.h"
35 #include "TH2.h"
36 #include "TProfile.h"
37 #include "TLorentzVector.h"
38 #include <TTree.h>
39 
40 // C++ Includes
41 #include <iostream>
42 #include <cstring>
43 #include <sys/stat.h>
44 
45 namespace simb{
46  class MCTruth;
47 }
48 
49 namespace sim{
50  class ParticleList;
51 }
52 
54 namespace larg4 {
55 
56  class LArG4Ana : public art::EDAnalyzer{
57  public:
58 
60  explicit LArG4Ana(fhicl::ParameterSet const& pset);
61  virtual ~LArG4Ana();
62 
63  void analyze (const art::Event& evt);
64  void beginJob();
65  void reconfigure(fhicl::ParameterSet const& pset);
66 
67  private:
68 
69  std::string fG4ModuleLabel;
70  std::string fTruthModuleLabel;
71 
72 
73  TH1D *fPDGCodes;
74  TH1D *fPi0Momentum;
75  TH1D *fnEnergy;
76  TH1D *fnDist;
77  TH1D *fnumChannels;
78  TProfile *fnumIDEs;
79  TH1D *fEventCharge;
80  TH1D *fEventEnergy;
81  TProfile *fChannelCharge;
82  TProfile *fChannelEnergy;
83 
84  // Int_t stringDim = 35;
85 
86  TTree *fTree;
87  Int_t fTEvt;
88  Int_t fTSub;
89  Int_t fTRun;
90  Int_t fTPdg;
91  Int_t fTID;
92  Int_t fTNds;
94  Int_t fTNds4;
95  Int_t *fTDID;
96  Int_t *fTDPdg;
97  Float_t *fTDWt;
98  Char_t fTProcess[35];
99  Char_t fTVolume[35];
100  Char_t fTTVolume[35]; // Termination Volume
101  Char_t fTMaterial[35];
102  Char_t fTDProcess[200][35];
103  Int_t fTParentID;
104  Int_t fTStatus;
105  Float_t fTWeight;
106  Float_t* fT4Origin;
107  Float_t* fT4DOrigin;
108  Float_t* fT4Termination; // Termination Coordinates
109  Float_t* fT4Momentum;
110  Float_t* fT4DMomentum;
111  };
112 
113 } // namespace larg4
114 
115 namespace larg4 {
116 
117  //-----------------------------------------------------------------------
118  // Constructor
119  LArG4Ana::LArG4Ana(fhicl::ParameterSet const& pset)
120  : EDAnalyzer(pset)
121  {
122  this->reconfigure(pset);
123  }
124 
125  //-----------------------------------------------------------------------
126  // Destructor
128  {
129  }
130 
131  //-----------------------------------------------------------------------
133  {
136 
137  fPDGCodes = tfs->make<TH1D>("pdgcodes", ";PDG Code;", 5000, -2500, 2500);
138  fPi0Momentum = tfs->make<TH1D>("pi0mom", ";#pi^{0} Momentum (GeV);", 1000, 0., 1000.);
139 
140  fTree = tfs->make<TTree>("MCTTree","MCTTree");
141  fnEnergy = tfs->make<TH1D>("nEnergy", ";n,#Lambda^{0},K^{0} Momentum (GeV);", 100, 0., 10.);
142  fnDist = tfs->make<TH1D>("nDistance", ";n,#Lambda^{0},K^{0} Distance (m);", 200, -30000.0, +30000.);
143 
144 
145  // Some histograms relating to drift electrons, active detector
146  // channels and charge/energy on channels
147  fnumChannels = tfs->make<TH1D>("fnumChannels",
148  "Active channels;Active channels;# events",
149  256, 0, geo->Nchannels());
150  fnumIDEs = tfs->make<TProfile>("fnumIDEs",
151  "Drift Electrons per channel;Channel;Drift electrons",
152  geo->Nchannels()+1, 0, geo->Nchannels(),
153  0, 1e4);
154  fEventCharge = tfs->make<TH1D>("fEventCharge",
155  "Charge in event;Total charge per event;# events",
156  100, 0, 2.5e8);
157  fEventEnergy = tfs->make<TH1D>("fEventEnergy",
158  "Energy in event;Total energy per event;# events",
159  100, 0, 1e4);
160  fChannelCharge = tfs->make<TProfile>("fChannelCharge",
161  "Charge on channel;Channel;Total charge per channel",
162  geo->Nchannels()+1,0,geo->Nchannels(),
163  0,1e5);
164  fChannelEnergy = tfs->make<TProfile>("fChannelEnergy",
165  "Energy on channel;Channel;Total energy per channel",
166  geo->Nchannels()+1,0,geo->Nchannels(),
167  0, 1e3);
168 
169 
170  fT4Origin = new Float_t[4];
171  fT4DOrigin = new Float_t[fTNds*4];
172  fT4Termination = new Float_t[4];
173  fT4Momentum = new Float_t[4];
174  fT4DMomentum = new Float_t[fTNds*4];
175  fTDID = new Int_t[fTNds];
176  fTDPdg = new Int_t[fTNds];
177  fTDWt = new Float_t[fTNds];
178  fTNds4 = fTNds*4; // TTree/Branch requirement to store this.
179 
180  fTree->Branch("MCEvt", &fTEvt, "MCEvt/I");
181  fTree->Branch("MCSub", &fTSub, "MCSub/I");
182  fTree->Branch("MCRun", &fTRun, "MCRun/I");
183  fTree->Branch("MCWt", &fTWeight, "MCWt/F");
184  fTree->Branch("MCPdg", &fTPdg, "MCPdg/I");
185  fTree->Branch("MCID", &fTID, "MCID/I");
186  fTree->Branch("MCParentID", &fTParentID, "MCParentID/I");
187  fTree->Branch("MCNumDs", &fTNds, "MCNumDs/I");
188  fTree->Branch("MCNumDs4", &fTNds4, "MCNumDs4/I");
189  fTree->Branch("MCDID", fTDID, "MCDID[MCNumDs]/I");
190  fTree->Branch("MCDPdg", fTDPdg, "MCDPdg[MCNumDs]/I");
191  fTree->Branch("MCDWt", fTDWt, "MCDWt[MCNumDs]/I");
192  fTree->Branch("MCProcess", fTProcess, "MCProcess/C");
193  fTree->Branch("MCVolume", fTVolume, "MCVolume/C");
194  fTree->Branch("MCTVolume", fTTVolume, "MCTVolume/C");
195  fTree->Branch("MCMaterial", fTMaterial, "MCMaterial/C");
196  fTree->Branch("MCDProcess", fTDProcess, "MCDProcess[MCNumDs]/C");
197  fTree->Branch("MCStatus", &fTStatus, "MCStatus/I");
198  fTree->Branch("MCOrigin", fT4Origin, "MCOrigin[4]/F");
199  fTree->Branch("MCDOrigin", fT4DOrigin, "MCDOrigin[MCNumDs4]/F");
200  fTree->Branch("MCTermination", fT4Termination, "MCTermination[4]/F");
201  fTree->Branch("MCMomentum", fT4Momentum, "MCMomentum[4]/F");
202  fTree->Branch("MCDMomentum", fT4DMomentum, "MCDMomentum[MCNumDs4]/F");
203 
204  }
205 
206  //-----------------------------------------------------------------------
208  {
209  fG4ModuleLabel = p.get< std::string >("GeantModuleLabel");
210  fTNdsOriginal = p.get< int >("Ndaughters" );
211  fTruthModuleLabel = p.get< std::string >("TruthModuleLabel");
213 
214  return;
215  }
216 
217  //-----------------------------------------------------------------------
219  {
220 
221  //get the list of particles from this event
223  const sim::ParticleList& plist = pi_serv->ParticleList();
225 
226  // loop over all sim::SimChannels in the event and make sure there are no
227  // sim::IDEs with trackID values that are not in the sim::ParticleList
228  std::vector<const sim::SimChannel*> sccol;
229  evt.getView(fG4ModuleLabel, sccol);
230 
231  double totalCharge=0.0;
232  double totalEnergy=0.0;
233  fnumChannels->Fill(sccol.size());
234  for(size_t sc = 0; sc < sccol.size(); ++sc){
235  double numIDEs=0.0;
236  double scCharge=0.0;
237  double scEnergy=0.0;
238  const auto & tdcidemap = sccol[sc]->TDCIDEMap();
239  for(auto mapitr = tdcidemap.begin(); mapitr != tdcidemap.end(); mapitr++){
240  const std::vector<sim::IDE> idevec = (*mapitr).second;
241  numIDEs += idevec.size();
242  for(size_t iv = 0; iv < idevec.size(); ++iv){
243  if(plist.find( idevec[iv].trackID ) == plist.end()
244  && idevec[iv].trackID != sim::NoParticleId)
245  mf::LogWarning("LArG4Ana") << idevec[iv].trackID << " is not in particle list";
246  totalCharge +=idevec[iv].numElectrons;
247  scCharge += idevec[iv].numElectrons;
248  totalEnergy +=idevec[iv].energy;
249  scEnergy += idevec[iv].energy;
250  }
251  }
252  fnumIDEs->Fill(sc,numIDEs);
253  fChannelCharge->Fill(sc,scCharge);
254  fChannelEnergy->Fill(sc,scEnergy);
255  }
256  fEventCharge->Fill(totalCharge);
257  fEventEnergy->Fill(totalEnergy);
258 
259  // get the particles from the back tracker
260  const sim::ParticleList& Particles = pi_serv->ParticleList();
261  std::vector<const simb::MCParticle*> pvec;
262  pvec.reserve(Particles.size());
263  for (const auto& PartPair: Particles) {
264  pvec.push_back(PartPair.second);
265  fPDGCodes->Fill(PartPair.second->PdgCode());
266  }
267 
268  // now look for pi0's that decay to 2 gammas
269  int pi0loc = -1;
270  int numpi0gamma = 0;
271  for(unsigned int i = 0; i < pvec.size(); ++i){
272  if(pvec[i]->PdgCode() == 111) pi0loc = i;
273  if(pvec[i]->Mother() == pi0loc+1 &&
274  pi0loc > 0 &&
275  pvec[i]->PdgCode() == 22){
276  mf::LogInfo("LArG4Ana") << pvec[i]->E() << " gamma energy ";
277  ++numpi0gamma;
278  }
279 
280  // n,Lambda,K0s,K0L,K0
281  if (pvec[i]->PdgCode() == 2112 ||
282  pvec[i]->PdgCode() == 3122 ||
283  pvec[i]->PdgCode() == 130 ||
284  pvec[i]->PdgCode() == 310 ||
285  pvec[i]->PdgCode() == 311 ) {
286  fnEnergy->Fill(pvec[i]->E(),pvec[i]->Weight());
287  fnDist->Fill(pvec[i]->Vx(),pvec[i]->Weight());
288  }
289 
290  fTPdg = pvec[i]->PdgCode();
291  fTID = pvec[i]->TrackId();
292  // 0 out strings, else there may be cruft in here from prev evt.
293  for (unsigned int s = 0; s < 35; ++s){
294  *(fTProcess+s) = 0;
295  *(fTProcess+s) = 0;
296  *(fTMaterial+s) = 0;
297  *(fTMaterial+s) = 0;
298  *(fTVolume+s) = 0;
299  *(fTVolume+s) = 0;
300  *(fTTVolume+s) = 0;
301  *(fTTVolume+s) = 0;
302  }
303 
304  for(unsigned int s = 0; s < pvec[i]->Process().length(); ++s) *(fTProcess+s) = pvec[i]->Process()[s];
305 
306  TVector3 dum = pvec[i]->Position().Vect();
307 
308  for (unsigned int s = 0; s < geom->MaterialName(pvec[i]->Position().Vect()).length(); ++s)
309  *(fTMaterial+s) = geom->MaterialName(pvec[i]->Position().Vect())[s];
310 
311  for (unsigned int s = 0; s < geom->VolumeName(pvec[i]->Position().Vect()).length(); ++s)
312  *(fTVolume+s) = geom->VolumeName(pvec[i]->Position().Vect())[s];
313 
314  for (unsigned int s = 0; s < geom->VolumeName(pvec[i]->EndPosition().Vect()).length(); ++s)
315  *(fTTVolume+s) = geom->VolumeName(pvec[i]->EndPosition().Vect())[s];
316 
317  fTEvt = evt.id().event();
318  fTSub = evt.subRun();
319  fTRun = evt.run();
320  fTParentID = pvec[i]->Mother();
321  fTStatus = pvec[i]->StatusCode();
322  int daughter = 9999;
323  fTNds = TMath::Min(pvec[i]->NumberDaughters(),fTNdsOriginal);
324  for( int d = 0; d < fTNds; d++ ){
325  daughter = pvec[i]->Daughter(d);
326  fTDID[d] = daughter;
327  // zero it out.
328  for (unsigned int s = 0; s < 35; ++s) *(fTDProcess[d]+s) = 0;
329 
330  for(unsigned int jj = i; jj < pvec.size(); ++jj){ // Don't look below i.
331 
332  if (fTDID[d] == pvec[jj]->TrackId()){
333  fTDPdg[d] = pvec[jj]->PdgCode(); // get the pointer,
334  fTDWt[d] = pvec[jj]->Weight();
335 
336  for (unsigned int s = 0; s < pvec[jj]->Process().length(); ++s)
337  *(fTDProcess[d]+s) = pvec[jj]->Process()[s];
338 
339  for (unsigned int kk = 0; kk < 4; ++kk){
340  fT4DOrigin[d*4+kk] = pvec[jj]->Position()[kk];
341  fT4DMomentum[d*4+kk] = pvec[jj]->Momentum()[kk];
342  }
343  break;
344  }
345  }
346  }//end loop over d
347 
348  for (unsigned int ii = 0; ii < 4; ++ii){
349  fT4Termination[ii] = 1e9;
350  fT4Origin[ii] = pvec[i]->Position()[ii];
351  if (ii!=3) fT4Termination[ii] = pvec[i]->EndPosition()[ii];
352  if (ii==4) fT4Termination[ii] = pvec[i]->Momentum()[ii]; // yes, odd
353  fT4Momentum[ii] = pvec[i]->Momentum()[ii];
354  }
355 
356  fTWeight = pvec[i]->Weight();
357  fTree->Fill();
358 
359  } // end loop on particles in list
360  if(numpi0gamma == 2 && pi0loc > 0){
361  mf::LogInfo("LArG4Ana") << pvec[pi0loc]->E();
362  fPi0Momentum->Fill(pvec[pi0loc]->E());
363  }
364 
365  return;
366  }
367 
368 
369 
370 } // namespace larg4
371 
372 namespace larg4 {
373 
375 
376 } // namespace LArG4
377 
378 #endif // LARG4_LARG4_H
379 
Float_t s
Definition: plot.C:23
SubRunNumber_t subRun() const
Definition: Event.h:72
TProfile * fChannelCharge
Charge per channel.
void analyze(const art::Event &evt)
std::string fG4ModuleLabel
module label for the Geant
Float_t E
Definition: plot.C:23
void reconfigure(fhicl::ParameterSet const &pset)
TH1D * fnumChannels
The number of channels recieving charge per event.
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Float_t * fT4Termination
Char_t fTTVolume[35]
Geant4 interface.
Char_t fTMaterial[35]
TH1D * fEventEnergy
Energy collected per event.
Char_t fTVolume[35]
unsigned int Nchannels() const
Returns the number of TPC readout channels in the detector.
iterator find(const key_type &key)
Definition: ParticleList.h:318
TProfile * fChannelEnergy
Energy per channel.
std::size_t getView(std::string const &moduleLabel, std::string const &productInstanceName, std::vector< ELEMENT const * > &result) const
Definition: DataViewImpl.h:474
Char_t fTDProcess[200][35]
Float_t * fT4DMomentum
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
void beginJob()
Definition: Breakpoints.cc:14
static const int NoParticleId
Definition: sim.h:28
T get(std::string const &key) const
Definition: ParameterSet.h:231
Float_t d
Definition: plot.C:237
Char_t fTProcess[35]
TH1D * fEventCharge
Charge collected per event.
Monte Carlo Simulation.
Framework includes.
T * make(ARGS...args) const
std::string fTruthModuleLabel
module label for the Geant
object containing MC truth information necessary for making RawDigits and doing back tracking ...
EventNumber_t event() const
Definition: EventID.h:117
std::string MaterialName(TVector3 const &point) const
Name of the deepest material containing the point xyz.
TCEvent evt
Definition: DataStructs.cxx:5
Particle list in DetSim contains Monte Carlo particle information.
RunNumber_t run() const
Definition: Event.h:77
Namespace collecting geometry-related classes utilities.
size_type size() const
Definition: ParticleList.h:313
Tools and modules for checking out the basics of the Monte Carlo.
std::string VolumeName(geo::Point_t const &point) const
Returns the name of the deepest volume containing specified point.
EventID id() const
Definition: Event.h:56
art framework interface to geometry description
TProfile * fnumIDEs
Number of drift electrons per channel.