41 #include "art/Framework/Core/EDProducer.h"
42 #include "art/Framework/Core/ModuleMacros.h"
43 #include "art/Framework/Principal/Event.h"
44 #include "art/Framework/Principal/Handle.h"
45 #include "art/Framework/Services/Registry/ServiceHandle.h"
46 #include "art_root_io/TFileDirectory.h"
47 #include "art_root_io/TFileService.h"
48 #include "canvas/Persistency/Common/Assns.h"
49 #include "canvas/Persistency/Common/FindManyP.h"
50 #include "canvas/Persistency/Common/Ptr.h"
51 #include "canvas/Persistency/Common/PtrVector.h"
52 #include "fhiclcpp/ParameterSet.h"
53 #include "messagefacility/MessageLogger/MessageLogger.h"
79 void produce(art::Event&
e)
override;
92 const std::vector<KalmanOutput>&
outputs,
94 std::vector<recob::Track>&
tracks,
95 std::vector<recob::SpacePoint>& spts,
96 art::Assns<recob::Track, recob::Hit>& th_assn,
97 art::Assns<recob::Track, recob::Hit, recob::TrackHitMeta>& thm_assn,
98 art::Assns<recob::Track, recob::SpacePoint>& tsp_assn,
99 art::Assns<recob::SpacePoint, recob::Hit>& sph_assn,
100 art::Assns<recob::PFParticle, recob::Track>& pfPartTrack_assns);
134 , fUseClusterHits(
false)
135 , fUsePFParticleHits(
false)
136 , fUsePFParticleSeeds(
false)
137 , fTKHAlg(pset.get<fhicl::ParameterSet>(
"Track3DKalmanHitAlg"))
138 , fSpacePointAlg(pset.get<fhicl::ParameterSet>(
"SpacePointAlg"))
143 fHist = pset.get<
bool>(
"Hist");
144 fUseClusterHits = pset.get<
bool>(
"UseClusterHits");
145 fUsePFParticleHits = pset.get<
bool>(
"UsePFParticleHits");
146 fUsePFParticleSeeds = pset.get<
bool>(
"UsePFParticleSeeds");
147 fHitModuleLabel = pset.get<std::string>(
"HitModuleLabel");
148 fClusterModuleLabel = pset.get<std::string>(
"ClusterModuleLabel");
149 fPFParticleModuleLabel = pset.get<std::string>(
"PFParticleModuleLabel");
150 if (fUseClusterHits && fUsePFParticleHits) {
151 throw cet::exception(
"Track3DKalmanHit")
152 <<
"Using input from both clustered and PFParticle hits.\n";
155 produces<std::vector<recob::Track>>();
156 produces<std::vector<recob::SpacePoint>>();
159 produces<art::Assns<recob::Track, recob::Hit, recob::TrackHitMeta>>();
160 produces<art::Assns<recob::Track, recob::SpacePoint>>();
161 produces<art::Assns<recob::SpacePoint, recob::Hit>>();
162 produces<art::Assns<recob::PFParticle, recob::Track>>();
166 mf::LogInfo(
"Track3DKalmanHit") <<
"Track3DKalmanHit configured with the following parameters:\n"
167 <<
" UseClusterHits = " << fUseClusterHits <<
"\n"
168 <<
" HitModuleLabel = " << fHitModuleLabel <<
"\n"
169 <<
" ClusterModuleLabel = " << fClusterModuleLabel;
179 art::ServiceHandle<art::TFileService>
tfs;
180 art::TFileDirectory
dir = tfs->mkdir(
"hitkalman",
"Track3DKalmanHit histograms");
182 fHIncChisq = dir.make<TH1F>(
"IncChisq",
"Incremental Chisquare", 100, 0., 20.);
183 fHPull = dir.make<TH1F>(
"Pull",
"Hit Pull", 100, -10., 10.);
194 auto tracks = std::make_unique<std::vector<recob::Track>>();
195 auto th_assn = std::make_unique<art::Assns<recob::Track, recob::Hit>>();
196 auto thm_assn = std::make_unique<art::Assns<recob::Track, recob::Hit, recob::TrackHitMeta>>();
197 auto tsp_assn = std::make_unique<art::Assns<recob::Track, recob::SpacePoint>>();
198 auto pfPartTrack_assns = std::make_unique<art::Assns<recob::PFParticle, recob::Track>>();
199 auto spts = std::make_unique<std::vector<recob::SpacePoint>>();
200 auto sph_assn = std::make_unique<art::Assns<recob::SpacePoint, recob::Hit>>();
203 auto inputs = getInput(evt);
205 auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
207 art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(evt, clockData);
210 if (fHist) { fillHistograms(
outputs); }
224 fSpacePointAlg.clearHitMap();
226 evt.put(std::move(
tracks));
227 evt.put(std::move(spts));
228 evt.put(std::move(th_assn));
229 evt.put(std::move(thm_assn));
230 evt.put(std::move(tsp_assn));
231 evt.put(std::move(sph_assn));
232 evt.put(std::move(pfPartTrack_assns));
240 mf::LogInfo(
"Track3DKalmanHit") <<
"Track3DKalmanHit statistics:\n"
241 <<
" Number of events = " << fNumEvent <<
"\n";
249 fSpacePointAlg.clearHitMap();
261 if (fUsePFParticleHits)
return getPFParticleStuff(evt);
271 art::Handle<std::vector<recob::Cluster>> clusterh;
272 evt.getByLabel(fClusterModuleLabel, clusterh);
274 if (!clusterh.isValid())
return hits;
276 art::FindManyP<recob::Hit> hitsbycluster(clusterh, evt, fClusterModuleLabel);
278 for (
size_t i = 0; i < clusterh->size(); ++i) {
279 std::vector<art::Ptr<recob::Hit>> clushits = hitsbycluster.at(i);
280 hits.insert(hits.end(), clushits.begin(), clushits.end());
291 art::Handle<std::vector<recob::Hit>> hith;
292 evt.getByLabel(fHitModuleLabel, hith);
293 if (!hith.isValid())
return hits;
294 size_t nhits = hith->size();
297 for (
size_t i = 0; i < nhits; ++i) {
298 hits.push_back(art::Ptr<recob::Hit>(hith, i));
314 art::Handle<std::vector<recob::PFParticle>> pfParticleHandle;
315 evt.getByLabel(fPFParticleModuleLabel, pfParticleHandle);
316 if (!pfParticleHandle.isValid())
return inputs;
318 art::Handle<std::vector<recob::Cluster>> clusterHandle;
319 evt.getByLabel(fClusterModuleLabel, clusterHandle);
322 if (!clusterHandle.isValid())
return inputs;
326 art::FindManyP<recob::Cluster> clusterAssns(pfParticleHandle, evt, fPFParticleModuleLabel);
329 art::FindManyP<recob::Seed> seedAssns(pfParticleHandle, evt, fPFParticleModuleLabel);
332 art::FindManyP<recob::Hit> clusterHitAssns(clusterHandle, evt, fClusterModuleLabel);
334 inputs.reserve(pfParticleHandle->size());
337 for (
size_t partIdx = 0; partIdx < pfParticleHandle->size(); partIdx++) {
340 inputs.emplace_back();
342 kalman_input.
pfPartPtr = art::Ptr<recob::PFParticle>(pfParticleHandle, partIdx);
347 std::vector<art::Ptr<recob::Cluster>> clusterVec = clusterAssns.at(partIdx);
349 for (
auto const&
cluster : clusterVec) {
350 std::vector<art::Ptr<recob::Hit>> hitVec = clusterHitAssns.at(
cluster.key());
351 hits.insert(hits.end(), hitVec.begin(), hitVec.end());
355 if (!fUsePFParticleSeeds)
continue;
356 art::PtrVector<recob::Seed>&
seeds = kalman_input.
seeds;
357 std::vector<art::Ptr<recob::Seed>> seedVec = seedAssns.at(partIdx);
358 seeds.insert(seeds.end(), seedVec.begin(), seedVec.end());
359 art::FindManyP<recob::Hit> seedHitAssns(seedVec, evt, fPFParticleModuleLabel);
361 for (
size_t seedIdx = 0; seedIdx < seedVec.size(); ++seedIdx) {
362 std::vector<art::Ptr<recob::Hit>> seedHitVec;
365 seedHitVec = seedHitAssns.at(seedIdx);
367 catch (art::Exception
const&) {
370 kalman_input.
seedhits.emplace_back();
372 seedhits.insert(seedhits.end(), seedHitVec.begin(), seedHitVec.end());
381 const art::Event&
evt,
383 std::vector<KalmanOutput>
const&
outputs,
385 std::vector<recob::Track>&
tracks,
386 std::vector<recob::SpacePoint>& spts,
387 art::Assns<recob::Track, recob::Hit>& th_assn,
388 art::Assns<recob::Track, recob::Hit, recob::TrackHitMeta>& thm_assn,
389 art::Assns<recob::Track, recob::SpacePoint>& tsp_assn,
390 art::Assns<recob::SpacePoint, recob::Hit>& sph_assn,
391 art::Assns<recob::PFParticle, recob::Track>& pfPartTrack_assns)
393 if (outputs.size() != inputs.size())
return;
395 size_t tracksSize(0);
396 for (
auto const& kalman_output : outputs) {
397 tracksSize += kalman_output.tracks.size();
399 tracks.reserve(tracksSize);
401 auto const tid = evt.getProductID<std::vector<recob::Track>>();
402 auto const tidgetter = evt.productGetter(tid);
404 auto const spacepointId = evt.getProductID<std::vector<recob::SpacePoint>>();
405 auto const getter = evt.productGetter(spacepointId);
407 for (
size_t i = 0; i < outputs.size(); ++i) {
409 const std::deque<KGTrack>& kalman_tracks = outputs[i].tracks;
411 for (
auto const& kalman_track : kalman_tracks) {
415 kalman_track.fillTrack(detProp, track, tracks.size());
418 tracks.emplace_back(std::move(track));
423 std::vector<unsigned int> hittpindex;
424 kalman_track.fillHits(trhits, hittpindex);
425 if (hittpindex.back() >= numtrajpts) {
426 throw cet::exception(
"Track3DKalmanHit")
427 <<
"Last hit corresponds to trajectory point index " << hittpindex.back()
428 <<
" while the number of trajectory points is " << numtrajpts <<
'\n';
432 auto nspt = spts.size();
433 fSpacePointAlg.fillSpacePoints(detProp, spts, kalman_track.TrackMap());
435 std::vector<art::Ptr<recob::SpacePoint>> sptvec;
436 for (
auto ispt = nspt; ispt < spts.size(); ++ispt) {
437 sptvec.emplace_back(spacepointId, ispt, getter);
440 auto const& sphits = fSpacePointAlg.getAssociatedHits((spts)[ispt]);
441 for (
auto const& sphit : sphits) {
442 sph_assn.addSingle(sptvec.back(), sphit);
446 art::Ptr<recob::Track> aptr(tid, tracks.size() - 1, tidgetter);
449 for (
size_t h = 0;
h < trhits.size(); ++
h) {
450 th_assn.addSingle(aptr, trhits[
h]);
452 thm_assn.addSingle(aptr, trhits[h], metadata);
456 for (
auto const& spt : sptvec) {
457 tsp_assn.addSingle(aptr, spt);
461 if (fUsePFParticleHits) { pfPartTrack_assns.addSingle(inputs[i].pfPartPtr, aptr); }
473 for (
auto const&
output : outputs) {
474 const std::deque<KGTrack>& kalman_tracks =
output.tracks;
475 for (
size_t i = 0; i < kalman_tracks.size(); ++i) {
476 const KGTrack& trg = kalman_tracks[i];
478 const std::multimap<double, KHitTrack>& trackmap = trg.
getTrackMap();
479 for (std::multimap<double, KHitTrack>::const_iterator ih = trackmap.begin();
480 ih != trackmap.end();
483 const std::shared_ptr<const KHitBase>&
hit = trh.
getHit();
484 double chisq = hit->getChisq();
485 fHIncChisq->Fill(chisq);
bool fUseClusterHits
Use clustered hits as input.
const std::shared_ptr< const KHitBase > & getHit() const
Measurement.
std::string fHitModuleLabel
Unclustered Hits.
Track3DKalmanHitAlg fTKHAlg
Track3DKalmanHit algorithm.
bool fUsePFParticleSeeds
Use PFParticle seeds.
ClusterModuleLabel join with tracks
KalmanInputs getInput(const art::Event &evt) const
void beginJob() override
Begin job method.
Declaration of signal hit object.
std::string fPFParticleModuleLabel
PFParticle label.
size_t NumberTrajectoryPoints() const
Various functions related to the presence and the number of (valid) points.
bool fHist
Make histograms.
std::string fClusterModuleLabel
Clustered Hits.
process_name use argoneut_mc_hitfinder track
Track3DKalmanHit(fhicl::ParameterSet const &pset)
const KVector< N >::type & getResVector() const
Residual vector.
Track3DKalmanHit Algorithm.
SpacePointAlg fSpacePointAlg
Space point algorithm.
const std::multimap< double, KHitTrack > & getTrackMap() const
KHitTrack collection, indexed by path distance.
const KSymMatrix< N >::type & getResError() const
Residual error matrix.
Hits getClusteredHits(const art::Event &evt) const
Fill a collection using clustered hits.
Kalman filter measurement class template.
Hits getAllHits(const art::Event &evt) const
If both UseClusteredHits and UsePFParticles is false use this method to fill in hits.
void fillHistograms(std::vector< KalmanOutput > &outputs)
Fill Histograms method.
Declaration of cluster object.
Provides recob::Track data product.
int fNumEvent
Number of events seen.
std::vector< TrajPoint > seeds
art::PtrVector< recob::Hit > Hits
bool fUsePFParticleHits
Use PFParticle hits as input.
void produce(art::Event &e) override
BEGIN_PROLOG sequence::SlidingWindowTriggerPatternsOppositeWindows END_PROLOG simSlidingORM6O6 effSlidingORW output
TrackCollectionProxyElement< TrackCollProxy > Track
Proxy to an element of a proxy collection of recob::Track objects.
2D representation of charge deposited in the TDC/wire plane
art::ServiceHandle< art::TFileService > tfs
Algorithm for generating space points from hits.
void createOutputs(const art::Event &evt, detinfo::DetectorPropertiesData const &detProp, const std::vector< KalmanOutput > &outputs, const KalmanInputs &inputs, std::vector< recob::Track > &tracks, std::vector< recob::SpacePoint > &spts, art::Assns< recob::Track, recob::Hit > &th_assn, art::Assns< recob::Track, recob::Hit, recob::TrackHitMeta > &thm_assn, art::Assns< recob::Track, recob::SpacePoint > &tsp_assn, art::Assns< recob::SpacePoint, recob::Hit > &sph_assn, art::Assns< recob::PFParticle, recob::Track > &pfPartTrack_assns)
std::vector< KalmanInput > KalmanInputs
void endJob() override
End job method.
TH1F * fHIncChisq
Incremental chisquare.
BEGIN_PROLOG sequence::SlidingWindowTriggerPatternsWindowPair END_PROLOG trigslidewindowOR6m output outputs
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track:
KalmanInputs getPFParticleStuff(const art::Event &evt) const
If UsePFParticles is true use this method to fill in hits.