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