HitFinderProtoDUNESP_module.cc
Go to the documentation of this file.
1 #ifndef HITFINDERPROTODUNESP_H
2 #define HITFINDERPROTODUNESP_H
3 
4 ////////////////////////////////////////////////////////////////////////
5 //
6 // HitFinderProtoDUNESP class
7 //
8 // trj@fnal.gov
9 //
10 // Runs the disambiguation algorithm for the single-phase ProtoDUNE detector
11 //
12 //
13 //
14 //
15 ////////////////////////////////////////////////////////////////////////
16 
17 #include <string>
18 #include <vector>
19 #include <utility> // std::move()
20 #include <memory> // std::unique_ptr<>
21 
22 // Framework includes
26 
27 
28 // LArSoft Includes
35 #include "DisambigAlgProtoDUNESP.h"
36 
37 // ROOT Includes
38 #include "TH1D.h"
39 #include "TDecompSVD.h"
40 #include "TMath.h"
41 #include "TF1.h"
42 #include "TTree.h"
43 #include "TVectorD.h"
44 #include "TVector2.h"
45 #include "TVector3.h"
46 
47 namespace dune{
49 
50  public:
51 
52  explicit HitFinderProtoDUNESP(fhicl::ParameterSet const& pset);
53 
54  void produce(art::Event& evt);
55  void beginJob();
56  void endJob();
57  void reconfigure(fhicl::ParameterSet const& p);
58 
59 
60  private:
61 
63 
65 
67  std::string fAlg; // which algorithm to use
68 
69  protected:
70 
71 
72  }; // class HitFinderProtoDUNESP
73 
74 
75  //-------------------------------------------------
76  //-------------------------------------------------
78  EDProducer(pset),
79  fDisambigAlg(pset.get< fhicl::ParameterSet >("DisambigAlg"))
80  {
81  this->reconfigure(pset);
82 
83  // let HitCollectionCreator declare that we are going to produce
84  // hits and associations with wires and raw digits
85  // (with no particular product label)
87  }
88 
89 
90  //-------------------------------------------------
91  //-------------------------------------------------
93  {
94 
95  fChanHitLabel = p.get< std::string >("ChanHitLabel");
96  fAlg = p.get < std::string >("Algorithm"); // switch on which algorithm to run
97 
98  }
99 
100  //-------------------------------------------------
101  //-------------------------------------------------
103  {
104 
105  }
106 
107  //-------------------------------------------------
108  //-------------------------------------------------
110  {
111 
112  }
113 
114 
115  //-------------------------------------------------
117  {
118 
119  auto ChannelHits = evt.getValidHandle<std::vector<recob::Hit>>(fChanHitLabel);
120 
121  // also get the associated wires and raw digits;
122  // we assume they have been created by the same module as the hits
123  art::FindOneP<raw::RawDigit> ChannelHitRawDigits
124  (ChannelHits, evt, fChanHitLabel);
125  art::FindOneP<recob::Wire> ChannelHitWires
126  (ChannelHits, evt, fChanHitLabel);
127 
128  // this object contains the hit collection
129  // and its associations to wires and raw digits:
131  /* doWireAssns */ ChannelHitWires.isValid(),
132  /* doRawDigitAssns */ ChannelHitRawDigits.isValid()
133  );
134 
135  // Make unambiguous collection hits
136  std::vector< art::Ptr<recob::Hit> > ChHits;
137  art::fill_ptr_vector(ChHits, ChannelHits);
138  for( size_t h = 0; h < ChHits.size(); h++ ){
139  if( ChHits[h]->View() != geo::kZ ) continue;
140 
141  art::Ptr<recob::Wire> wire = ChannelHitWires.at(h);
142  art::Ptr<raw::RawDigit> rawdigits = ChannelHitRawDigits.at(h);
143 
144  // just copy it
145  hcol.emplace_back(*ChHits[h], wire, rawdigits);
146 
147  } // for
148 
149  // Run alg on all APAs -- check fAlg if we have more than one algorithm defined (future expansion)
150  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(evt);
151  fDisambigAlg.RunDisambig(detProp, ChHits);
152 
153  for( size_t t=0; t < fDisambigAlg.fDisambigHits.size(); t++ ){
155  geo::WireID wid = fDisambigAlg.fDisambigHits[t].second;
156 
157  // create a new hit copy of the original one, but with new wire ID
158  recob::HitCreator disambiguous_hit(*hit, wid);
159 
160  // get the objects associated with the original hit;
161  // since hit comes from ChannelHits, its key is the index in that collection
162  // and also the index for the query of associated objects
163  art::Ptr<recob::Hit>::key_type hit_index = hit.key();
164  art::Ptr<recob::Wire> wire = ChannelHitWires.at(hit_index);
165  art::Ptr<raw::RawDigit> rawdigits = ChannelHitRawDigits.at(hit_index);
166 
167  hcol.emplace_back(disambiguous_hit.move(), wire, rawdigits);
168  } // for
169 
170 
171  // put the hit collection and associations into the event
172  hcol.put_into(evt);
173 
174  }
175 
176 
178 
179 } // end of dune namespace
180 #endif // HITFINDERPROTODUNESP_H
AdcChannelData::View View
std::string string
Definition: nybbler.cc:12
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
Planes which measure Z direction.
Definition: geo_types.h:132
std::size_t key_type
Definition: Ptr.h:79
void RunDisambig(detinfo::DetectorPropertiesData const &detProp, const std::vector< art::Ptr< recob::Hit > > &OrigHits)
Run disambiguation as currently configured.
static void declare_products(art::ProducesCollector &collector, std::string instance_name="", bool doWireAssns=true, bool doRawDigitAssns=true)
Declares the hit products we are going to fill.
Definition: HitCreator.cxx:248
art framework interface to geometry description
art::ServiceHandle< geo::Geometry > fGeom
Class managing the creation of a new recob::Hit object.
Definition: HitCreator.h:83
std::vector< std::pair< art::Ptr< recob::Hit >, geo::WireID > > fDisambigHits
The final list of hits to pass back to be made.
Helper functions to create a hit.
A class handling a collection of hits and its associations.
Definition: HitCreator.h:508
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
key_type key() const noexcept
Definition: Ptr.h:216
T get(std::string const &key) const
Definition: ParameterSet.h:271
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
Definition: DataViewImpl.h:441
p
Definition: test.py:223
Detector simulation of raw signals on wires.
ProducesCollector & producesCollector() noexcept
Declaration of signal hit object.
void reconfigure(fhicl::ParameterSet const &p)
TCEvent evt
Definition: DataStructs.cxx:7
auto const & get(AssnsNode< L, R, D > const &r)
Definition: AssnsNode.h:115
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
HitFinderProtoDUNESP(fhicl::ParameterSet const &pset)
recob::Hit && move()
Prepares the constructed hit to be moved away.
Definition: HitCreator.h:343