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