dEdx_module.cc
Go to the documentation of this file.
1 #ifndef dEdx_Module
2 #define dEdx_Module
3 
4 // LArSoft includes
11 
12 // Framework includes
17 #include "art_root_io/TFileService.h"
19 #include "canvas/Persistency/Common/FindManyP.h"
21 #include "fhiclcpp/ParameterSet.h"
22 
23 // ROOT includes.
24 #include "TF1.h"
25 #include "TH1.h"
26 #include "TH2.h"
27 #include "TH3.h"
28 #include "TTree.h"
29 #include "TLorentzVector.h"
30 #include "TVector3.h"
31 #include "THStack.h"
32 
33 // C++ Includes
34 #include <map>
35 #include <vector>
36 #include <algorithm>
37 #include <iostream>
38 #include <string>
39 #include <cmath>
40 
41 ///////////////////////////////////////////////////////////////////////
42 // This module is used to calculate dE/dx from simulated (cheated) //
43 // information. It was intended for a e/gamma separation study, //
44 // but can be adapted for similar purposes if desired. //
45 // //
46 // Kevin R. Wood -- krwood214@gmail.com (07/23/2014) //
47 ///////////////////////////////////////////////////////////////////////
48 
49 
50 namespace AnalysisExample
51 {
52 
53  class dEdx : public art::EDAnalyzer
54  {
55  public:
56 
57  explicit dEdx(fhicl::ParameterSet const& pset);
58  virtual ~dEdx();
59  void beginJob();
60  void beginRun(const art::Run& run);
61  void reconfigure(fhicl::ParameterSet const& pset);
62  void analyze(const art::Event& evt);
63  void endJob();
64  void findDepositionVertex(std::vector<sim::SimChannel> SCHandle);
65  void checkDepositionVertex(std::vector<sim::SimChannel> SCHandle, double x, double y, double z);
66 
67  private:
68 
72 
73  // dE/dx "cheated" information
74  TH1D* fMCdEdxHist;
75  std::vector<double> fMCdEdxVec;
76  double fMaxMCdEdx = 0.;
77  double dE;
78 
79  // distance of energy deposition from conversion point
80  double trackLength;
81 
82  // energy deposition 3d positions for conversion points, ~2.5cm of "track"
83  // upstream of conversion point (for dE/dx), and shower respectively
87 
88  // list of event numbers corresponding to unusually low dEdx measurement
89  std::vector<unsigned int> fNullList;
90 
91  // position of conversion point
92  double xi = 0.;
93  double yi = 0.;
94  double zi = 1000.;
95 
96  // for finding conversion point
97  double ziFalse = -500.;
98  unsigned int nearDeps = 0;
99  bool goodDepVertex = false;
100 
101  // for conversion distance
102  double fConvDist;
104  double Vxi; //
105  double Vyi; // primary vertex coordinates
106  double Vzi; //
107 
108  };
109 
110  //-----------------------------------------------------------------------
111 
112  dEdx::dEdx(fhicl::ParameterSet const& parameterSet)
113  : EDAnalyzer(parameterSet)
114  {
115  this->reconfigure(parameterSet);
116  }
117 
118  //-----------------------------------------------------------------------
119 
121  {
122  }
123 
124  //-----------------------------------------------------------------------
125 
127  {
128 
129  // hard-coded fd4apa cryostat dimensions
130  fStartEdepHist = tfs->make<TH3D>("fStartEdepHist","fStartEdepHist",50,-240.,240.,50,-750.,800.,50,-50.,560.);
131  fTrackEdepHist = tfs->make<TH3D>("fTrackEdepHist","fTrackEdepHist",50,-240.,240.,50,-750.,800.,50,-50.,560.);
132  fShowerEdepHist = tfs->make<TH3D>("fShowerEdepHist","fShowerEdepHist",50,-240.,240.,50,-750.,800.,50,-50.,560.);
133 
134  fStartEdepHist->GetXaxis()->SetTitle("x");
135  fStartEdepHist->GetYaxis()->SetTitle("y");
136  fStartEdepHist->GetZaxis()->SetTitle("z");
137 
138  fShowerEdepHist->GetXaxis()->SetTitle("x");
139  fShowerEdepHist->GetYaxis()->SetTitle("y");
140  fShowerEdepHist->GetZaxis()->SetTitle("z");
141 
142  fTrackEdepHist->GetXaxis()->SetTitle("x");
143  fTrackEdepHist->GetYaxis()->SetTitle("y");
144  fTrackEdepHist->GetZaxis()->SetTitle("z");
145 
146  fConvDistHist = tfs->make<TH1D>("fConvDistHist","fConvDistHist",100,0.,50.);
147 
148  }
149 
150  //-----------------------------------------------------------------------
151 
153  {
154  }
155 
156  //-----------------------------------------------------------------------
157 
159  {
160  fSimulationProducerLabel = p.get<std::string>("SimulationLabel");
161  }
162 
163  //-----------------------------------------------------------------------
164 
166  {
167 
168  auto simChanHandle = event.getHandle<std::vector<sim::SimChannel>>(fSimulationProducerLabel);
169 
170  // initialize at the start of every event
171  dE = 0.;
172  trackLength = 0.;
173  goodDepVertex = false;
174  unsigned int nTries = 0;
175  ziFalse = -500.;
176 
177  // find conversion point
178  while(!goodDepVertex && nTries < 300)
179  {
180  this->findDepositionVertex( (*simChanHandle) );
181  this->checkDepositionVertex( (*simChanHandle), xi, yi, zi );
182  nTries++;
183  }
184 
185  std::cout << "Number of attemps at finding deposition vertex = " << nTries << "." << std::endl;
186 
187  // find energy deposited into first 2.5 cm of "track"
188  if(nTries < 300 && goodDepVertex)
189  {
190  fStartEdepHist->Fill(xi,yi,zi);
191  for(auto const& channel : (*simChanHandle))
192  {
193  if(fGeom->SignalType(channel.Channel()) == geo::kCollection)
194  {
195  auto const& timeSlices = channel.TDCIDEMap();
196  for(auto const& t : timeSlices)
197  {
198  auto const& eDeps = t.second;
199  for(auto const& eDep : eDeps)
200  {
201  fShowerEdepHist->Fill(eDep.x,eDep.y,eDep.z);
202  trackLength = std::sqrt( (eDep.x-xi)*(eDep.x-xi) + (eDep.y-yi)*(eDep.y-yi) + (eDep.z-zi)*(eDep.z-zi) );
203  if(trackLength < 2.5 /* && eDep.z >= zi*/)
204  {
205  dE += eDep.energy;
206  fTrackEdepHist->Fill(eDep.x,eDep.y,eDep.z);
207  }
208  } // every energy deposition
209  } // in every time slice
210  } // for the collection plane
211  } // and for every wire
212  dE /= 2.5;
213  fMCdEdxVec.push_back(dE);
214  if(dE > fMaxMCdEdx) fMaxMCdEdx = dE;
215  if(dE < 0.5) fNullList.push_back(event.id().event());
216  }
217 
218  auto particleHandle = event.getHandle< std::vector<simb::MCParticle> >(fSimulationProducerLabel);
219 
220  for(auto const& particle : (*particleHandle))
221  {
222  if(particle.Process() == "primary")
223  {
224  Vxi = particle.Vx();
225  Vyi = particle.Vy();
226  Vzi = particle.Vz();
227  }
228  }
229  fConvDist = std::sqrt( (Vxi-xi)*(Vxi-xi) + (Vyi-yi)*(Vyi-yi) + (Vzi-zi)*(Vzi-zi) );
230  fConvDistHist->Fill(fConvDist);
231  }
232 
233  //-----------------------------------------------------------------------
234 
236  {
237  fMCdEdxHist = tfs->make<TH1D>("fMCdEdxHist",";dE/dx (MeV/cm)",200,0.0,10.);
238 
239  for(unsigned int e = 0; e < fMCdEdxVec.size(); e++)
240  {
241  fMCdEdxHist->Fill(fMCdEdxVec[e]);
242  }
243 
244  std::cout << "Events with < 0.5 dE/dx entry: ";
245  for(unsigned int a = 0; a < fNullList.size(); a++) std::cout << fNullList[a] << ", ";
246 
247  }
248 
249  void dEdx::findDepositionVertex(std::vector<sim::SimChannel> SCHandle)
250  {
251  zi = 1000.;
252  std::cout << "Finding Deposition Vertex... " << std::endl;
253  for(auto const& channel : SCHandle)
254  {
255  if(fGeom->SignalType(channel.Channel()) == geo::kCollection)
256  {
257  auto const& timeSlices = channel.TDCIDEMap();
258  for(auto const& t : timeSlices)
259  {
260  auto const& eDeps = t.second;
261  for(auto const& eDep : eDeps)
262  {
263  if( (eDep.z < zi) && (eDep.energy > 0.) && (eDep.z > ziFalse) )
264  {
265  xi = eDep.x;
266  yi = eDep.y;
267  zi = eDep.z;
268  } // if z position is a "valid" minimum
269  } // energy deposition loop
270  } // time slice loop
271  } // if collection wire
272  } // sim channel loop
273  } // findDepositionVertex
274 
275  void dEdx::checkDepositionVertex(std::vector<sim::SimChannel> SCHandle, double x,double y, double z)
276  {
277  std::cout << "Checking for Deposition Vertex" << std::endl;
278  nearDeps = 0;
279  double track = 0;
280  for(auto const& channel : SCHandle)
281  {
282  if(fGeom->SignalType(channel.Channel()) == geo::kCollection)
283  {
284  auto const& timeSlices = channel.TDCIDEMap();
285  for(auto const& t : timeSlices)
286  {
287  auto const& eDeps = t.second;
288  for(auto const& eDep : eDeps)
289  {
290  track = std::sqrt( (eDep.x-x)*(eDep.x-x) + (eDep.y-y)*(eDep.y-y) + (eDep.z-z)*(eDep.z-z) );
291  if(track < 2.5 && eDep.energy > 0.)
292  {
293  nearDeps++;
294  } // if z position is minimum
295  } // energy deposition loop
296  } // time slice loop
297  } // if collection wire
298  } // sim channel loop
299  if( nearDeps > 30) goodDepVertex = true;
300  else ziFalse = z + 0.05;
301  } // checkDepositionVertex
302 
304 
305 } // namespace AnalysisExample
306 
307 #endif // dEdx_Module
art::ServiceHandle< geo::Geometry > fGeom
Definition: dEdx_module.cc:70
Store parameters for running LArG4.
dEdx(fhicl::ParameterSet const &pset)
Definition: dEdx_module.cc:112
std::string string
Definition: nybbler.cc:12
unsigned int nearDeps
Definition: dEdx_module.cc:98
std::string fSimulationProducerLabel
Definition: dEdx_module.cc:71
uint8_t channel
Definition: CRTFragment.hh:201
SigType_t SignalType(geo::PlaneID const &pid) const
Returns the type of signal on the channels of specified TPC plane.
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
void checkDepositionVertex(std::vector< sim::SimChannel > SCHandle, double x, double y, double z)
Definition: dEdx_module.cc:275
Particle class.
art framework interface to geometry description
Definition: Run.h:17
void reconfigure(fhicl::ParameterSet const &pset)
Definition: dEdx_module.cc:158
const double e
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
void findDepositionVertex(std::vector< sim::SimChannel > SCHandle)
Definition: dEdx_module.cc:249
art::ServiceHandle< art::TFileService > tfs
Definition: dEdx_module.cc:69
const double a
T get(std::string const &key) const
Definition: ParameterSet.h:271
void analyze(const art::Event &evt)
Definition: dEdx_module.cc:165
p
Definition: test.py:223
Definition of data types for geometry description.
void beginRun(const art::Run &run)
Definition: dEdx_module.cc:152
EventNumber_t event() const
Definition: EventID.h:116
list x
Definition: train.py:276
std::vector< unsigned int > fNullList
Definition: dEdx_module.cc:89
TCEvent evt
Definition: DataStructs.cxx:7
EventID id() const
Definition: Event.cc:34
QTextStream & endl(QTextStream &s)
Event finding and building.
std::vector< double > fMCdEdxVec
Definition: dEdx_module.cc:75
Signal from collection planes.
Definition: geo_types.h:146