All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
shwf::ShowerReco Class Reference
Inheritance diagram for shwf::ShowerReco:

Public Member Functions

 ShowerReco (fhicl::ParameterSet const &pset)
 

Private Member Functions

void beginJob ()
 
void beginRun (art::Run &run)
 
void produce (art::Event &evt)
 
void GetVertexAndAnglesFromCluster (art::Ptr< recob::Cluster > clust, unsigned int plane)
 
void LongTransEnergy (geo::GeometryCore const *geom, detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, unsigned int set, std::vector< art::Ptr< recob::Hit >> hitlist)
 
void ClearandResizeVectors (unsigned int nPlanes)
 

Private Attributes

int fRun
 
int fEvent
 
int fSubRun
 
float slope [3]
 
float angle [3]
 
std::string fClusterModuleLabel
 
float ftimetick
 
double fMean_wire_pitch
 
fhicl::ParameterSet fCaloPSet
 
std::vector< double > fRMS_2cm
 
std::vector< int > fNpoints_2cm
 
std::vector< double > fCorr_MeV_2cm
 
std::vector< double > fCorr_Charge_2cm
 
std::vector< int > fNpoints_corr_ADC_2cm
 
std::vector< int > fNpoints_corr_MeV_2cm
 
std::vector< double > fTotChargeADC
 
std::vector< double > fTotChargeMeV
 
std::vector< double > fTotChargeMeV_MIPs
 
std::vector< double > fChargeADC_2cm
 
std::vector< double > fChargeMeV_2cm
 
std::vector< double > fChargeMeV_2cm_refined
 
std::vector< double > fChargeMeV_2cm_axsum
 
std::vector< std::vector
< double > > 
fDistribChargeADC
 
std::vector< std::vector
< double > > 
fDistribChargeMeV
 
std::vector< std::vector
< double > > 
fDistribHalfChargeMeV
 
std::vector< std::vector
< double > > 
fDistribChargeposition
 
std::vector< std::vector
< double > > 
fSingleEvtAngle
 
std::vector< std::vector
< double > > 
fSingleEvtAngleVal
 
std::vector< unsigned int > fWire_vertex
 
std::vector< double > fTime_vertex
 
std::vector< double > fWire_vertexError
 
std::vector< double > fTime_vertexError
 
std::vector< unsigned int > fWire_last
 
std::vector< double > fTime_last
 
std::vector< double > xyz_vertex_fit
 
std::vector< std::vector
< double > > 
fNPitch
 
float Kin_En
 
std::vector< float > vdEdx
 
std::vector< float > vresRange
 
std::vector< float > vdQdx
 
std::vector< float > deadwire
 
float Trk_Length
 
float fTrkPitchC
 
float fdEdxlength
 
float fcalodEdxlength
 
bool fUseArea
 
double xphi
 
double xtheta
 
unsigned int fNPlanes
 
unsigned int fNAngles
 
TTree * ftree_shwf
 
double fWirePitch
 
double fTimeTick
 
double fDriftVelocity
 
double fWireTimetoCmCm
 
std::vector< int > fNhitsperplane
 
std::vector< double > fTotADCperplane
 

Detailed Description

Definition at line 58 of file ShowerReco_module.cc.

Constructor & Destructor Documentation

shwf::ShowerReco::ShowerReco ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 162 of file ShowerReco_module.cc.

162  : EDProducer{pset}
163  {
164  fClusterModuleLabel = pset.get<std::string>("ClusterModuleLabel");
165  fCaloPSet = pset.get<fhicl::ParameterSet>("CaloAlg");
166 
167  fdEdxlength =
168  pset.get<double>("dEdxlength"); // distance that gets used to determine e/gamma separation
170  pset.get<double>("calodEdxlength"); // cutoff distance for hits saved to the calo object.
171  fUseArea = pset.get<bool>("UseArea");
172 
173  produces<std::vector<recob::Shower>>();
174  produces<art::Assns<recob::Shower, recob::Cluster>>();
175  produces<art::Assns<recob::Shower, recob::Hit>>();
176  produces<std::vector<anab::Calorimetry>>();
177  produces<art::Assns<recob::Shower, anab::Calorimetry>>();
178  }
fhicl::ParameterSet fCaloPSet
std::string fClusterModuleLabel

Member Function Documentation

void shwf::ShowerReco::beginJob ( )
private
Todo:
the call to geo->Nplanes() assumes this is a single cryostat and single TPC detector; need to generalize to multiple cryostats and TPCs

Get TFileService and define output Histograms

All-knowing tree with reconstruction information

Definition at line 182 of file ShowerReco_module.cc.

183  {
184  art::ServiceHandle<geo::Geometry const> geo;
185 
186  /// \todo the call to geo->Nplanes() assumes this is a single cryostat and
187  /// single TPC detector; need to generalize to multiple cryostats and
188  /// TPCs
189  fNPlanes = geo->Nplanes();
190  fMean_wire_pitch = geo->WirePitch(); // wire pitch in cm
191 
192  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataForJob();
193  ftimetick = sampling_rate(clockData) / 1000.;
194 
195  /**Get TFileService and define output Histograms*/
196  art::ServiceHandle<art::TFileService const> tfs;
197 
198  ftree_shwf = tfs->make<TTree>("ShowerReco",
199  "Results"); /**All-knowing tree with reconstruction information*/
200  ftree_shwf->Branch("run", &fRun, "run/I");
201  ftree_shwf->Branch("subrun", &fSubRun, "subrun/I");
202  ftree_shwf->Branch("event", &fEvent, "event/I");
203  ftree_shwf->Branch("nplanes", &fNPlanes, "nplanes/I");
204  ftree_shwf->Branch("nangles", &fNAngles, "nangles/I");
205  ftree_shwf->Branch("xtheta", &xtheta, "xtheta/D");
206  ftree_shwf->Branch("xphi", &xphi, "xphi/D");
207  ftree_shwf->Branch("ftotChargeADC", "std::vector<double>", &fTotChargeADC);
208  ftree_shwf->Branch("ftotChargeMeV", "std::vector<double>", &fTotChargeMeV);
209  ftree_shwf->Branch("fTotChargeMeV_MIPs", "std::vector<double>", &fTotChargeMeV_MIPs);
210  ftree_shwf->Branch("NPitch", "std::vector< std::vector<double> >", &fNPitch);
211 
212  // this should be temporary - until the omega is sorted out.
213  ftree_shwf->Branch("RMS_2cm", "std::vector<double>", &fRMS_2cm);
214  ftree_shwf->Branch("Npoints_2cm", "std::vector<int>", &fNpoints_2cm);
215  ftree_shwf->Branch("ChargeADC_2cm", "std::vector<double>", &fChargeADC_2cm);
216  ftree_shwf->Branch("ChargeMeV_2cm", "std::vector<double>", &fChargeMeV_2cm);
217  ftree_shwf->Branch("ChargeMeV_2cm_refined", "std::vector<double>", &fChargeMeV_2cm_refined);
218  ftree_shwf->Branch("ChargeMeV_2cm_axsum", "std::vector<double>", &fChargeMeV_2cm_axsum);
219  ftree_shwf->Branch("fNhitsperplane", "std::vector<int>", &fNhitsperplane);
220  ftree_shwf->Branch("fTotADCperplane", "std::vector<double>", &fTotADCperplane);
221  ftree_shwf->Branch(
222  "ChargedistributionADC", "std::vector<std::vector<double>>", &fDistribChargeADC);
223  ftree_shwf->Branch(
224  "ChargedistributionMeV", "std::vector<std::vector<double>>", &fDistribChargeMeV);
225  ftree_shwf->Branch(
226  "DistribHalfChargeMeV", "std::vector<std::vector<double>>", &fDistribHalfChargeMeV);
227  ftree_shwf->Branch(
228  "ChargedistributionPosition", "std::vector<std::vector<double>>", &fDistribChargeposition);
229  ftree_shwf->Branch("xyz_vertex_fit", "std::vector<double>", &xyz_vertex_fit);
230  }
std::vector< double > fTotChargeADC
std::vector< int > fNhitsperplane
std::vector< std::vector< double > > fNPitch
std::vector< std::vector< double > > fDistribChargeMeV
std::vector< double > fTotADCperplane
std::vector< std::vector< double > > fDistribHalfChargeMeV
std::vector< double > fChargeMeV_2cm_refined
std::vector< int > fNpoints_2cm
std::vector< double > fRMS_2cm
std::vector< double > xyz_vertex_fit
std::vector< double > fChargeMeV_2cm
std::vector< double > fTotChargeMeV_MIPs
std::vector< double > fChargeADC_2cm
std::vector< std::vector< double > > fDistribChargeADC
art::ServiceHandle< art::TFileService > tfs
std::vector< double > fChargeMeV_2cm_axsum
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
std::vector< std::vector< double > > fDistribChargeposition
std::vector< double > fTotChargeMeV
void shwf::ShowerReco::beginRun ( art::Run &  run)
private

Definition at line 233 of file ShowerReco_module.cc.

234  {
235  auto const* geom = lar::providerFrom<geo::Geometry>();
236  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataForJob();
237  auto const detProp =
238  art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataForJob();
239 
240  fWirePitch = geom->WirePitch(); // wire pitch in cm
241  fTimeTick = sampling_rate(clockData) / 1000.;
242  fDriftVelocity = detProp.DriftVelocity(detProp.Efield(), detProp.Temperature());
244  }
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
auto const detProp
void shwf::ShowerReco::ClearandResizeVectors ( unsigned int  nPlanes)
private
void shwf::ShowerReco::GetVertexAndAnglesFromCluster ( art::Ptr< recob::Cluster clust,
unsigned int  plane 
)
private

Actual routine that reconstruct the shower

Definition at line 815 of file ShowerReco_module.cc.

817  {
818  // convert to cm/cm units needed in the calculation
819  angle[plane] = clust->StartAngle();
820  slope[plane] = std::tan(clust->StartAngle());
821  fWire_vertex[plane] = clust->StartWire();
822  fTime_vertex[plane] = clust->StartTick();
823 
824  fWire_vertexError[plane] = clust->SigmaStartWire(); // wire coordinate of vertex for each plane
825  fTime_vertexError[plane] = clust->SigmaStartTick(); // time coordinate of vertex for each plane
826 
827  fWire_last[plane] = clust->EndWire(); // wire coordinate of last point for each plane
828  fTime_last[plane] = clust->EndTick();
829  }
std::vector< double > fTime_last
std::vector< double > fWire_vertexError
std::vector< double > fTime_vertex
std::vector< unsigned int > fWire_last
std::vector< double > fTime_vertexError
std::vector< unsigned int > fWire_vertex
void shwf::ShowerReco::LongTransEnergy ( geo::GeometryCore const *  geom,
detinfo::DetectorClocksData const &  clockData,
detinfo::DetectorPropertiesData const &  detProp,
unsigned int  set,
std::vector< art::Ptr< recob::Hit >>  hitlist 
)
private

third loop to get only points inside of 1RMS of value.

Definition at line 614 of file ShowerReco_module.cc.

619  {
620  // alogorithm for energy vs dx of the shower (roto-translation) COLLECTION
621  // VIEW
622 
624 
625  double totCnrg = 0,
626  totCnrg_corr = 0; // tot enegry of the shower in collection
627 
628  double time;
629  unsigned int wire = 0, plane = fNPlanes - 1;
630 
631  double mevav2cm = 0.;
632  double sum = 0.;
633  double npoints_calo = 0;
634 
635  int direction = -1;
636 
637  //override direction if phi (XZ angle) is less than 90 degrees
638  if (fabs(xphi) < 90) direction = 1;
639 
640  //variables to check whether a hit is close to the shower axis.
641  double ortdist, linedist;
642  double wire_on_line, time_on_line;
643 
644  //get effective pitch using 3D angles
645  util::GeometryUtilities const gser{*geom, clockData, detProp};
646  double newpitch = gser.PitchInView(plane, xphi, xtheta);
647 
648  using lar::to_element;
650  for (auto const& hit : hitlist | transform(to_element)) {
651  time = hit.PeakTime();
652  wire = hit.WireID().Wire;
653  plane = hit.WireID().Plane;
654 
655  double dEdx_new;
656 
657  if (fUseArea) { dEdx_new = calalg.dEdx_AREA(clockData, detProp, hit, newpitch); }
658  else // this will hopefully go away, once all of the calibration factors
659  // are calculated.
660  {
661  dEdx_new = calalg.dEdx_AMP(clockData, detProp, hit, newpitch);
662  }
663 
664  //calculate total energy.
665  totCnrg_corr += dEdx_new;
666 
667  // calculate the wire,time coordinates of the hit projection on to the 2D shower axis
668  gser.GetPointOnLine(slope[plane] / fWireTimetoCmCm,
669  fWire_vertex[plane],
670  fTime_vertex[plane],
671  wire,
672  time,
673  wire_on_line,
674  time_on_line);
675  linedist =
676  gser.Get2DDistance(wire_on_line, time_on_line, fWire_vertex[plane], fTime_vertex[plane]);
677  ortdist = gser.Get2DDistance(wire_on_line, time_on_line, wire, time);
678 
679  //calculate the distance from the vertex using the effective pitch metric
680  double wdist = (((double)wire - (double)fWire_vertex[plane]) * newpitch) *
681  direction; //wdist is always positive
682 
683  if ((wdist < fcalodEdxlength) && (wdist > 0.2)) {
684 
685  vdEdx.push_back(dEdx_new);
686  vresRange.push_back(fabs(wdist));
687  vdQdx.push_back(hit.PeakAmplitude() / newpitch);
688  Trk_Length = wdist;
689  fTrkPitchC = fNPitch[set][plane];
690  Kin_En += dEdx_new * newpitch;
691  npoints_calo++;
692  sum += dEdx_new;
693 
694  if (wdist < fdEdxlength &&
695  ((direction == 1 && wire > fWire_vertex[plane]) // take no hits before vertex
696  // (depending on direction)
697  || (direction == -1 && wire < fWire_vertex[plane])) &&
698  ortdist < 4.5 && linedist < fdEdxlength) {
699  fChargeMeV_2cm[set] += dEdx_new;
700  fNpoints_2cm[set]++;
701  }
702 
703  // fill out for 4cm preshower
704 
705  fDistribChargeMeV[set].push_back(dEdx_new); // vector with the first De/Dx points
706  fDistribChargeposition[set].push_back(
707  wdist); // vector with the first De/Dx points' positions
708 
709  } // end inside range if statement
710 
711  } // end first loop on hits.
712 
713  auto const signalType =
714  hitlist.empty() ? geo::kMysteryType : geom->SignalType(hitlist.front()->WireID());
715 
716  if (signalType == geo::kCollection) {
717  fTotChargeADC[set] = totCnrg * newpitch;
718  fTotChargeMeV[set] = totCnrg_corr * newpitch;
719  }
720 
721  // calculate average dE/dx
722  if (fNpoints_2cm[set] > 0) { mevav2cm = fChargeMeV_2cm[set] / fNpoints_2cm[set]; }
723 
724  // second loop to calculate RMS
725  for (auto const& hit : hitlist | transform(to_element)) {
726  time = hit.PeakTime();
727  wire = hit.WireID().Wire;
728  plane = hit.WireID().Plane;
729  double dEdx = 0;
730 
731  if (fUseArea) { dEdx = calalg.dEdx_AREA(clockData, detProp, hit, newpitch); }
732  else // this will hopefully go away, once all of the calibration factors
733  // are calculated.
734  {
735  dEdx = calalg.dEdx_AMP(clockData, detProp, hit, newpitch);
736  }
737 
738  gser.GetPointOnLine(slope[plane] / fWireTimetoCmCm,
739  fWire_vertex[plane],
740  fTime_vertex[plane],
741  wire,
742  time,
743  wire_on_line,
744  time_on_line);
745  linedist =
746  gser.Get2DDistance(wire_on_line, time_on_line, fWire_vertex[plane], fTime_vertex[plane]);
747  ortdist = gser.Get2DDistance(wire_on_line, time_on_line, wire, time);
748 
749  double wdist = (((double)wire - (double)fWire_vertex[plane]) * newpitch) * direction;
750 
751  if ((wdist < fcalodEdxlength) && (wdist > 0.2)) {
752  if (wdist < fdEdxlength &&
753  ((direction == 1 && wire > fWire_vertex[plane]) ||
754  (direction == -1 && wire < fWire_vertex[plane])) &&
755  ortdist < 4.5 && linedist < fdEdxlength) {
756  fRMS_2cm[set] += (dEdx - mevav2cm) * (dEdx - mevav2cm);
757  }
758 
759  } // end if on correct hits.
760  } // end RMS_calculating loop.
761 
762  if (fNpoints_2cm[set] > 0) { fRMS_2cm[set] = TMath::Sqrt(fRMS_2cm[set] / fNpoints_2cm[set]); }
763 
764  /// third loop to get only points inside of 1RMS of value.
765 
766  for (auto const& hit : hitlist | transform(to_element)) {
767  time = hit.PeakTime();
768  wire = hit.WireID().Wire;
769  plane = hit.WireID().Plane;
770 
771  double dEdx = 0;
772  if (fUseArea) { dEdx = calalg.dEdx_AREA(clockData, detProp, hit, newpitch); }
773  else // this will hopefully go away, once all of the calibration factors
774  // are calculated.
775  {
776  dEdx = calalg.dEdx_AMP(clockData, detProp, hit, newpitch);
777  }
778 
779  gser.GetPointOnLine(slope[plane] / fWireTimetoCmCm,
780  fWire_vertex[plane],
781  fTime_vertex[plane],
782  wire,
783  time,
784  wire_on_line,
785  time_on_line);
786  linedist =
787  gser.Get2DDistance(wire_on_line, time_on_line, fWire_vertex[plane], fTime_vertex[plane]);
788  ortdist = gser.Get2DDistance(wire_on_line, time_on_line, wire, time);
789 
790  double wdist = (((double)wire - (double)fWire_vertex[plane]) * newpitch) * direction;
791 
792  if ((wdist < fcalodEdxlength) && (wdist > 0.2 &&
793  ((direction == 1 && wire > fWire_vertex[plane]) ||
794  (direction == -1 && wire < fWire_vertex[plane])) &&
795  ortdist < 4.5 && linedist < fdEdxlength)) {
796  if (wdist < fdEdxlength) {
797  if (((dEdx > (mevav2cm - fRMS_2cm[set])) && (dEdx < (mevav2cm + fRMS_2cm[set]))) ||
798  (newpitch > 0.3 * fdEdxlength)) {
799  fCorr_MeV_2cm[set] += dEdx;
800  fNpoints_corr_MeV_2cm[set]++;
801  }
802 
803  } // end if on good hits
804  }
805  } // end of third loop on hits
806 
807  if (fNpoints_corr_MeV_2cm[set] > 0) {
810  }
811  }
fhicl::ParameterSet fCaloPSet
std::vector< double > fTotChargeADC
constexpr to_element_t to_element
Definition: ToElement.h:24
Who knows?
Definition: geo_types.h:147
static constexpr Sample_t transform(Sample_t sample)
std::vector< std::vector< double > > fNPitch
std::vector< std::vector< double > > fDistribChargeMeV
std::vector< int > fNpoints_corr_MeV_2cm
std::vector< double > fTime_vertex
process_name hit
Definition: cheaterreco.fcl:51
std::vector< double > fChargeMeV_2cm_refined
std::vector< float > vdEdx
std::vector< float > vdQdx
std::vector< int > fNpoints_2cm
std::vector< float > vresRange
std::vector< double > fRMS_2cm
float dEdx(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, const TCSlice &slc, TP3D &tp3d)
Definition: PFPUtils.cxx:2687
std::vector< double > fCorr_MeV_2cm
std::vector< double > fChargeMeV_2cm
Var Sqrt(const Var &v)
Use to take sqrt of a var.
Definition: Var.cxx:172
std::vector< std::vector< double > > fDistribChargeposition
std::vector< unsigned int > fWire_vertex
auto const detProp
std::vector< double > fTotChargeMeV
Signal from collection planes.
Definition: geo_types.h:146
void shwf::ShowerReco::produce ( art::Event &  evt)
private

Get Clusters

Todo:
really need to determine the values of the arguments of the recob::Shower ctor

Definition at line 352 of file ShowerReco_module.cc.

353  {
354  auto const* geom = lar::providerFrom<geo::Geometry>();
355  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
356  auto const detProp =
357  art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(evt, clockData);
358 
359  util::GeometryUtilities const gser{*geom, clockData, detProp};
360  fNPlanes = geom->Nplanes();
361  auto Shower3DVector = std::make_unique<std::vector<recob::Shower>>();
362  auto cassn = std::make_unique<art::Assns<recob::Shower, recob::Cluster>>();
363  auto hassn = std::make_unique<art::Assns<recob::Shower, recob::Hit>>();
364  auto calorimetrycol = std::make_unique<std::vector<anab::Calorimetry>>();
365  auto calassn = std::make_unique<art::Assns<anab::Calorimetry, recob::Shower>>();
366 
367  /**Get Clusters*/
368 
369  art::Handle<std::vector<recob::Cluster>> clusterListHandle;
370  evt.getByLabel(fClusterModuleLabel, clusterListHandle);
371 
372  art::Handle<std::vector<art::PtrVector<recob::Cluster>>> clusterAssociationHandle;
373  evt.getByLabel(fClusterModuleLabel, clusterAssociationHandle);
374 
375  art::FindManyP<recob::Hit> fmh(clusterListHandle, evt, fClusterModuleLabel);
376 
377  fRun = evt.id().run();
378  fSubRun = evt.id().subRun();
379  fEvent = evt.id().event();
380 
381  // find all the hits associated to all the clusters (once and for all);
382  // the index of the query matches the index of the cluster in the collection
383  // (conveniently carried around in its art pointer)
384  art::FindManyP<recob::Hit> ClusterHits(clusterListHandle, evt, fClusterModuleLabel);
385 
386  std::vector<art::PtrVector<recob::Cluster>>::const_iterator clusterSet =
387  clusterAssociationHandle->begin();
388  // loop over vector of vectors (each size of NPlanes) and reconstruct showers from each of those
389  for (size_t iClustSet = 0; iClustSet < clusterAssociationHandle->size(); iClustSet++) {
390 
391  const art::PtrVector<recob::Cluster> CurrentClusters = (*(clusterSet++));
392 
393  // do some error checking - i.e. are the clusters themselves present.
394  if (clusterListHandle->size() < 2 || CurrentClusters.size() < 2) {
395  ftree_shwf->Fill();
396  return;
397  }
398 
400 
401  std::vector<std::vector<art::Ptr<recob::Hit>>> hitlist_all;
402  hitlist_all.resize(fNPlanes);
403 
404  for (size_t iClust = 0; iClust < CurrentClusters.size(); iClust++) {
405  art::Ptr<recob::Cluster> const& pclust = CurrentClusters[iClust];
406 
407  // get all the hits for this cluster;
408  // pclust is a art::Ptr to the original cluster collection stored in the event;
409  // its key corresponds to its index in the collection
410  // (and therefore to the index in the query)
411  std::vector<art::Ptr<recob::Hit>> const& hitlist = ClusterHits.at(pclust.key());
412 
413  unsigned int p(0); //c=channel, p=plane, w=wire
414 
415  if (hitlist.size() == 0) continue;
416 
417  p = (*hitlist.begin())->WireID().Plane;
418  // get vertex position and slope information to start with - ii is the
419  // posistion of the correct cluster:
421 
422  double ADCcharge = 0;
423  //loop over cluster hits
424  for (art::Ptr<recob::Hit> const& hit : hitlist) {
425  p = hit->WireID().Plane;
426  hitlist_all[p].push_back(hit);
427  ADCcharge += hit->PeakAmplitude();
428  }
429  fNhitsperplane[p] = hitlist_all[p].size();
430  fTotADCperplane[p] = ADCcharge;
431  } // End loop on clusters.
432  // Now I have the Hitlists and the relevent clusters parameters saved.
433 
434  // find best set:
435  unsigned int bp1 = 0, bp2 = 0;
436  double minerror1 = 99999999, minerror2 = 9999999;
437  for (unsigned int ii = 0; ii < fNPlanes; ++ii) {
438  double locerror =
440  fTime_vertexError[ii] * fTime_vertexError[ii]; // time coordinate of vertex for each plane
441 
442  if (minerror1 >= locerror) // the >= sign is to favorize collection
443  {
444  minerror1 = locerror;
445  bp1 = ii;
446  }
447  }
448  for (unsigned int ij = 0; ij < fNPlanes; ++ij) {
449  double locerror =
451  fTime_vertexError[ij] * fTime_vertexError[ij]; // time coordinate of vertex for each plane
452 
453  if (minerror2 >= locerror && ij != bp1) {
454  minerror2 = locerror;
455  bp2 = ij;
456  }
457  }
458 
459  gser.Get3DaxisN(bp1, bp2, angle[bp1], angle[bp2], xphi, xtheta);
460 
461  ///////////////////////////////////////////////////////////
462  const double origin[3] = {0.};
463  std::vector<std::vector<double>> position;
464  double fTimeTick = sampling_rate(clockData) / 1000.;
465  double fDriftVelocity = detProp.DriftVelocity(detProp.Efield(), detProp.Temperature());
466  // get starting positions for all planes
467  for (unsigned int xx = 0; xx < fNPlanes; xx++) {
468  double pos1[3];
469  geom->Plane(xx).LocalToWorld(origin, pos1);
470  std::vector<double> pos2;
471  pos2.push_back(pos1[0]);
472  pos2.push_back(pos1[1]);
473  pos2.push_back(pos1[2]);
474  position.push_back(pos2);
475  }
476  // Assuming there is no problem ( and we found the best pair that comes
477  // close in time ) we try to get the Y and Z coordinates for the start of
478  // the shower.
479  try {
480  int chan1 = geom->PlaneWireToChannel(bp1, fWire_vertex[bp1], 0);
481  int chan2 = geom->PlaneWireToChannel(bp2, fWire_vertex[bp2], 0);
482 
483  double y, z;
484  geom->ChannelsIntersect(chan1, chan2, y, z);
485 
486  xyz_vertex_fit[1] = y;
487  xyz_vertex_fit[2] = z;
488  xyz_vertex_fit[0] =
489  (fTime_vertex[bp1] - trigger_offset(clockData)) * fDriftVelocity * fTimeTick +
490  position[0][0];
491  }
492  catch (cet::exception const& e) {
493  mf::LogWarning("ShowerReco") << "caught exception \n" << e;
494  xyz_vertex_fit[1] = 0;
495  xyz_vertex_fit[2] = 0;
496  xyz_vertex_fit[0] = 0;
497  }
498 
499  // if collection is not best plane, project starting point from that
500  if (bp1 != fNPlanes - 1 && bp2 != fNPlanes - 1) {
501  double pos[3];
502  unsigned int wirevertex;
503 
504  geom->Plane(fNPlanes - 1).LocalToWorld(origin, pos);
505 
506  pos[1] = xyz_vertex_fit[1];
507  pos[2] = xyz_vertex_fit[2];
508  wirevertex = geom->NearestWire(pos, fNPlanes - 1);
509 
510  double drifttick =
511  (xyz_vertex_fit[0] / detProp.DriftVelocity(detProp.Efield(), detProp.Temperature())) *
512  (1. / fTimeTick);
513  fWire_vertex[fNPlanes - 1] = wirevertex; // wire coordinate of vertex for each plane
514  fTime_vertex[fNPlanes - 1] =
515  drifttick -
516  (pos[0] / detProp.DriftVelocity(detProp.Efield(), detProp.Temperature())) *
517  (1. / fTimeTick) +
518  trigger_offset(clockData);
519  }
520 
521  if (fabs(xphi) < 5.) {
522  xtheta = gser.Get3DSpecialCaseTheta(
523  bp1, bp2, fWire_last[bp1] - fWire_vertex[bp1], fWire_last[bp2] - fWire_vertex[bp2]);
524  }
525 
526  // zero the arrays just to make sure
527  for (unsigned int i = 0; i < fNAngles; ++i) {
528  fTotChargeADC[i] = 0;
529  fTotChargeMeV[i] = 0;
530  fTotChargeMeV_MIPs[i] = 0;
531  fNpoints_corr_ADC_2cm[i] = 0;
532  fNpoints_corr_MeV_2cm[i] = 0;
533 
534  fRMS_2cm[i] = 0;
535  fNpoints_2cm[i] = 0;
536 
537  fCorr_MeV_2cm[i] = 0;
538  fCorr_Charge_2cm[i] = 0;
539 
540  fChargeADC_2cm[i] = 0; //Initial charge in ADC/cm for each plane angle calculation 4cm
541  fChargeMeV_2cm[i] = 0; //initial charge in MeV/cm for each angle calculation first 4cm
542  }
543 
544  // do loop and choose Collection. With ful calorimetry can do all.
545  if (!(fabs(xphi) > 89 && fabs(xphi) < 91)) // do not calculate pitch for extreme angles
546  LongTransEnergy(geom,
547  clockData,
548  detProp,
549  0,
550  hitlist_all[fNPlanes - 1]); // temporary only plane 2.
551 
552  //////create spacepoints, and direction cosines for Shower creation
553 
554  // make an art::PtrVector of the clusters
555  art::PtrVector<recob::Cluster> prodvec;
556  for (unsigned int i = 0; i < clusterListHandle->size(); ++i) {
557  art::Ptr<recob::Cluster> prod(clusterListHandle, i);
558  prodvec.push_back(prod);
559  }
560 
561  //create a singleSpacePoint at vertex.
562  std::vector<recob::SpacePoint> spcpts;
563 
564  //get direction cosines and set them for the shower
565  // TBD determine which angle to use for the actual shower
566  double fPhi = xphi;
567  double fTheta = xtheta;
568 
569  TVector3 dcosVtx(
570  TMath::Cos(fPhi * TMath::Pi() / 180) * TMath::Sin(fTheta * TMath::Pi() / 180),
571  TMath::Cos(fTheta * TMath::Pi() / 180),
572  TMath::Sin(fPhi * TMath::Pi() / 180) * TMath::Sin(fTheta * TMath::Pi() / 180));
573  /// \todo really need to determine the values of the arguments of the
574  /// recob::Shower ctor
575  // fill with bogus values for now
576  TVector3 dcosVtxErr(util::kBogusD, util::kBogusD, util::kBogusD);
577  recob::Shower singShower;
578  singShower.set_direction(dcosVtx);
579  singShower.set_direction_err(dcosVtxErr);
580 
581  Shower3DVector->push_back(singShower);
582  // associate the shower with its clusters
583  util::CreateAssn(evt, *Shower3DVector, prodvec, *cassn);
584 
585  // get the hits associated with each cluster and associate those with the shower
586  for (size_t p = 0; p < prodvec.size(); ++p) {
587  std::vector<art::Ptr<recob::Hit>> hits = fmh.at(p);
588  util::CreateAssn(evt, *Shower3DVector, hits, *hassn);
589  }
590 
591  geo::PlaneID planeID(0, 0, fNPlanes - 1);
592  calorimetrycol->emplace_back(
594 
595  art::PtrVector<recob::Shower> ssvec;
596 
597  art::ProductID aid = evt.getProductID<std::vector<recob::Shower>>();
598  art::Ptr<recob::Shower> aptr(aid, 0, evt.productGetter(aid));
599  ssvec.push_back(aptr);
600 
601  util::CreateAssn(evt, *calorimetrycol, ssvec, *calassn);
602  ftree_shwf->Fill();
603  } // end loop on Vectors of "Associated clusters"
604 
605  evt.put(std::move(Shower3DVector));
606  evt.put(std::move(cassn));
607  evt.put(std::move(hassn));
608  evt.put(std::move(calorimetrycol));
609  evt.put(std::move(calassn));
610  }
process_name opflash particleana ie ie ie z
std::vector< double > fWire_vertexError
std::vector< double > fTotChargeADC
std::vector< double > fCorr_Charge_2cm
void set_direction_err(const TVector3 &dir_e)
Definition: Shower.h:136
unsigned int event
Definition: DataStructs.h:634
unsigned int run
Definition: DataStructs.h:635
std::vector< int > fNhitsperplane
pdgs p
Definition: selectors.fcl:22
The data type to uniquely identify a Plane.
Definition: geo_types.h:472
std::vector< int > fNpoints_corr_MeV_2cm
std::vector< double > fTime_vertex
std::vector< double > fTotADCperplane
std::vector< unsigned int > fWire_last
process_name hit
Definition: cheaterreco.fcl:51
std::vector< float > vdEdx
void GetVertexAndAnglesFromCluster(art::Ptr< recob::Cluster > clust, unsigned int plane)
std::vector< float > vdQdx
std::vector< int > fNpoints_2cm
std::vector< double > fTime_vertexError
void set_direction(const TVector3 &dir)
Definition: Shower.h:135
process_name opflash particleana ie ie y
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
std::vector< float > vresRange
std::vector< double > fRMS_2cm
std::vector< double > fCorr_MeV_2cm
std::string fClusterModuleLabel
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.
std::vector< double > xyz_vertex_fit
std::vector< double > fChargeMeV_2cm
std::vector< double > fTotChargeMeV_MIPs
std::vector< double > fChargeADC_2cm
do i e
std::vector< float > deadwire
int trigger_offset(DetectorClocksData const &data)
void LongTransEnergy(geo::GeometryCore const *geom, detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, unsigned int set, std::vector< art::Ptr< recob::Hit >> hitlist)
unsigned int subRun
Definition: DataStructs.h:636
constexpr double kBogusD
obviously bogus double value
void ClearandResizeVectors(unsigned int nPlanes)
TCEvent evt
Definition: DataStructs.cxx:8
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
std::vector< int > fNpoints_corr_ADC_2cm
std::vector< unsigned int > fWire_vertex
auto const detProp
constexpr Point origin()
Returns a origin position with a point of the specified type.
Definition: geo_vectors.h:227
std::vector< double > fTotChargeMeV

Member Data Documentation

float shwf::ShowerReco::angle[3]
private

Definition at line 81 of file ShowerReco_module.cc.

std::vector<float> shwf::ShowerReco::deadwire
private

Definition at line 138 of file ShowerReco_module.cc.

float shwf::ShowerReco::fcalodEdxlength
private

Definition at line 143 of file ShowerReco_module.cc.

fhicl::ParameterSet shwf::ShowerReco::fCaloPSet
private

Definition at line 88 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fChargeADC_2cm
private

Definition at line 102 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fChargeMeV_2cm
private

Definition at line 103 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fChargeMeV_2cm_axsum
private

Definition at line 106 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fChargeMeV_2cm_refined
private

Definition at line 105 of file ShowerReco_module.cc.

std::string shwf::ShowerReco::fClusterModuleLabel
private

Definition at line 83 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fCorr_Charge_2cm
private

Definition at line 93 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fCorr_MeV_2cm
private

Definition at line 92 of file ShowerReco_module.cc.

float shwf::ShowerReco::fdEdxlength
private

Definition at line 141 of file ShowerReco_module.cc.

std::vector<std::vector<double> > shwf::ShowerReco::fDistribChargeADC
private

Definition at line 108 of file ShowerReco_module.cc.

std::vector<std::vector<double> > shwf::ShowerReco::fDistribChargeMeV
private

Definition at line 110 of file ShowerReco_module.cc.

std::vector<std::vector<double> > shwf::ShowerReco::fDistribChargeposition
private

Definition at line 113 of file ShowerReco_module.cc.

std::vector<std::vector<double> > shwf::ShowerReco::fDistribHalfChargeMeV
private

Definition at line 111 of file ShowerReco_module.cc.

double shwf::ShowerReco::fDriftVelocity
private

Definition at line 154 of file ShowerReco_module.cc.

int shwf::ShowerReco::fEvent
private

Definition at line 78 of file ShowerReco_module.cc.

double shwf::ShowerReco::fMean_wire_pitch
private

Definition at line 87 of file ShowerReco_module.cc.

unsigned int shwf::ShowerReco::fNAngles
private

Definition at line 148 of file ShowerReco_module.cc.

std::vector<int> shwf::ShowerReco::fNhitsperplane
private

Definition at line 157 of file ShowerReco_module.cc.

std::vector<std::vector<double> > shwf::ShowerReco::fNPitch
private

Definition at line 131 of file ShowerReco_module.cc.

unsigned int shwf::ShowerReco::fNPlanes
private

Definition at line 147 of file ShowerReco_module.cc.

std::vector<int> shwf::ShowerReco::fNpoints_2cm
private

Definition at line 91 of file ShowerReco_module.cc.

std::vector<int> shwf::ShowerReco::fNpoints_corr_ADC_2cm
private

Definition at line 95 of file ShowerReco_module.cc.

std::vector<int> shwf::ShowerReco::fNpoints_corr_MeV_2cm
private

Definition at line 96 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fRMS_2cm
private

Definition at line 90 of file ShowerReco_module.cc.

int shwf::ShowerReco::fRun
private

Definition at line 78 of file ShowerReco_module.cc.

std::vector<std::vector<double> > shwf::ShowerReco::fSingleEvtAngle
private

Definition at line 115 of file ShowerReco_module.cc.

std::vector<std::vector<double> > shwf::ShowerReco::fSingleEvtAngleVal
private

Definition at line 116 of file ShowerReco_module.cc.

int shwf::ShowerReco::fSubRun
private

Definition at line 78 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fTime_last
private

Definition at line 125 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fTime_vertex
private

Definition at line 119 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fTime_vertexError
private

Definition at line 122 of file ShowerReco_module.cc.

float shwf::ShowerReco::ftimetick
private

Definition at line 85 of file ShowerReco_module.cc.

double shwf::ShowerReco::fTimeTick
private

Definition at line 153 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fTotADCperplane
private

Definition at line 158 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fTotChargeADC
private

Definition at line 98 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fTotChargeMeV
private

Definition at line 99 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fTotChargeMeV_MIPs
private

Definition at line 100 of file ShowerReco_module.cc.

TTree* shwf::ShowerReco::ftree_shwf
private

Definition at line 149 of file ShowerReco_module.cc.

float shwf::ShowerReco::fTrkPitchC
private

Definition at line 140 of file ShowerReco_module.cc.

bool shwf::ShowerReco::fUseArea
private

Definition at line 144 of file ShowerReco_module.cc.

std::vector<unsigned int> shwf::ShowerReco::fWire_last
private

Definition at line 124 of file ShowerReco_module.cc.

std::vector<unsigned int> shwf::ShowerReco::fWire_vertex
private

Definition at line 118 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::fWire_vertexError
private

Definition at line 121 of file ShowerReco_module.cc.

double shwf::ShowerReco::fWirePitch
private

Definition at line 152 of file ShowerReco_module.cc.

double shwf::ShowerReco::fWireTimetoCmCm
private

Definition at line 155 of file ShowerReco_module.cc.

float shwf::ShowerReco::Kin_En
private

Definition at line 134 of file ShowerReco_module.cc.

float shwf::ShowerReco::slope[3]
private

Definition at line 80 of file ShowerReco_module.cc.

float shwf::ShowerReco::Trk_Length
private

Definition at line 139 of file ShowerReco_module.cc.

std::vector<float> shwf::ShowerReco::vdEdx
private

Definition at line 135 of file ShowerReco_module.cc.

std::vector<float> shwf::ShowerReco::vdQdx
private

Definition at line 137 of file ShowerReco_module.cc.

std::vector<float> shwf::ShowerReco::vresRange
private

Definition at line 136 of file ShowerReco_module.cc.

double shwf::ShowerReco::xphi
private

Definition at line 146 of file ShowerReco_module.cc.

double shwf::ShowerReco::xtheta
private

Definition at line 146 of file ShowerReco_module.cc.

std::vector<double> shwf::ShowerReco::xyz_vertex_fit
private

Definition at line 128 of file ShowerReco_module.cc.


The documentation for this class was generated from the following file: