All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MCRecoEdep.cxx
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // MCRecoEdep source
4 //
5 ////////////////////////////////////////////////////////////////////////
6 
7 #include "art/Framework/Services/Registry/ServiceHandle.h"
8 #include "messagefacility/MessageLogger/MessageLogger.h"
9 #include "fhiclcpp/ParameterSet.h"
10 
15 
16 #include "MCRecoEdep.h"
17 
18 #include <iostream>
19 #include <map>
20 #include <utility>
21 #include <vector>
22 
23 namespace sim {
24 
25  namespace details {
26  std::map<geo::PlaneID, size_t> createPlaneIndexMap(){
27  art::ServiceHandle<geo::Geometry const> geom;
28  std::map<geo::PlaneID, size_t> m;
29  size_t i = 0;
30  for(auto const& pid : geom->IteratePlaneIDs()){
31  m[pid] = i;
32  i++;
33  }
34  return m;
35  }
36  }
37  //const unsigned short MCEdepHit::kINVALID_USHORT = std::numeric_limits<unsigned short>::max();
38 
39  //const short MCEdep::kINVALID_SHORT = std::numeric_limits<short>::max();
40 
41  //##################################################################
42  MCRecoEdep::MCRecoEdep(fhicl::ParameterSet const& pset)
43  //##################################################################
44  {
45  _debug_mode = pset.get<bool>("DebugMode");
46  _save_mchit = pset.get<bool>("SaveMCHit");
47  }
48 
49  const std::vector<sim::MCEdep>& MCRecoEdep::GetEdepArrayAt(size_t edep_index) const
50  {
51  if(edep_index >= _mc_edeps.size())
52  throw cet::exception(__FUNCTION__) << Form("Track ID %zu not found!",edep_index);
53  return _mc_edeps.at(edep_index);
54  }
55 
56  std::vector<sim::MCEdep>& MCRecoEdep::__GetEdepArray__(unsigned int track_id)
57  {
58  if(ExistTrack(track_id)) return _mc_edeps.at((*_track_index.find(track_id)).second);
59  _track_index.insert(std::pair<unsigned int,size_t>(track_id,_mc_edeps.size()));
60  _mc_edeps.push_back(std::vector<sim::MCEdep>());
61  return (*(_mc_edeps.rbegin()));
62  }
63 
64  void MCRecoEdep::MakeMCEdep(const std::vector<sim::SimChannel>& schArray)
65  {
66  _mc_edeps.clear();
67  _track_index.clear();
68 
69  art::ServiceHandle<geo::Geometry const> geom;
70 // const detinfo::DetectorProperties* detp = lar::providerFrom<detinfo::DetectorPropertiesService>();
71 
72  // Key map to identify a unique particle energy deposition point
73  std::map<std::pair<UniquePosition, unsigned int>, int> hit_index_m;
74 
75  auto pindex = details::createPlaneIndexMap();
76 
77  if(_debug_mode) std::cout<<"Processing "<<schArray.size()<<" channels..."<<std::endl;
78  // Loop over channels
79  for(size_t i=0; i<schArray.size(); ++i) {
80 
81  // Get data to loop over
82  auto const& sch = schArray[i];
83  const auto &sch_map(sch.TDCIDEMap());
84  // Channel
85  UInt_t ch = sch.Channel();
86  // Loop over ticks
87  for(auto tdc_iter = sch_map.begin(); tdc_iter!=sch_map.end(); ++tdc_iter) {
88  // for c2: hit_time is unused
89  //unsigned short hit_time = (*tdc_iter).first;
90  // Loop over IDEs
91  for(auto const &ide : (*tdc_iter).second) {
92 
93  int track_id = ide.trackID;
94  if(track_id < 0) track_id = track_id * (-1);
95  unsigned int real_track_id = track_id;
96 
97  UniquePosition pos(ide.x, ide.y, ide.z);
98 
99  int hit_index = -1;
100  auto key = std::make_pair(pos, real_track_id);
101  auto hit_index_track_iter = hit_index_m.find(key);
102  if(hit_index_track_iter == hit_index_m.end()) {
103  int new_hit_index = this->__GetEdepArray__(real_track_id).size();
104  hit_index_m[key]= new_hit_index;
105  }
106  else {
107  hit_index = (*hit_index_track_iter).second;
108  }
109  auto const pid = geom->ChannelToWire(ch)[0].planeID();
110  auto const channel_id = pindex[pid];
111  double charge = ide.numElectrons;
112  if(hit_index < 0) {
113  // This particle energy deposition is never recorded so far. Create a new Edep
114  //float charge = ide.numElectrons * detp->ElectronsToADC();
115  this->__GetEdepArray__(real_track_id).emplace_back(pos, pid, pindex.size(), ide.energy, charge, channel_id);
116  } else {
117  // Append charge to the relevant edep (@ hit_index)
118  //float charge = ide.numElectrons * detp->ElectronsToADC();
119  MCEdep &edep = this->__GetEdepArray__(real_track_id).at(hit_index);
120  edep.deps[channel_id].charge += charge;
121  edep.deps[channel_id].energy += ide.energy;
122  }
123  } // end looping over ides in this tick
124  } // end looping over ticks in this channel
125  }// end looping over channels
126 
127  if(_debug_mode) {
128  std::cout<< Form(" Collected %zu particles' energy depositions...",_mc_edeps.size()) << std::endl;
129  // for c2: disable the entire loop instead of just the print statement
130  //for(auto const& track_id_index : _track_index ) {
131  //auto track_id = track_id_index.first;
132  //auto edep_index = track_id_index.second;
133  // std::cout<< Form(" Track ID: %d ... %zu Edep!", track_id, edep_index) << std::endl;
134  //}
135  std::cout<<std::endl;
136  }
137  }
138 
139  void MCRecoEdep::MakeMCEdep(const std::vector<sim::SimEnergyDeposit>& sedArray)
140  {
141  _mc_edeps.clear();
142  _track_index.clear();
143 
144  art::ServiceHandle<geo::Geometry const> geom;
145 // const detinfo::DetectorProperties* detp = lar::providerFrom<detinfo::DetectorPropertiesService>();
146 
147  // Key map to identify a unique particle energy deposition point
148  std::map<std::pair<UniquePosition, unsigned int>, int> hit_index_m;
149 
150  auto pindex = details::createPlaneIndexMap();
151 
152  if(_debug_mode) std::cout<<"Processing "<<sedArray.size()<<" energy deposits..."<<std::endl;
153  // Loop over channels
154  for(size_t i=0; i<sedArray.size(); ++i) {
155 
156  // Get data to loop over
157  auto const& sed = sedArray[i];
158 
159  // David Caratelli:
160  // much of the code below is taken from the module:
161  // https://cdcvs.fnal.gov/redmine/projects/larsim/repository/revisions/develop/entry/larsim/ElectronDrift/SimDriftElectrons_module.cc
162 
163  // given this SimEnergyDeposit, find the TPC channel information
164  // "xyz" is the position of the energy deposit in world
165  // coordinates. Note that the units of distance in
166  // sim::SimEnergyDeposit are supposed to be cm.
167  auto const mp = sed.MidPoint();
168  double const xyz[3] = { mp.X(), mp.Y(), mp.Z() };
169  // From the position in world coordinates, determine the
170  // cryostat and tpc. If somehow the step is outside a tpc
171  // (e.g., cosmic rays in rock) just move on to the next one.
172  unsigned int cryostat = 0;
173  try {
174  geom->PositionToCryostat(xyz, cryostat);
175  }
176  catch(cet::exception &e){
177  mf::LogWarning("SimDriftElectrons") << "step "// << energyDeposit << "\n"
178  << "cannot be found in a cryostat\n"
179  << e;
180  continue;
181  }
182  unsigned int tpc = 0;
183  try {
184  geom->PositionToTPC(xyz, tpc, cryostat);
185  }
186  catch(cet::exception &e){
187  mf::LogWarning("SimDriftElectrons") << "step "// << energyDeposit << "\n"
188  << "cannot be found in a TPC\n"
189  << e;
190  continue;
191  }
192  const geo::TPCGeo& tpcGeo = geom->TPC(tpc, cryostat);
193 
194  //Define charge drift direction: driftcoordinate (x, y or z) and driftsign (positive or negative). Also define coordinates perpendicular to drift direction.
195  // unused int driftcoordinate = std::abs(tpcGeo.DetectDriftDirection())-1; //x:0, y:1, z:2
196 
197  // make a collection of electrons for each plane
198  for(size_t p = 0; p < tpcGeo.Nplanes(); ++p){
199 
200  // grab the nearest channel to the fDriftClusterPos position
201  // David Caratelli, comment begin:
202  // NOTE: the below code works only when the drift coordinate is indeed in x (i.e. 0th coordinate)
203  // see code linked above for a much more careful treatment of the coordinate system
204  // David Caratelli, comment end.
205  raw::ChannelID_t ch;
206  try {
207  ch = geom->NearestChannel(xyz, p, tpc, cryostat);
208  }
209  catch(cet::exception &e){
210  mf::LogWarning("SimDriftElectrons") << "step "// << energyDeposit << "\n"
211  << "nearest wire not in TPC\n"
212  << e;
213  continue;
214  }
215 
216  int track_id = sed.TrackID();
217 
218  if(track_id < 0) track_id = track_id * (-1);
219  unsigned int real_track_id = track_id;
220 
221  UniquePosition pos(xyz[0], xyz[1], xyz[2]);
222 
223  int hit_index = -1;
224  auto key = std::make_pair(pos, real_track_id);
225  auto hit_index_track_iter = hit_index_m.find(key);
226  if(hit_index_track_iter == hit_index_m.end()) {
227  int new_hit_index = this->__GetEdepArray__(real_track_id).size();
228  hit_index_m[key]= new_hit_index;
229  }
230  else {
231  hit_index = (*hit_index_track_iter).second;
232  }
233  auto const pid = geom->ChannelToWire(ch)[0].planeID();
234  auto const channel_id = pindex[pid];
235  double charge = sed.NumElectrons();
236  if(hit_index < 0) {
237  // This particle energy deposition is never recorded so far. Create a new Edep
238  //float charge = ide.numElectrons * detp->ElectronsToADC();
239  this->__GetEdepArray__(real_track_id).emplace_back(pos, pid, pindex.size(), sed.Energy(), charge, channel_id);
240  } else {
241  // Append charge to the relevant edep (@ hit_index)
242  //float charge = ide.numElectrons * detp->ElectronsToADC();
243  MCEdep &edep = this->__GetEdepArray__(real_track_id).at(hit_index);
244  edep.deps[channel_id].charge += charge;
245  edep.deps[channel_id].energy += sed.Energy();
246  }
247  } // end looping over planes
248  }// end looping over SimEnergyDeposits
249 
250  if(_debug_mode) {
251  std::cout<< Form(" Collected %zu particles' energy depositions...",_mc_edeps.size()) << std::endl;
252  // for c2: disable the entire loop instead of just the print statement
253  //for(auto const& track_id_index : _track_index ) {
254  //auto track_id = track_id_index.first;
255  //auto edep_index = track_id_index.second;
256  // std::cout<< Form(" Track ID: %d ... %zu Edep!", track_id, edep_index) << std::endl;
257  //}
258  std::cout<<std::endl;
259  }
260 
261  return;
262  }
263 
264  void MCRecoEdep::MakeMCEdep(const std::vector<sim::SimEnergyDepositLite>& sedArray) {
265  // Create a substitue array of sim::SimEnergyDeposit to avoid duplicating code...
266  // Note that this makes use of the explicit conversion operator
267  // defined in SimEnergyDepositLite class. Information will be partial.
268  // Most notably for MakeMCEdep, charge (numElectrons) will be 0.
269  std::vector<sim::SimEnergyDeposit> new_sedArray(sedArray.begin(), sedArray.end());
270  MakeMCEdep(new_sedArray);
271  }
272 
273 
274 }
void MakeMCEdep(const std::vector< sim::SimChannel > &schArray)
Definition: MCRecoEdep.cxx:64
unsigned int Nplanes() const
Number of planes in this tpc.
Definition: TPCGeo.h:165
pdgs p
Definition: selectors.fcl:22
Geometry information for a single TPC.
Definition: TPCGeo.h:38
std::map< geo::PlaneID, size_t > createPlaneIndexMap()
Definition: MCRecoEdep.cxx:26
const std::vector< sim::MCEdep > & GetEdepArrayAt(size_t edep_index) const
Returns a vector of MCEdep object at the given index.
Definition: MCRecoEdep.cxx:49
tuple m
now if test mode generate materials, CRT shell, world, gdml header else just generate CRT shell for u...
Access the description of detector geometry.
MCRecoEdep(fhicl::ParameterSet const &pset)
Default constructor with fhicl parameters.
Definition: MCRecoEdep.cxx:42
bool ExistTrack(const unsigned int track_id) const
Definition: MCRecoEdep.h:109
std::vector< std::vector< sim::MCEdep > > _mc_edeps
Definition: MCRecoEdep.h:139
std::vector< sim::MCEdep > & __GetEdepArray__(unsigned int track_id)
Definition: MCRecoEdep.cxx:56
do i e
std::vector< deposit > deps
Definition: MCRecoEdep.h:83
std::map< unsigned int, size_t > _track_index
Definition: MCRecoEdep.h:138
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:28
art framework interface to geometry description
BEGIN_PROLOG could also be cout
Encapsulate the construction of a single detector plane.