9 #include "art/Framework/Core/EDAnalyzer.h"
10 #include "art/Framework/Core/ModuleMacros.h"
11 #include "art/Framework/Principal/Event.h"
12 #include "art/Framework/Principal/Handle.h"
13 #include "art/Framework/Services/Registry/ServiceHandle.h"
14 #include "art_root_io/TFileService.h"
15 #include "canvas/Persistency/Common/FindManyP.h"
16 #include "fhiclcpp/ParameterSet.h"
25 #include "nusimdata/SimulationBase/MCTruth.h"
46 const simb::MCParticle*& MCparticle,
88 , fHitModuleLabel(pset.
get<
std::string>(
"HitModuleLabel",
"trajcluster"))
89 , fShowerModuleLabel(pset.
get<
std::string>(
"ShowerModuleLabel",
"tcshower"))
90 , fGenieGenModuleLabel(pset.
get<
std::string>(
"GenieGenModuleLabel",
"generator"))
91 , fDigitModuleLabel(pset.
get<
std::string>(
"DigitModuleLabel",
"largeant"))
92 , fCalorimetryAlg(pset.
get<fhicl::ParameterSet>(
"CalorimetryAlg"))
101 art::ServiceHandle<art::TFileService const>
tfs;
103 fTree = tfs->make<TTree>(
"tcshowerana",
"tcshowerana");
104 fTree->Branch(
"run", &run,
"run/I");
105 fTree->Branch(
"subrun", &subrun,
"subrun/I");
106 fTree->Branch(
"event", &event,
"event/I");
108 fTree->Branch(
"nuPDG_truth", &nuPDG_truth,
"nuPDG_truth/I");
109 fTree->Branch(
"ccnc_truth", &ccnc_truth,
"ccnc_truth/I");
110 fTree->Branch(
"mode_truth", &mode_truth,
"mode_truth/I");
112 fTree->Branch(
"nshws", &nshws,
"nshws/I");
113 fTree->Branch(
"shwid", shwid,
"shwid[nshws]/I");
114 fTree->Branch(
"shwdcosx", shwdcosx,
"shwdcosx[nshws]/F");
115 fTree->Branch(
"shwdcosy", shwdcosy,
"shwdcosy[nshws]/F");
116 fTree->Branch(
"shwdcosz", shwdcosz,
"shwdcosz[nshws]/F");
117 fTree->Branch(
"shwstartx", shwstartx,
"shwstartx[nshws]/F");
118 fTree->Branch(
"shwstarty", shwstarty,
"shwstarty[nshws]/F");
119 fTree->Branch(
"shwstartz", shwstartz,
"shwstartz[nshws]/F");
120 fTree->Branch(
"shwdedx", shwdedx,
"shwdedx[nshws][2]/D");
121 fTree->Branch(
"shwbestplane", shwbestplane,
"shwbestplane[nshws]/I");
123 fTree->Branch(
"highestHitsPDG", &highestHitsPDG,
"highestHitsPDG/I");
124 fTree->Branch(
"highestHitsFrac", &highestHitsFrac,
"highestHitsFrac/D");
135 subrun = evt.subRun();
136 event = evt.id().event();
138 art::Handle<std::vector<recob::Hit>> hitListHandle;
139 std::vector<art::Ptr<recob::Hit>> hitlist;
140 if (evt.getByLabel(fHitModuleLabel, hitListHandle)) art::fill_ptr_vector(hitlist, hitListHandle);
142 art::Handle<std::vector<sim::SimChannel>> scListHandle;
143 std::vector<art::Ptr<sim::SimChannel>> simchanlist;
144 if (evt.getByLabel(fDigitModuleLabel, scListHandle))
145 art::fill_ptr_vector(simchanlist, scListHandle);
147 art::Handle<std::vector<recob::Shower>> showerListHandle;
148 std::vector<art::Ptr<recob::Shower>> showerlist;
149 if (evt.getByLabel(fShowerModuleLabel, showerListHandle))
150 art::fill_ptr_vector(showerlist, showerListHandle);
152 art::Handle<std::vector<simb::MCTruth>> mctruthListHandle;
153 std::vector<art::Ptr<simb::MCTruth>> mclist;
154 if (evt.getByLabel(fGenieGenModuleLabel, mctruthListHandle))
155 art::fill_ptr_vector(mclist, mctruthListHandle);
157 art::FindManyP<recob::Hit> shwfm(showerListHandle, evt, fShowerModuleLabel);
160 if (showerListHandle.isValid()) {
161 nshws = showerlist.size();
163 for (
int i = 0; i < std::min(
int(showerlist.size()),
kMaxShowers); ++i) {
164 shwid[i] = showerlist[i]->ID();
165 shwdcosx[i] = showerlist[i]->Direction().X();
166 shwdcosy[i] = showerlist[i]->Direction().Y();
167 shwdcosz[i] = showerlist[i]->Direction().Z();
168 shwstartx[i] = showerlist[i]->ShowerStart().X();
169 shwstarty[i] = showerlist[i]->ShowerStart().Y();
170 shwstartz[i] = showerlist[i]->ShowerStart().Z();
171 for (
size_t j = 0; j < (showerlist[i]->dEdx()).
size(); ++j) {
172 shwdedx[i][j] = showerlist[i]->dEdx()[j];
174 shwbestplane[i] = showerlist[i]->best_plane();
179 auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
182 art::Ptr<simb::MCTruth> mctruth = mclist[0];
183 if (mctruth->NeutrinoSet()) {
185 nuPDG_truth = mctruth->GetNeutrino().Nu().PdgCode();
186 ccnc_truth = mctruth->GetNeutrino().CCNC();
187 mode_truth = mctruth->GetNeutrino().Mode();
189 if (showerlist.size()) {
191 std::vector<art::Ptr<recob::Hit>> showerhits = shwfm.at(0);
193 const simb::MCParticle* particle;
194 double tmpEfrac = 0.0;
195 double tmpEcomplet = 0;
196 truthMatcher(clockData, hitlist, showerhits, particle, tmpEfrac, tmpEcomplet);
198 std::cout <<
"shower pdg: " << particle->PdgCode() <<
" efrac " << tmpEfrac << std::endl;
200 highestHitsPDG = particle->PdgCode();
201 highestHitsFrac = tmpEfrac;
221 nuPDG_truth = -99999;
228 shwdcosx[i] = -99999;
229 shwdcosy[i] = -99999;
230 shwdcosz[i] = -99999;
231 shwstartx[i] = -99999;
232 shwstarty[i] = -99999;
233 shwstartz[i] = -99999;
234 for (
int j = 0; j < 2; ++j) {
235 shwdedx[i][j] = -99999;
237 shwbestplane[i] = -99999;
240 highestHitsPDG = -99999;
241 highestHitsFrac = -99999;
251 const simb::MCParticle*& MCparticle,
260 art::ServiceHandle<cheat::BackTrackerService const> bt_serv;
261 art::ServiceHandle<cheat::ParticleInventoryService const> pi_serv;
262 std::map<int, double> trkID_E;
263 for (
size_t j = 0; j < shower_hits.size(); ++j) {
264 art::Ptr<recob::Hit>
hit = shower_hits[j];
266 if (hit->View() != 1)
continue;
267 std::vector<sim::TrackIDE>
TrackIDs = bt_serv->HitToEveTrackIDEs(clockData, hit);
268 for (
size_t k = 0;
k < TrackIDs.size();
k++) {
269 if (trkID_E.find(
std::abs(TrackIDs[
k].trackID)) == trkID_E.end())
270 trkID_E[
std::abs(TrackIDs[k].trackID)] = 0;
271 trkID_E[
std::abs(TrackIDs[k].trackID)] += TrackIDs[
k].energy;
274 double max_E = -999.0;
275 double total_E = 0.0;
277 double partial_E = 0.0;
279 if (!trkID_E.size())
return;
280 for (std::map<int, double>::iterator ii = trkID_E.begin(); ii != trkID_E.end(); ++ii) {
281 total_E += ii->second;
282 if ((ii->second) > max_E) {
283 partial_E = ii->second;
288 MCparticle = pi_serv->TrackIdToParticle_P(TrackID);
290 Efrac = partial_E / total_E;
293 double totenergy = 0;
294 for (
size_t k = 0;
k < all_hits.size(); ++
k) {
295 art::Ptr<recob::Hit>
hit = all_hits[
k];
296 std::vector<sim::TrackIDE>
TrackIDs = bt_serv->HitToEveTrackIDEs(clockData, hit);
297 for (
size_t l = 0; l < TrackIDs.size(); ++l) {
298 if (
std::abs(TrackIDs[l].trackID) ==
TrackID) { totenergy += TrackIDs[l].energy; }
301 Ecomplet = partial_E / totenergy;
double std(const std::vector< short > &wf, const double ped_mean, size_t start, size_t nsample)
float shwstartx[kMaxShowers]
double shwdedx[kMaxShowers][2]
TCShowerAnalysis(fhicl::ParameterSet const &pset)
std::vector< TrackID > TrackIDs
void truthMatcher(detinfo::DetectorClocksData const &clockData, std::vector< art::Ptr< recob::Hit >> all_hits, std::vector< art::Ptr< recob::Hit >> shower_hits, const simb::MCParticle *&MCparticle, double &Efrac, double &Ecomplet)
Declaration of signal hit object.
int shwbestplane[kMaxShowers]
std::size_t size(FixedBins< T, C > const &) noexcept
calo::CalorimetryAlg fCalorimetryAlg
float shwdcosz[kMaxShowers]
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
float shwstarty[kMaxShowers]
constexpr int kMaxShowers
std::string fHitModuleLabel
std::string fDigitModuleLabel
float shwdcosy[kMaxShowers]
std::string fShowerModuleLabel
float shwdcosx[kMaxShowers]
float shwstartz[kMaxShowers]
std::string fGenieGenModuleLabel
void analyze(const art::Event &evt) override
Contains all timing reference information for the detector.
object containing MC truth information necessary for making RawDigits and doing back tracking ...
art::ServiceHandle< art::TFileService > tfs
BEGIN_PROLOG could also be cout