24 #include "art/Framework/Core/EDProducer.h"
25 #include "art/Framework/Core/ModuleMacros.h"
26 #include "art/Framework/Principal/Event.h"
27 #include "art_root_io/TFileService.h"
28 #include "canvas/Persistency/Common/FindManyP.h"
29 #include "cetlib_except/exception.h"
30 #include "messagefacility/MessageLogger/MessageLogger.h"
49 #include "nusimdata/SimulationBase/MCParticle.h"
55 accepted_particle(
int apdg)
71 void produce(art::Event&
e)
override;
96 DEFINE_ART_MODULE(TrackKalmanCheater)
109 , fHist(pset.get<
bool>(
"Hist"))
110 , fKFAlg(pset.get<fhicl::ParameterSet>(
"KalmanFilterAlg"))
111 , fSpacePointAlg(pset.get<fhicl::ParameterSet>(
"SpacePointAlg"))
112 , fUseClusterHits(pset.get<
bool>(
"UseClusterHits"))
113 , fHitModuleLabel{pset.get<std::string>(
"HitModuleLabel")}
114 , fClusterModuleLabel{pset.get<std::string>(
"ClusterModuleLabel")}
115 , fMaxTcut(pset.get<
double>(
"MaxTcut"))
117 produces<std::vector<recob::Track>>();
118 produces<std::vector<recob::SpacePoint>>();
119 produces<art::Assns<recob::Track, recob::Hit>>();
120 produces<art::Assns<recob::Track, recob::SpacePoint>>();
121 produces<art::Assns<recob::SpacePoint, recob::Hit>>();
125 mf::LogInfo(
"TrackKalmanCheater")
126 <<
"TrackKalmanCheater configured with the following parameters:\n"
127 <<
" UseClusterHits = " << fUseClusterHits <<
"\n"
128 <<
" HitModuleLabel = " << fHitModuleLabel <<
"\n"
129 <<
" ClusterModuleLabel = " << fClusterModuleLabel;
141 art::ServiceHandle<art::TFileService>
tfs;
142 art::TFileDirectory
dir = tfs->mkdir(
"hitkalman",
"TrackKalmanCheater histograms");
144 fHIncChisq = dir.make<TH1F>(
"IncChisq",
"Incremental Chisquare", 100, 0., 20.);
145 fHPull = dir.make<TH1F>(
"Pull",
"Hit Pull", 100, -10., 10.);
167 std::unique_ptr<std::vector<recob::Track>>
tracks(
new std::vector<recob::Track>);
168 std::unique_ptr<art::Assns<recob::Track, recob::Hit>> th_assn(
169 new art::Assns<recob::Track, recob::Hit>);
170 std::unique_ptr<art::Assns<recob::Track, recob::SpacePoint>> tsp_assn(
171 new art::Assns<recob::Track, recob::SpacePoint>);
176 std::unique_ptr<std::vector<recob::SpacePoint>> spts(
new std::vector<recob::SpacePoint>);
177 std::unique_ptr<art::Assns<recob::SpacePoint, recob::Hit>> sph_assn(
178 new art::Assns<recob::SpacePoint, recob::Hit>);
182 std::deque<KGTrack> kalman_tracks;
186 art::ServiceHandle<geo::Geometry const> geom;
187 art::ServiceHandle<cheat::BackTrackerService const> bt_serv;
188 art::ServiceHandle<cheat::ParticleInventoryService const> pi_serv;
189 auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
193 fSpacePointAlg.clearHitMap();
197 art::PtrVector<recob::Hit> hits;
199 if (fUseClusterHits) {
203 art::Handle<std::vector<recob::Cluster>> clusterh;
204 evt.getByLabel(fClusterModuleLabel, clusterh);
207 art::FindManyP<recob::Hit> fm(clusterh, evt, fClusterModuleLabel);
209 if (clusterh.isValid()) {
210 int nclus = clusterh->size();
212 for (
int i = 0; i < nclus; ++i) {
213 art::Ptr<recob::Cluster> pclus(clusterh, i);
214 std::vector<art::Ptr<recob::Hit>> clushits = fm.at(i);
215 int nhits = clushits.size();
216 hits.reserve(hits.size() + nhits);
218 for (
std::vector<art::Ptr<recob::Hit>>::const_iterator ihit = clushits.begin();
219 ihit != clushits.end();
221 hits.push_back(*ihit);
230 art::Handle<std::vector<recob::Hit>> hith;
231 evt.getByLabel(fHitModuleLabel, hith);
232 if (hith.isValid()) {
233 int nhits = hith->size();
236 for (
int i = 0; i < nhits; ++i)
237 hits.push_back(art::Ptr<recob::Hit>(hith, i));
243 std::map<int, art::PtrVector<recob::Hit>> hitmap;
247 for (art::PtrVector<recob::Hit>::const_iterator ihit = hits.begin(); ihit != hits.end(); ++ihit) {
252 std::vector<sim::TrackIDE> tids = bt_serv->HitToTrackIDEs(clockData, *ihit);
256 for (std::vector<sim::TrackIDE>::const_iterator itid = tids.begin(); itid != tids.end();
258 int trackID = itid->trackID;
262 hitmap[trackID].push_back(*ihit);
268 sim::ParticleList
const& plist = pi_serv->ParticleList();
270 auto const det_prop =
271 art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(evt);
280 mf::LogDebug log(
"TrackKalmanCheater");
281 for (sim::ParticleList::const_iterator ipart = plist.begin(); ipart != plist.end(); ++ipart) {
282 const simb::MCParticle* part = (*ipart).second;
283 if (!part)
throw cet::exception(
"TrackKalmanCheater") <<
"no particle!\n";
284 int pdg = part->PdgCode();
289 if (not accepted_particle(apdg)) {
continue; }
291 int trackid = part->TrackId();
292 int stat = part->StatusCode();
294 if (hitmap.count(trackid) != 0) nhit = hitmap[trackid].
size();
295 const TLorentzVector& pos = part->Position();
296 const TLorentzVector& mom = part->Momentum();
298 if (mf::isDebugEnabled()) {
299 log <<
"Trackid=" << trackid <<
", pdgid=" << pdg <<
", status=" << stat
300 <<
", NHit=" << nhit <<
"\n"
301 <<
" x = " << pos.X() <<
", y = " << pos.Y() <<
", z = " << pos.Z() <<
"\n"
302 <<
" px= " << mom.Px() <<
", py= " << mom.Py() <<
", pz= " << mom.Pz() <<
"\n";
308 double px = mom.Px();
309 double py = mom.Py();
310 double pz = mom.Pz();
311 double p = std::hypot(px, py, pz);
313 if (nhit > 0 && pz != 0.) {
314 const art::PtrVector<recob::Hit>& trackhits = hitmap[trackid];
315 if (trackhits.empty())
throw cet::exception(
"TrackKalmanCheater") <<
"No hits in track\n";
319 std::shared_ptr<const Surface> psurf(
new SurfYZPlane(0., y, z, 0.));
327 KTrack trk(psurf, vec, dir, pdg);
332 cont.
fill(det_prop, trackhits, -1);
337 std::vector<unsigned int> planehits(3, 0);
338 for (art::PtrVector<recob::Hit>::const_iterator ih = trackhits.begin();
339 ih != trackhits.end();
344 if (plane >= planehits.size()) {
345 throw cet::exception(
"TrackKalmanCheater") <<
"plane " << plane <<
"...\n";
349 unsigned int prefplane = 0;
350 for (
unsigned int i = 0; i < planehits.size(); ++i) {
351 if (planehits[i] > planehits[prefplane]) prefplane = i;
353 if (mf::isDebugEnabled()) log <<
"Preferred plane = " << prefplane <<
"\n";
358 fKFAlg.setPlane(prefplane);
360 if (fKFAlg.smoothTrackIter(5, trg, prop)) {
362 if (fKFAlg.fitMomentum(trg, prop, tremom)) { fKFAlg.updateMomentum(tremom, prop, trg); }
364 kalman_tracks.push_back(trg);
377 for (
auto const& trg : kalman_tracks) {
381 const std::multimap<double, KHitTrack>& trackmap = trg.getTrackMap();
382 for (std::multimap<double, KHitTrack>::const_iterator ih = trackmap.begin();
383 ih != trackmap.end();
386 const std::shared_ptr<const KHitBase>&
hit = trh.
getHit();
387 double chisq = hit->getChisq();
388 fHIncChisq->Fill(chisq);
400 tracks->reserve(kalman_tracks.size());
401 for (
auto const& kalman_track : kalman_tracks) {
403 kalman_track.fillTrack(det_prop, tracks->back(), tracks->size() - 1);
407 art::PtrVector<recob::Hit> trhits;
408 std::vector<unsigned int> hittpindex;
409 kalman_track.fillHits(hits, hittpindex);
414 int nspt = spts->size();
415 fSpacePointAlg.fillSpacePoints(det_prop, *spts, kalman_track.TrackMap());
420 art::PtrVector<recob::SpacePoint> sptvec;
424 for (
unsigned int ispt = nspt; ispt < spts->size(); ++ispt) {
426 art::ProductID sptid = evt.getProductID<std::vector<recob::SpacePoint>>();
427 art::Ptr<recob::SpacePoint> sptptr(sptid, ispt, evt.productGetter(sptid));
428 sptvec.push_back(sptptr);
432 const art::PtrVector<recob::Hit>& sphits = fSpacePointAlg.getAssociatedHits(spt);
443 evt.put(std::move(tracks));
444 evt.put(std::move(spts));
445 evt.put(std::move(th_assn));
446 evt.put(std::move(tsp_assn));
447 evt.put(std::move(sph_assn));
455 mf::LogInfo(
"TrackKalmanCheater") <<
"TrackKalmanCheater statistics:\n"
456 <<
" Number of events = " << fNumEvent <<
"\n"
457 <<
" Number of tracks created = " << fNumTrack;
void endJob() override
End job method.
process_name opflash particleana ie ie ie z
TrackDirection
Track direction enum.
Basic Kalman filter track class, plus one measurement on same surface.
const std::shared_ptr< const KHitBase > & getHit() const
Measurement.
Propagate to PropYZPlane surface.
void fill(const detinfo::DetectorPropertiesData &clock_data, const art::PtrVector< recob::Hit > &hits, int only_plane) override
ClusterModuleLabel join with tracks
A KHitContainer for KHitWireX type measurements.
process_name opflash particleana ie x
std::string fClusterModuleLabel
Clustered Hits.
geo::WireID WireID() const
Planar surface parallel to x-axis.
Declaration of signal hit object.
std::size_t size(FixedBins< T, C > const &) noexcept
KalmanFilterAlg fKFAlg
Kalman filter algorithm.
const KVector< N >::type & getResVector() const
Residual vector.
int fNumEvent
Number of events seen.
std::string fHitModuleLabel
Unclustered Hits.
void produce(art::Event &e) override
bool fUseClusterHits
Use cluster hits or all hits?
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
TrackKalmanCheater(fhicl::ParameterSet const &pset)
const KSymMatrix< N >::type & getResError() const
Residual error matrix.
process_name opflash particleana ie ie y
Kalman filter measurement class template.
void beginJob() override
Begin job method.
KVector< 5 >::type TrackVector
Track state vector, dimension 5.
PlaneID_t Plane
Index of the plane within its TPC.
Declaration of cluster object.
Provides recob::Track data product.
bool CreateAssn(art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t index=UINT_MAX)
Creates a single one-to-one association.
int fNumTrack
Number of tracks produced.
process_name physics producers generator physics producers generator physics producers generator py
bool fHist
Make histograms.
SpacePointAlg fSpacePointAlg
Space point algorithm.
TH1F * fHIncChisq
Incremental chisquare.
A collection of KHitTracks.
2D representation of charge deposited in the TDC/wire plane
art::ServiceHandle< art::TFileService > tfs
double fMaxTcut
Maximum delta ray energy in MeV for restricted dE/dx.
Algorithm for generating space points from hits.
art framework interface to geometry description
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track: