Public Member Functions | Private Attributes | List of all members
ex::HelloAuxDet Class Reference
Inheritance diagram for ex::HelloAuxDet:
art::EDAnalyzer art::detail::Analyzer art::detail::LegacyModule art::Observer art::ModuleBase

Public Member Functions

 HelloAuxDet (fhicl::ParameterSet const &p)
 
 HelloAuxDet (HelloAuxDet const &)=delete
 
 HelloAuxDet (HelloAuxDet &&)=delete
 
HelloAuxDetoperator= (HelloAuxDet const &)=delete
 
HelloAuxDetoperator= (HelloAuxDet &&)=delete
 
void analyze (art::Event const &e) override
 
void beginJob () override
 
- Public Member Functions inherited from art::EDAnalyzer
 EDAnalyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDAnalyzer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Analyzer
virtual ~Analyzer () noexcept
 
 Analyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 Analyzer (Table< Config > const &config)
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Observer
 ~Observer () noexcept
 
 Observer (Observer const &)=delete
 
 Observer (Observer &&)=delete
 
Observeroperator= (Observer const &)=delete
 
Observeroperator= (Observer &&)=delete
 
void registerProducts (ProductDescriptions &, ModuleDescription const &)
 
void fillDescriptions (ModuleDescription const &)
 
fhicl::ParameterSetID selectorConfig () const
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Private Attributes

std::vector< art::InputTagfDigitLabels
 
std::vector< art::InputTagfSimLabels
 
art::InputTag fPartLabel
 
TH1D * fXDistToTrueTrackFront
 
TH1D * fYDistToTrueTrackFront
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
std::string const & processName () const
 
bool wantAllEvents () const noexcept
 
bool wantEvent (ScheduleID id, Event const &e) const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) const
 
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &select_paths, std::vector< std::string > const &reject_paths, fhicl::ParameterSet const &config)
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 99 of file HelloAuxDet_module.cc.

Constructor & Destructor Documentation

ex::HelloAuxDet::HelloAuxDet ( fhicl::ParameterSet const &  p)
explicit

Definition at line 130 of file HelloAuxDet_module.cc.

131  :
132  EDAnalyzer(p) // ,
133  // More initializers here.
134 {
135  fDigitLabels = p.get<std::vector<art::InputTag>>("DigitLabels", {});
136  fSimLabels = p.get<std::vector<art::InputTag>>("SimLabels", {});
137  fPartLabel = p.get<art::InputTag>("PartLabel", "largeant");
138 }
std::vector< art::InputTag > fDigitLabels
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
p
Definition: test.py:223
art::InputTag fPartLabel
std::vector< art::InputTag > fSimLabels
ex::HelloAuxDet::HelloAuxDet ( HelloAuxDet const &  )
delete
ex::HelloAuxDet::HelloAuxDet ( HelloAuxDet &&  )
delete

Member Function Documentation

void ex::HelloAuxDet::analyze ( art::Event const &  e)
overridevirtual

Implements art::EDAnalyzer.

Definition at line 140 of file HelloAuxDet_module.cc.

141 {
142  //Print() AuxDetDigits
143  //::PrintAllIfEmpty<raw::AuxDetDigit>(fDigitLabels, e);
144 
145  //Print() AuxDetSimChannels
146  //::PrintAllIfEmpty<sim::AuxDetSimChannel>(fSimLabels, e);
147 
148  //Produce a map from TrackId to MCParticle
149  const auto parts = e.getValidHandle<std::vector<simb::MCParticle>>(fPartLabel);
150  std::map<int, simb::MCParticle> idToPart; //Map of Geant TrackId to particle
151  for(const auto& part: *parts) idToPart[part.TrackId()] = part;
152 
153  idToPart[-1] = idToPart[0]; //TODO: Remove this hack that makes sure the primary particle is recognized
154  mf::LogInfo("Number of Particles") << "There are " << parts->size() << " simb::MCParticles in this event.\n";
155 
156  //Make histograms showing that channel mapping works for AuxDetSimChannels
158  for(const auto& label: fSimLabels)
159  {
160  const auto channels = e.getValidHandle<std::vector<sim::AuxDetSimChannel>>(label);
161  for(const auto& channel: *channels)
162  {
163  const auto& auxDet = geom->AuxDet(channel.AuxDetID());
164  const auto& sens = auxDet.SensitiveVolume(channel.AuxDetSensitiveID());
165  const auto center = sens.GetCenter();
166  TLorentzVector centerV(center.X(), center.Y(), center.Z(), 0.); //TODO: Use GenVector when I have time to figure it out
167  for(const auto& ide: channel.AuxDetIDEs())
168  {
169  const auto found = idToPart.find(ide.trackID);
170  if(found != idToPart.end())
171  {
172  const auto& part = found->second;
173  const auto closestPt = std::min_element(part.Trajectory().begin(), part.Trajectory().end(),
174  [&centerV](const auto& first, const auto& second)
175  {
176  //SimulationBase coordinates are in mm, but Geometry service coordinates are in cm!
177  //Convert to cm.
178  return (first.first*0.1-centerV).Vect().Mag() < (second.first*0.1-centerV).Vect().Mag();
179  });
180 
181  if(closestPt != part.Trajectory().end())
182  {
183  //Figure out the x and y positions where the CRT was hit
184  //TODO: This is a hack to partially get around my lack of knowledge about the CRT channel mapping.
185  // I'd like to figure this out by identifying the AuxDet.
186  if(sens.toWorldCoords(geo::AuxDetSensitiveGeo::LocalPoint_t(0, 0, sens.HalfLength())).Y() - center.Y()
187  > sens.toWorldCoords(geo::AuxDetSensitiveGeo::LocalPoint_t(0, 0, sens.HalfLength())).X() - center.X())
188  //If a CRT hit that gives x position
189  {
190  fXDistToTrueTrackFront->Fill(centerV.X() - closestPt->first.X());
191  }
192  else //This is a CRT hit that gives y position
193  {
194  fYDistToTrueTrackFront->Fill(centerV.Y() - closestPt->first.Y());
195  }
196  }
197  else mf::LogWarning("No Closest Point") << "No point found that is considered \"closest\" to AuxDetSensitive!\n";
198  }
199  else mf::LogWarning("No Such TrackID") << "No map entry for TrackID " << ide.trackID << "\n";
200  }
201  }
202  }
203 }
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
AuxDetSensitiveGeo const & SensitiveVolume(size_t sv) const
Definition: AuxDetGeo.h:171
uint8_t channel
Definition: CRTFragment.hh:201
const double e
art::InputTag fPartLabel
AuxDetGeo const & AuxDet(unsigned int const ad=0) const
Returns the specified auxiliary detector.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
def center(depos, point)
Definition: depos.py:117
geo::Point3DBase_t< AuxDetSensitiveGeoCoordinatesTag > LocalPoint_t
Type of points in the local GDML auxiliary detector frame.
second_as<> second
Type of time stored in seconds, in double precision.
Definition: spacetime.h:85
std::vector< art::InputTag > fSimLabels
void ex::HelloAuxDet::beginJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 205 of file HelloAuxDet_module.cc.

206 {
207  // Print all of the AuxDetGeos that the Geometry service knows about.
208  const auto geom = lar::providerFrom<geo::Geometry>();
209 
210  //Get AuxDetGeos and print their names.
211  for(size_t det = 0; det < geom->NAuxDets(); ++det)
212  {
213  const auto& geo = geom->AuxDet(det);
214  mf::LogInfo("AuxDetGeometry") << "AuxDetGeo number " << det << ", " << geo.Name() << ", is centered at " << geo.GetCenter()
215  << " and has sensitive volumes:\n";
216  for(size_t sens = 0; sens < geo.NSensitiveVolume(); ++sens)
217  {
218  const auto& sensGeo = geo.SensitiveVolume(sens);
219  const auto center = sensGeo.GetCenter();
220  const bool xPos = (sensGeo.toWorldCoords(geo::AuxDetSensitiveGeo::LocalPoint_t(0, 0, sensGeo.HalfLength())).Y() - center.Y()
221  > sensGeo.toWorldCoords(geo::AuxDetSensitiveGeo::LocalPoint_t(0, 0, sensGeo.HalfLength())).X() - center.X());
222  mf::LogInfo("AuxDetGeometry") << "Sensitive volume " << sens << ": " << sensGeo.TotalVolume()->GetName() << " with center "
223  << sensGeo.GetCenter() << " which gives information about " << std::string(xPos?"X":"Y") << " coordinates.\n";
224  }
225  }
226 
227  //Tell the framework that I intend to consume AuxDetDigits and AuxDetSimChannels
228  for(const auto& label: fDigitLabels) consumes<std::vector<raw::AuxDetDigit>>(label);
229  for(const auto& label: fSimLabels) consumes<std::vector<sim::AuxDetSimChannel>>(label);
230  consumes<std::vector<simb::MCParticle>>(fPartLabel);
231 
232  //Create histograms owned by "The Framework" that I will fill
234  fXDistToTrueTrackFront = tfs->make<TH1D>("XDistToTrueTrackFront", "X Distance from True Trajectory to CRT Hits"
235  ";X[mm];Events", 150, -30, 30);
236  fYDistToTrueTrackFront = tfs->make<TH1D>("YDistToTrueTrackFront", "Y Distance from True Trajectory to CRT Hits"
237  ";Y[mm];Events", 150, -30, 30);
238 }
std::string string
Definition: nybbler.cc:12
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
std::vector< art::InputTag > fDigitLabels
art::InputTag fPartLabel
def center(depos, point)
Definition: depos.py:117
geo::Point3DBase_t< AuxDetSensitiveGeoCoordinatesTag > LocalPoint_t
Type of points in the local GDML auxiliary detector frame.
LArSoft geometry interface.
Definition: ChannelGeo.h:16
std::vector< art::InputTag > fSimLabels
HelloAuxDet& ex::HelloAuxDet::operator= ( HelloAuxDet const &  )
delete
HelloAuxDet& ex::HelloAuxDet::operator= ( HelloAuxDet &&  )
delete

Member Data Documentation

std::vector<art::InputTag> ex::HelloAuxDet::fDigitLabels
private

Definition at line 120 of file HelloAuxDet_module.cc.

art::InputTag ex::HelloAuxDet::fPartLabel
private

Definition at line 122 of file HelloAuxDet_module.cc.

std::vector<art::InputTag> ex::HelloAuxDet::fSimLabels
private

Definition at line 121 of file HelloAuxDet_module.cc.

TH1D* ex::HelloAuxDet::fXDistToTrueTrackFront
private

Definition at line 125 of file HelloAuxDet_module.cc.

TH1D* ex::HelloAuxDet::fYDistToTrueTrackFront
private

Definition at line 126 of file HelloAuxDet_module.cc.


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