SpacePointSolver_module.cc
Go to the documentation of this file.
1 // Christopher Backhouse - bckhouse@fnal.gov
2 
3 // Test file at Caltech: /nfs/raid11/dunesam/prodgenie_nu_dune10kt_1x2x6_mcc7.0/prodgenie_nu_dune10kt_1x2x6_63_20160811T171439_merged.root
4 
5 // C/C++ standard libraries
6 #include <iostream>
7 #include <string>
8 
9 // framework libraries
18 #include "cetlib/pow.h"
19 #include "fhiclcpp/ParameterSet.h"
20 
21 // LArSoft libraries
23 #include "larcoreobj/SimpleTypesAndConstants/RawTypes.h" // raw::ChannelID_t
27 
30 
32 
33 #include "Solver.h"
34 #include "TripletFinder.h"
35 
36 namespace reco3d {
37 
38  /// The position in the drift direction that an induction wire's time implies
39  /// depends on which TPC we assume the hit originated in.
42 
43  bool
45  {
46  return xpos < w.xpos;
47  }
48 
50  double xpos;
51  };
52 
54  public:
55  explicit SpacePointSolver(const fhicl::ParameterSet& pset);
56 
57  private:
58  void produce(art::Event& evt) override;
59  void beginJob() override;
60 
61  void AddNeighbours(const std::vector<SpaceCharge*>& spaceCharges) const;
62 
63  typedef std::map<const WireHit*, const recob::Hit*> HitMap_t;
64 
65  void BuildSystem(const std::vector<HitTriplet>& triplets,
66  std::vector<CollectionWireHit*>& cwires,
67  std::vector<InductionWireHit*>& iwires,
68  std::vector<SpaceCharge*>& orphanSCs,
69  bool incNei,
70  HitMap_t& hitmap) const;
71 
72  void Minimize(const std::vector<CollectionWireHit*>& cwires,
73  const std::vector<SpaceCharge*>& orphanSCs,
74  double alpha,
75  int maxiterations);
76 
77  /// return whether the point was inserted (only happens when it has charge)
78  bool AddSpacePoint(const SpaceCharge& sc,
79  int id,
81 
82  void FillSystemToSpacePoints(const std::vector<CollectionWireHit*>& cwires,
83  const std::vector<SpaceCharge*>& orphanSCs,
85 
86  void FillSystemToSpacePointsAndAssns(const std::vector<art::Ptr<recob::Hit>>& hitlist,
87  const std::vector<CollectionWireHit*>& cwires,
88  const std::vector<SpaceCharge*>& orphanSCs,
89  const HitMap_t& hitmap,
92 
94 
95  bool fFit;
96  bool fAllowBadInductionHit, fAllowBadCollectionHit;
97 
98  double fAlpha;
99 
100  double fDistThresh;
102 
105 
106  double fXHitOffset;
107 
109  std::unique_ptr<reco3d::IHitReader> fHitReader; ///< Expt specific tool for reading hits
110  };
111 
113 
114  // ---------------------------------------------------------------------------
116  : EDProducer{pset}
117  , fHitLabel(pset.get<std::string>("HitLabel"))
118  , fFit(pset.get<bool>("Fit"))
119  , fAllowBadInductionHit(pset.get<bool>("AllowBadInductionHit"))
120  , fAllowBadCollectionHit(pset.get<bool>("AllowBadCollectionHit"))
121  , fAlpha(pset.get<double>("Alpha"))
122  , fDistThresh(pset.get<double>("WireIntersectThreshold"))
123  , fDistThreshDrift(pset.get<double>("WireIntersectThresholdDriftDir"))
124  , fMaxIterationsNoReg(pset.get<int>("MaxIterationsNoReg"))
125  , fMaxIterationsReg(pset.get<int>("MaxIterationsReg"))
126  , fXHitOffset(pset.get<double>("XHitOffset"))
127  {
128  recob::ChargedSpacePointCollectionCreator::produces(producesCollector(), "pre");
129  if (fFit) {
131  produces<art::Assns<recob::SpacePoint, recob::Hit>>();
132  recob::ChargedSpacePointCollectionCreator::produces(producesCollector(), "noreg");
133  }
134 
135  fHitReader = art::make_tool<reco3d::IHitReader>(pset.get<fhicl::ParameterSet>("HitReaderTool"));
136  }
137 
138  // ---------------------------------------------------------------------------
139  void
141  {
143  }
144 
145  // ---------------------------------------------------------------------------
146  void
147  SpacePointSolver::AddNeighbours(const std::vector<SpaceCharge*>& spaceCharges) const
148  {
149  static const double kCritDist = 5;
150 
151  // Could use a QuadTree or VPTree etc, but seems like overkill
152  class IntCoord {
153  public:
154  IntCoord(const SpaceCharge& sc)
155  : fX(sc.fX / kCritDist), fY(sc.fY / kCritDist), fZ(sc.fZ / kCritDist)
156  {}
157 
158  bool
159  operator<(const IntCoord& i) const
160  {
161  return std::make_tuple(fX, fY, fZ) < std::make_tuple(i.fX, i.fY, i.fZ);
162  }
163 
164  std::vector<IntCoord>
165  Neighbours() const
166  {
167  std::vector<IntCoord> ret;
168  for (int dx = -1; dx <= +1; ++dx) {
169  for (int dy = -1; dy <= +1; ++dy) {
170  for (int dz = -1; dz <= +1; ++dz) {
171  ret.push_back(IntCoord(fX + dx, fY + dy, fZ + dz));
172  }
173  }
174  }
175  return ret;
176  }
177 
178  protected:
179  IntCoord(int x, int y, int z) : fX(x), fY(y), fZ(z) {}
180 
181  int fX, fY, fZ;
182  };
183 
184  std::map<IntCoord, std::vector<SpaceCharge*>> scMap;
185  for (SpaceCharge* sc : spaceCharges) {
186  scMap[IntCoord(*sc)].push_back(sc);
187  }
188 
189  std::cout << "Neighbour search..." << std::endl;
190 
191  // Now that we know all the space charges, can go through and assign neighbours
192 
193  int Ntests = 0;
194  int Nnei = 0;
195  for (SpaceCharge* sc1 : spaceCharges) {
196  IntCoord ic(*sc1);
197  for (IntCoord icn : ic.Neighbours()) {
198  for (SpaceCharge* sc2 : scMap[icn]) {
199 
200  ++Ntests;
201 
202  if (sc1 == sc2) continue;
203  double dist2 =
204  cet::sum_of_squares(sc1->fX - sc2->fX, sc1->fY - sc2->fY, sc1->fZ - sc2->fZ);
205 
206  if (dist2 > cet::square(kCritDist)) continue;
207 
208  if (dist2 == 0) {
209  std::cout << "ZERO DISTANCE SOMEHOW?" << std::endl;
210  std::cout << sc1->fCWire << " " << sc1->fWire1 << " " << sc1->fWire2 << std::endl;
211  std::cout << sc2->fCWire << " " << sc2->fWire1 << " " << sc2->fWire2 << std::endl;
212  std::cout << dist2 << " " << sc1->fX << " " << sc2->fX << " " << sc1->fY << " "
213  << sc2->fY << " " << sc1->fZ << " " << sc2->fZ << std::endl;
214  continue;
215  dist2 = cet::square(kCritDist);
216  }
217 
218  ++Nnei;
219 
220  // This is a pretty random guess
221  const double coupling = exp(-sqrt(dist2) / 2);
222  sc1->fNeighbours.emplace_back(sc2, coupling);
223 
224  if (isnan(1 / sqrt(dist2)) || isinf(1 / sqrt(dist2))) {
225  std::cout << dist2 << " " << sc1->fX << " " << sc2->fX << " " << sc1->fY << " "
226  << sc2->fY << " " << sc1->fZ << " " << sc2->fZ << std::endl;
227  abort();
228  }
229  } // end for sc2
230  } // end for icn
231 
232  // The neighbours lists use the most memory, so be careful to trim
233  sc1->fNeighbours.shrink_to_fit();
234  } // end for sc1
235 
236  for (SpaceCharge* sc : spaceCharges) {
237  for (Neighbour& nei : sc->fNeighbours) {
238  sc->fNeiPotential += nei.fCoupling * nei.fSC->fPred;
239  }
240  }
241 
242  std::cout << Ntests << " tests to find " << Nnei << " neighbours" << std::endl;
243  }
244 
245  // ---------------------------------------------------------------------------
246  void
247  SpacePointSolver::BuildSystem(const std::vector<HitTriplet>& triplets,
248  std::vector<CollectionWireHit*>& cwires,
249  std::vector<InductionWireHit*>& iwires,
250  std::vector<SpaceCharge*>& orphanSCs,
251  bool incNei,
252  HitMap_t& hitmap) const
253  {
254  std::set<const recob::Hit*> ihits;
255  std::set<const recob::Hit*> chits;
256  for (const HitTriplet& trip : triplets) {
257  if (trip.x) chits.insert(trip.x);
258  if (trip.u) ihits.insert(trip.u);
259  if (trip.v) ihits.insert(trip.v);
260  }
261 
262  std::map<const recob::Hit*, InductionWireHit*> inductionMap;
263  for (const recob::Hit* hit : ihits) {
264  InductionWireHit* iwire = new InductionWireHit(hit->Channel(), hit->Integral());
265  inductionMap[hit] = iwire;
266  iwires.emplace_back(iwire);
267  hitmap[iwire] = hit;
268  }
269 
270  std::map<const recob::Hit*, std::vector<SpaceCharge*>> collectionMap;
271  std::map<const recob::Hit*, std::vector<SpaceCharge*>> collectionMapBad;
272 
273  std::set<InductionWireHit*> satisfiedInduction;
274 
275  for (const HitTriplet& trip : triplets) {
276  // Don't have a cwire object yet, set it later
277  SpaceCharge* sc = new SpaceCharge(
278  trip.pt.x, trip.pt.y, trip.pt.z, 0, inductionMap[trip.u], inductionMap[trip.v]);
279 
280  if (trip.u && trip.v) {
281  collectionMap[trip.x].push_back(sc);
282  if (trip.x) {
283  satisfiedInduction.insert(inductionMap[trip.u]);
284  satisfiedInduction.insert(inductionMap[trip.v]);
285  }
286  }
287  else {
288  collectionMapBad[trip.x].push_back(sc);
289  }
290  }
291 
292  std::vector<SpaceCharge*> spaceCharges;
293 
294  for (const recob::Hit* hit : chits) {
295  // Find the space charges associated with this hit
296  std::vector<SpaceCharge*>& scs = collectionMap[hit];
297  if (scs.empty()) {
298  // If there are no full triplets try the triplets with one bad channel
299  scs = collectionMapBad[hit];
300  }
301  else {
302  // If there were good triplets, delete the bad hit ones
303  for (SpaceCharge* sc : collectionMapBad[hit])
304  delete sc;
305  }
306  // Still no space points, don't bother making a wire
307  if (scs.empty()) continue;
308 
309  CollectionWireHit* cwire = new CollectionWireHit(hit->Channel(), hit->Integral(), scs);
310  hitmap[cwire] = hit;
311  cwires.push_back(cwire);
312  spaceCharges.insert(spaceCharges.end(), scs.begin(), scs.end());
313  for (SpaceCharge* sc : scs)
314  sc->fCWire = cwire;
315  } // end for hit
316 
317  // Space charges whose collection wire is bad, which we have no other way of
318  // addressing.
319  for (SpaceCharge* sc : collectionMap[0]) {
320  // Only count orphans where an induction wire has no other explanation
321  if (satisfiedInduction.count(sc->fWire1) == 0 || satisfiedInduction.count(sc->fWire2) == 0) {
322  orphanSCs.push_back(sc);
323  }
324  else {
325  delete sc;
326  }
327  }
328  spaceCharges.insert(spaceCharges.end(), orphanSCs.begin(), orphanSCs.end());
329 
330  std::cout << cwires.size() << " collection wire objects" << std::endl;
331  std::cout << spaceCharges.size() << " potential space points" << std::endl;
332 
333  if (incNei) AddNeighbours(spaceCharges);
334  }
335 
336  // ---------------------------------------------------------------------------
337  bool
339  int id,
341  {
342  static const double err[6] = {
343  0,
344  };
345 
346  const float charge = sc.fPred;
347  if (charge == 0) return false;
348 
349  const double xyz[3] = {sc.fX, sc.fY, sc.fZ};
350  points.add({xyz, err, 0.0, id}, charge);
351 
352  return true;
353  }
354 
355  // ---------------------------------------------------------------------------
356  void
357  SpacePointSolver::FillSystemToSpacePoints(const std::vector<CollectionWireHit*>& cwires,
358  const std::vector<SpaceCharge*>& orphanSCs,
360  {
361  int iPoint = 0;
362  for (const CollectionWireHit* cwire : cwires) {
363  for (const SpaceCharge* sc : cwire->fCrossings) {
364  AddSpacePoint(*sc, iPoint++, points);
365  } // for sc
366  } // for cwire
367 
368  for (const SpaceCharge* sc : orphanSCs)
369  AddSpacePoint(*sc, iPoint++, points);
370  }
371 
372  // ---------------------------------------------------------------------------
373  void
375  const std::vector<art::Ptr<recob::Hit>>& hitlist,
376  const std::vector<CollectionWireHit*>& cwires,
377  const std::vector<SpaceCharge*>& orphanSCs,
378  const HitMap_t& hitmap,
381  {
382  std::map<const recob::Hit*, art::Ptr<recob::Hit>> ptrmap;
383  for (art::Ptr<recob::Hit> hit : hitlist)
384  ptrmap[hit.get()] = hit;
385 
386  std::vector<const SpaceCharge*> scs;
387  for (const SpaceCharge* sc : orphanSCs)
388  scs.push_back(sc);
389  for (const CollectionWireHit* cwire : cwires)
390  for (const SpaceCharge* sc : cwire->fCrossings)
391  scs.push_back(sc);
392 
393  int iPoint = 0;
394 
395  for (const SpaceCharge* sc : scs) {
396  if (!AddSpacePoint(*sc, iPoint++, points)) continue;
397  const auto& spsPtr = points.lastSpacePointPtr();
398 
399  if (sc->fCWire) { assn.addSingle(spsPtr, ptrmap[hitmap.at(sc->fCWire)]); }
400  if (sc->fWire1) { assn.addSingle(spsPtr, ptrmap[hitmap.at(sc->fWire1)]); }
401  if (sc->fWire2) { assn.addSingle(spsPtr, ptrmap[hitmap.at(sc->fWire2)]); }
402  }
403  }
404 
405  // ---------------------------------------------------------------------------
406  void
407  SpacePointSolver::Minimize(const std::vector<CollectionWireHit*>& cwires,
408  const std::vector<SpaceCharge*>& orphanSCs,
409  double alpha,
410  int maxiterations)
411  {
412  double prevMetric = Metric(cwires, alpha);
413  std::cout << "Begin: " << prevMetric << std::endl;
414  for (int i = 0; i < maxiterations; ++i) {
415  Iterate(cwires, orphanSCs, alpha);
416  const double metric = Metric(cwires, alpha);
417  std::cout << i << " " << metric << std::endl;
418  if (metric > prevMetric) {
419  std::cout << "Warning: metric increased" << std::endl;
420  return;
421  }
422  if (fabs(metric - prevMetric) < 1e-3 * fabs(prevMetric)) return;
423  prevMetric = metric;
424  }
425  }
426 
427  // ---------------------------------------------------------------------------
428  void
430  {
432  std::vector<art::Ptr<recob::Hit>> hitlist;
433  if (evt.getByLabel(fHitLabel, hits)) art::fill_ptr_vector(hitlist, hits);
434 
435  auto spcol_pre = recob::ChargedSpacePointCollectionCreator::forPtrs(evt, "pre");
436  auto spcol_noreg = recob::ChargedSpacePointCollectionCreator::forPtrs(evt, "noreg");
438  auto assns = std::make_unique<art::Assns<recob::SpacePoint, recob::Hit>>();
439 
440  // Skip very small events
441  if (hits->size() < 20) {
442  spcol_pre.put();
443  if (fFit) {
444  spcol.put();
445  evt.put(std::move(assns));
446  spcol_noreg.put();
447  }
448  return;
449  }
450 
452 
453  std::vector<art::Ptr<recob::Hit>> xhits, uhits, vhits;
454  bool is2view = fHitReader->readHits(hitlist, xhits, uhits, vhits);
455 
456  std::vector<raw::ChannelID_t> xbadchans, ubadchans, vbadchans;
457  if (fAllowBadInductionHit || fAllowBadCollectionHit) {
458  for (raw::ChannelID_t cid :
459  art::ServiceHandle<lariov::ChannelStatusService const>()->GetProvider().BadChannels()) {
460  if (geom->SignalType(cid) == geo::kCollection) {
461  if (fAllowBadCollectionHit && geom->View(cid) == geo::kZ) { xbadchans.push_back(cid); }
462  }
463  else {
464  if (fAllowBadInductionHit) {
465  if (geom->View(cid) == geo::kU) ubadchans.push_back(cid);
466  if (geom->View(cid) == geo::kV) vbadchans.push_back(cid);
467  }
468  }
469  }
470  }
471  std::cout << xbadchans.size() << " X, " << ubadchans.size() << " U, " << vbadchans.size()
472  << " V bad channels" << std::endl;
473 
474  std::vector<CollectionWireHit*> cwires;
475  // So we can find them all to free the memory
476  std::vector<InductionWireHit*> iwires;
477  // Nodes with a bad collection wire that we otherwise can't address
478  std::vector<SpaceCharge*> orphanSCs;
479 
480  auto const detProp =
482 
483  HitMap_t hitmap;
484  if (is2view) {
485  std::cout << "Finding 2-view coincidences..." << std::endl;
486  TripletFinder tf(detProp,
487  xhits,
488  uhits,
489  {},
490  xbadchans,
491  ubadchans,
492  {},
493  fDistThresh,
494  fDistThreshDrift,
495  fXHitOffset);
496  BuildSystem(tf.TripletsTwoView(), cwires, iwires, orphanSCs, fAlpha != 0, hitmap);
497  }
498  else {
499  std::cout << "Finding XUV coincidences..." << std::endl;
500  TripletFinder tf(detProp,
501  xhits,
502  uhits,
503  vhits,
504  xbadchans,
505  ubadchans,
506  vbadchans,
507  fDistThresh,
508  fDistThreshDrift,
509  fXHitOffset);
510  BuildSystem(tf.Triplets(), cwires, iwires, orphanSCs, fAlpha != 0, hitmap);
511  }
512 
513  FillSystemToSpacePoints(cwires, orphanSCs, spcol_pre);
514  spcol_pre.put();
515 
516  if (fFit) {
517  std::cout << "Iterating with no regularization..." << std::endl;
518  Minimize(cwires, orphanSCs, 0, fMaxIterationsNoReg);
519 
520  FillSystemToSpacePoints(cwires, orphanSCs, spcol_noreg);
521  spcol_noreg.put();
522 
523  std::cout << "Now with regularization..." << std::endl;
524  Minimize(cwires, orphanSCs, fAlpha, fMaxIterationsReg);
525 
526  FillSystemToSpacePointsAndAssns(hitlist, cwires, orphanSCs, hitmap, spcol, *assns);
527  spcol.put();
528  evt.put(std::move(assns));
529  } // end if fFit
530 
531  for (InductionWireHit* i : iwires)
532  delete i;
533  for (CollectionWireHit* c : cwires)
534  delete c;
535  for (SpaceCharge* s : orphanSCs)
536  delete s;
537  }
538 
539 } // end namespace reco3d
const geo::GeometryCore * geom
std::string string
Definition: nybbler.cc:12
Planes which measure V.
Definition: geo_types.h:130
art::Ptr< recob::SpacePoint > lastSpacePointPtr() const
Returns an art pointer to the last inserted space point (no check!).
void Minimize(const std::vector< CollectionWireHit * > &cwires, const std::vector< SpaceCharge * > &orphanSCs, double alpha, int maxiterations)
constexpr T sum_of_squares(T x, T y)
Definition: pow.h:139
static void produces(art::ProducesCollector &producesCollector, std::string const &instanceName={})
Declares the data products being produced.
double fCoupling
Definition: Solver.h:37
struct vector vector
This provides an art tool interface definition for reading hits into the SpacePointSolver universe...
std::unique_ptr< reco3d::IHitReader > fHitReader
Expt specific tool for reading hits.
virtual const provider_type * provider() const override
void AddNeighbours(const std::vector< SpaceCharge * > &spaceCharges) const
Planes which measure Z direction.
Definition: geo_types.h:132
InductionWireWithXPos(InductionWireHit *w, double x)
SigType_t SignalType(geo::PlaneID const &pid) const
Returns the type of signal on the channels of specified TPC plane.
Definition: tf_graph.h:23
constexpr T square(T x)
Definition: pow.h:21
art framework interface to geometry description
void produce(art::Event &evt) override
Planes which measure U.
Definition: geo_types.h:129
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Definition: DataViewImpl.h:633
double fZ
Definition: Solver.h:52
const double e
std::map< const WireHit *, const recob::Hit * > HitMap_t
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
void beginJob()
Definition: Breakpoints.cc:14
void Iterate(CollectionWireHit *cwire, double alpha)
Definition: Solver.cxx:263
InductionWireHit * fWire2
Definition: Solver.h:54
double fY
Definition: Solver.h:52
def move(depos, offset)
Definition: depos.py:107
bool AddSpacePoint(const SpaceCharge &sc, int id, recob::ChargedSpacePointCollectionCreator &points) const
return whether the point was inserted (only happens when it has charge)
double alpha
Definition: doAna.cpp:15
InductionWireHit * fWire1
Definition: Solver.h:54
void FillSystemToSpacePointsAndAssns(const std::vector< art::Ptr< recob::Hit >> &hitlist, const std::vector< CollectionWireHit * > &cwires, const std::vector< SpaceCharge * > &orphanSCs, const HitMap_t &hitmap, recob::ChargedSpacePointCollectionCreator &points, art::Assns< recob::SpacePoint, recob::Hit > &assn) const
double fX
Definition: Solver.h:52
double Metric(double q, double p)
Definition: Solver.cxx:71
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
bool operator<(const InductionWireWithXPos &w) const
Description of geometry of one entire detector.
View_t View(geo::PlaneID const &pid) const
Returns the view (wire orientation) on the channels of specified TPC plane.
void err(const char *fmt,...)
Definition: message.cpp:226
Detector simulation of raw signals on wires.
static ChargedSpacePointCollectionCreator forPtrs(art::Event &event, std::string const &instanceName={})
Static function binding a new object to a specific art event.
void BuildSystem(const std::vector< HitTriplet > &triplets, std::vector< CollectionWireHit * > &cwires, std::vector< InductionWireHit * > &iwires, std::vector< SpaceCharge * > &orphanSCs, bool incNei, HitMap_t &hitmap) const
Creates a collection of space points with associated charge.
CollectionWireHit * fCWire
Definition: Solver.h:53
void addSingle(Ptr< left_t > const &left, Ptr< right_t > const &right, data_t const &data)
Definition: Assns.h:546
std::vector< HitTriplet > Triplets()
void FillSystemToSpacePoints(const std::vector< CollectionWireHit * > &cwires, const std::vector< SpaceCharge * > &orphanSCs, recob::ChargedSpacePointCollectionCreator &pts) const
Interface for experiment-specific channel quality info provider.
list x
Definition: train.py:276
Helpers to create space points with associated charge.
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:48
TCEvent evt
Definition: DataStructs.cxx:7
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:28
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
Interface for experiment-specific service for channel quality info.
double fPred
Definition: Solver.h:58
void add(recob::SpacePoint const &spacePoint, recob::PointCharge const &charge)
Inserts the specified space point and associated data into the collection.
static QCString * s
Definition: config.cpp:1042
QTextStream & endl(QTextStream &s)
Signal from collection planes.
Definition: geo_types.h:146
SpaceCharge * fSC
Definition: Solver.h:36