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

Public Member Functions

 protonanalysis (fhicl::ParameterSet const &pset)
 
virtual ~protonanalysis ()
 
void beginJob ()
 
void endJob ()
 
void beginRun (const art::Run &run)
 
void analyze (const art::Event &evt)
 
void reset ()
 
- 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

TH1F * pdg_beamparticles
 
TH1F * calibration_constants
 
TH2F * electric_field_YZ
 
TH2F * electric_field_YZ_3
 
TH1F * EFieldX
 
TH3F * EField3DX
 
TTree * fEventTree
 
Int_t run
 
Int_t subrun
 
Int_t event
 
Double_t Total_energy [kMaxTracks]
 
Double_t EndE [kMaxTracks]
 
Int_t pdg [kMaxTracks]
 
Double_t evttime
 
Int_t year_month_date
 
Int_t hour_min_sec
 
Int_t beam_proton
 
Int_t tot_trks
 
Int_t total_trks
 
Int_t pdg_all_beam [kMaxTracks]
 
Float_t xprojectedlen [kMaxTracks]
 
Float_t trackthetaxz [kMaxTracks]
 
Float_t trackthetayz [kMaxTracks]
 
Float_t trkstartx [kMaxTracks]
 
Float_t trkstarty [kMaxTracks]
 
Float_t trkstartz [kMaxTracks]
 
Float_t trkendx [kMaxTracks]
 
Float_t trkendy [kMaxTracks]
 
Float_t trkendz [kMaxTracks]
 
Float_t trklen [kMaxTracks]
 
Int_t TrkID [kMaxTracks]
 
Int_t TrackId_geant [kMaxTracks]
 
Int_t origin [kMaxTracks]
 
Float_t EndPointx [kMaxTracks]
 
Float_t EndPointy [kMaxTracks]
 
Float_t EndPointz [kMaxTracks]
 
Float_t EndPointx_corrected [kMaxTracks]
 
Float_t EndPointy_corrected [kMaxTracks]
 
Float_t EndPointz_corrected [kMaxTracks]
 
Float_t StartPointx [kMaxTracks]
 
Float_t StartPointy [kMaxTracks]
 
Float_t StartPointz [kMaxTracks]
 
Float_t trkstartcosxyz [kMaxTracks][3]
 
Float_t trkendcosxyz [kMaxTracks][3]
 
Int_t ntrkhits [kMaxTracks][3]
 
Float_t trkdqdx [kMaxTracks][3][3000]
 
Float_t trkdedx [kMaxTracks][3][3000]
 
Float_t trkresrange [kMaxTracks][3][3000]
 
Float_t trkhitx [kMaxTracks][3][3000]
 
Float_t trkhity [kMaxTracks][3][3000]
 
Float_t trkhitz [kMaxTracks][3][3000]
 
Float_t trkpitch [kMaxTracks][3][3000]
 
Float_t Efield [kMaxTracks][3][3000]
 
Float_t peakT_max [kMaxTracks]
 
Float_t peakT_min [kMaxTracks]
 
Double_t momentum [kMaxTracks]
 
std::string fHitsModuleLabel
 
std::string fTrackModuleLabel
 
std::string fCalorimetryModuleLabel
 
bool fSaveCaloInfo
 
bool fSaveTrackInfo
 

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 84 of file protonanalysis_module.cc.

Constructor & Destructor Documentation

dune::protonanalysis::protonanalysis ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 163 of file protonanalysis_module.cc.

163  :
164  EDAnalyzer(pset),
165  fHitsModuleLabel (pset.get< std::string >("HitsModuleLabel","") ),
166  fTrackModuleLabel (pset.get< std::string >("TrackModuleLabel","") ),
167  fCalorimetryModuleLabel (pset.get< std::string >("CalorimetryModuleLabel","") ),
168  fSaveCaloInfo (pset.get< bool>("SaveCaloInfo",false)),
169  fSaveTrackInfo (pset.get< bool>("SaveTrackInfo",false))
170  {
171  if (fSaveTrackInfo == false) fSaveCaloInfo = false;
172  }
std::string string
Definition: nybbler.cc:12
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
dune::protonanalysis::~protonanalysis ( )
virtual

Definition at line 175 of file protonanalysis_module.cc.

175  {
176  }

Member Function Documentation

void dune::protonanalysis::analyze ( const art::Event evt)

Definition at line 262 of file protonanalysis_module.cc.

262  {
263  reset();
264 
265 
266 
267  std::vector<art::Ptr<recob::Track> > tracklist;
268  auto trackListHandle = evt.getHandle< std::vector<recob::Track> >(fTrackModuleLabel);
269  if (trackListHandle) art::fill_ptr_vector(tracklist, trackListHandle);
270 
271  std::vector<art::Ptr<recob::Hit>> hitlist;
272  auto hitListHandle = evt.getHandle< std::vector<recob::Hit> >(fHitsModuleLabel); // to get information about the hits
273  if (hitListHandle) art::fill_ptr_vector(hitlist, hitListHandle);
274  art::FindManyP<recob::Hit> fmtht(trackListHandle, evt, fTrackModuleLabel); // to associate tracks and hits
275  art::FindManyP<anab::Calorimetry> fmcal(trackListHandle, evt, fCalorimetryModuleLabel);
276 
277  run = evt.run();
278  subrun = evt.subRun();
279  event = evt.id().event();
280  art::Timestamp ts = evt.time();
281  TTimeStamp tts(ts.timeHigh(), ts.timeLow());
282  evttime=tts.AsDouble();
283 
284  UInt_t year=0;
285  UInt_t month=0;
286  UInt_t day=0;
287 
288  year_month_date=tts.GetDate(kTRUE,0,&year,&month,&day);
289 
290  UInt_t hour=0;
291  UInt_t min=0;
292  UInt_t sec=0;
293 
294  hour_min_sec=tts.GetTime(kTRUE,0,&hour,&min,&sec);
295 
296  beam_proton=0;
297  total_trks=0;
298  tot_trks=0;
299 
300 
301 
302  auto const* SCE = lar::providerFrom<spacecharge::SpaceChargeService>();
303  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService>()->DataFor(evt);
304  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService>()->DataFor(evt, clockData);
305  std::cout<<detProp.Temperature()<<endl;
306  std::cout<<"drift velocity "<<detProp.DriftVelocity(0.50,87.0)<<std::endl;
307  double efield=detProp.Efield();
308  double ycord=300;
309  double zcord=350;
310  for(int i=0;i<120;i++){
311  double xcord=-360+i*6;
312  geo::Vector_t fEfieldOffsets = SCE->GetEfieldOffsets(geo::Point_t{xcord,ycord,zcord});
313  double Ef=efield+efield*fEfieldOffsets.X();
314  EFieldX->SetBinContent(i+1, Ef);
315  }
316 
317  for(int i=0;i<144;i++){
318  for(int j=0;j<120;j++){
319  for(int k=0;k<139;k++){
320  double x1=i*5.0-357.5;
321  double y1=j*5.0+2.5;
322  double z1=k*5.0+2.5;
323  geo::Vector_t fEfieldOffsets=SCE->GetEfieldOffsets(geo::Point_t{x1,y1,z1});
324  double EfX=efield+efield*fEfieldOffsets.X();
325  EField3DX->SetBinContent(i+1,j+1,k+1,EfX);
326  }
327  }
328  }
329 
330 
331 
332 
333 
334  double EField1=efield;
335  Float_t x1=-300;
336  for(int i=0;i<139;i++){
337  Float_t z=5*i+2.5;
338  for(int j=0;j<120;j++){
339  Float_t y=5*j+2.5;
340  geo::Vector_t fEfieldOffsets = SCE->GetEfieldOffsets(geo::Point_t{x1,y,z});
341  EField1 = std::sqrt((efield + efield*fEfieldOffsets.X())*(efield + efield*fEfieldOffsets.X())+(efield*fEfieldOffsets.Y())*(efield*fEfieldOffsets.Y())+(efield*fEfieldOffsets.Z())*(efield*fEfieldOffsets.Z()));
342  electric_field_YZ->SetBinContent(i+1,j+1,EField1);
343  }
344  }
345  Float_t x4=-25;
346  for(int i=0;i<139;i++){
347  Float_t z=5*i+2.5;
348  for(int j=0;j<120;j++){
349  Float_t y=5*j+2.5;
350  geo::Vector_t fEfieldOffsets = SCE->GetEfieldOffsets(geo::Point_t{x4,y,z});
351  EField1 = std::sqrt((efield + efield*fEfieldOffsets.X())*(efield + efield*fEfieldOffsets.X())+(efield*fEfieldOffsets.Y()*efield*fEfieldOffsets.Y())+(efield*fEfieldOffsets.Z()*efield*fEfieldOffsets.Z()));
352  electric_field_YZ_3->SetBinContent(i+1,j+1,EField1);
353  }
354  }
355 
356 
357  size_t NTracks = tracklist.size();
358  for(size_t i=0; i<NTracks;++i){
359  art::Ptr<recob::Track> ptrack(trackListHandle, i);
360  std::vector<art::Ptr<anab::Calorimetry>> calos=fmcal.at(i);
361  const recob::Track& track = *ptrack;
362  // TVector3 pos, dir_start, dir_end, end;
363  auto pos = track.Vertex();
364  auto dir_start = track.VertexDirection();
365  auto dir_end = track.EndDirection();
366  auto end = track.End();
367  double theta_xz = std::atan2(dir_start.X(), dir_start.Z());
368  double theta_yz = std::atan2(dir_start.Y(), dir_start.Z());
369 
370 
371 
372  /*********************storing hits for each tracks*******************************/
373 
374  std::vector<art::Ptr<recob::Hit>> allHits=fmtht.at(i); //storing hits for ith track
377  int trackid=-1;
378  std::map<int,double> trkide;
379  std::map<int,double> trknumelec;
380  std::vector<float> hitpeakT;
381  //selection based on hitpeakT
382  for(size_t h1=0;h1<allHits.size();h1++){
383  hitpeakT.push_back(allHits[h1]->PeakTime());
384  }
385  float max=-999999;
386  float min=-999999;
387  min=*min_element(hitpeakT.begin(),hitpeakT.end());
388  max=*max_element(hitpeakT.begin(),hitpeakT.end());
389  hitpeakT.clear();
390 
391  for(size_t h=0; h<allHits.size();h++){
392  art::Ptr<recob::Hit> hit=allHits[h];
393  std::vector<sim::TrackIDE> eveIDs = bt_serv->HitToTrackIDEs(clockData, hit);
394  for(size_t e=0;e<eveIDs.size(); ++e){
395  trkide[eveIDs[e].trackID] += eveIDs[e].energy;
396  }
397  }
398  double maxe = -1;
399  double tote = 0;
400  for(std::map<int,double>::iterator ii = trkide.begin(); ii!=trkide.end(); ++ii){
401  // cout<<" trkid "<<ii->first<<" energy deposited = "<<ii->second<<endl;
402  tote += ii->second;
403  if((ii->second)>maxe){
404  maxe = ii->second;
405  trackid = ii->first;
406 
407  }
408  }
409 
410  //*******************new section************************//
411  float total_energy=0.0;
412  for(size_t h=0; h<allHits.size();h++){
413  art::Ptr<recob::Hit> hit=allHits[h];
414  if (hit->WireID().Plane!=2) continue;
415  std::vector<sim::TrackIDE> eveIDs = bt_serv->HitToTrackIDEs(clockData, hit);
416  for(size_t e=0;e<eveIDs.size(); ++e){
417  if (eveIDs[e].trackID == trackid) total_energy += eveIDs[e].energy;
418  }
419  }
420 
421  int pdg1=-99999;
422 
423  const art::Ptr<simb::MCTruth> mc=pi_serv->TrackIdToMCTruth_P(trackid);
424  const simb::MCParticle *particle = pi_serv->TrackIdToParticle_P(trackid);
425  if(!particle) continue;
426  pdg1=particle->PdgCode();
427  int origin_type=mc->Origin();
428  // int no_daughters=particle->NumberDaughters();
429  if((particle->Process()=="primary") && origin_type==4) pdg_beamparticles->Fill(pdg1);
430 
431  if(particle->Process()=="primary" && origin_type==4){
432  tot_trks++;
433  pdg_all_beam[tot_trks-1]=pdg1;
434  }
435 
436  if(!(pdg1==2212 && (particle->Process()=="primary") && origin_type==4/* && no_daughters==0*/)) continue;//non interacting beam protons
437 
438  for(size_t h=0; h<allHits.size();h++){
439  art::Ptr<recob::Hit> hit=allHits[h];
440  std::vector<const sim::IDE*> eventIDs = bt_serv->HitToSimIDEs_Ps(clockData, hit);
441  float charge=0;
442  int number_electrons=0;
443  if(hit->WireID().Plane!=2) continue;
444  charge = hit->Integral();
445  for(size_t e=0;e<eventIDs.size(); ++e){
446  // trknumelec[eventIDs[e].trackID]+=eventIDs[e].numElectrons;
447  //std::cout<<eventIDs[e]->numElectrons<<" "<<eventIDs[e]->energy<<std::endl;
448  number_electrons+=eventIDs[e]->numElectrons;
449  }
450  calibration_constants->Fill(charge/float(number_electrons));
451  }
452 
453 
454 
455 
456 
457  total_trks++;
458 
459  beam_proton++;
460 
461 
462 
463 
464 
465  /* auto const* SCE = lar::providerFrom<spacecharge::SpaceChargeService>();
466 
467 
468  ****************************************************************************************************************/
469 
470  EndPointx_corrected[beam_proton-1]=particle->EndPosition()[0]-SCE->GetPosOffsets(geo::Point_t(particle->EndPosition()[0],particle->EndPosition()[1],particle->EndPosition()[2])).X();
471  EndPointy_corrected[beam_proton-1]=particle->EndPosition()[1]+SCE->GetPosOffsets(geo::Point_t(particle->EndPosition()[0],particle->EndPosition()[1],particle->EndPosition()[2])).Y();
472  EndPointz_corrected[beam_proton-1]=particle->EndPosition()[2]+SCE->GetPosOffsets(geo::Point_t(particle->EndPosition()[0],particle->EndPosition()[1],particle->EndPosition()[2])).Z();
473  origin[beam_proton-1]=mc->Origin();
474  EndE[beam_proton-1]=particle->EndE();
475  momentum[beam_proton-1]=particle->Momentum().Vect().Mag();
476  pdg[beam_proton-1]=particle->PdgCode();
477  TrackId_geant[beam_proton-1]=particle->TrackId();
478  EndPointx[beam_proton-1]=particle->EndPosition()[0];
479  EndPointy[beam_proton-1]=particle->EndPosition()[1];
480  EndPointz[beam_proton-1]=particle->EndPosition()[2];
481  StartPointx[beam_proton-1]=particle->Vx();
482  StartPointy[beam_proton-1]=particle->Vy();
483  StartPointz[beam_proton-1]=particle->Vz();
484  Total_energy[beam_proton-1]=total_energy;
485  trackthetaxz[beam_proton-1]=theta_xz;
486  trackthetayz[beam_proton-1]=theta_yz;
487  trkstartx[beam_proton-1]=pos.X();
488  trkstarty[beam_proton-1]=pos.Y();
489  trkstartz[beam_proton-1]=pos.Z();
490  trkendx[beam_proton-1]=end.X();
491  trkendy[beam_proton-1]=end.Y();
492  trkendz[beam_proton-1]=end.Z();
493  trklen[beam_proton-1]=track.Length();
496  TrkID[beam_proton-1]=track.ID();
497  trkstartcosxyz[beam_proton-1][0]=dir_start.X();
498  trkstartcosxyz[beam_proton-1][1]=dir_start.Y();
499  trkstartcosxyz[beam_proton-1][2]=dir_start.Z();
500  trkendcosxyz[beam_proton-1][0]=dir_end.X();
501  trkendcosxyz[beam_proton-1][1]=dir_end.Y();
502  trkendcosxyz[beam_proton-1][2]=dir_end.Z();
503  for(size_t ical = 0; ical<calos.size(); ++ical){
504  if(!calos[ical]) continue;
505  if(!calos[ical]->PlaneID().isValid) continue;
506  int planenum = calos[ical]->PlaneID().Plane;
507  if(planenum<0||planenum>2) continue;
508  const size_t NHits = calos[ical] -> dEdx().size();
509  ntrkhits[beam_proton-1][planenum]=int(NHits);
510  for(size_t iHit = 0; iHit < NHits; ++iHit){
511  const auto& TrkPos = (calos[ical] -> XYZ())[iHit];
512  trkdqdx[beam_proton-1][planenum][iHit]=(calos[ical] -> dQdx())[iHit];
513  trkdedx[beam_proton-1][planenum][iHit]=(calos[ical] -> dEdx())[iHit];
514  trkresrange[beam_proton-1][planenum][iHit]=(calos[ical]->ResidualRange())[iHit];
515  trkhitx[beam_proton-1][planenum][iHit]=TrkPos.X();
516  trkhity[beam_proton-1][planenum][iHit]=TrkPos.Y();
517  trkhitz[beam_proton-1][planenum][iHit]=TrkPos.Z();
518  trkpitch[beam_proton-1][planenum][iHit]=(calos[ical]->TrkPitchVec())[iHit];
519 
520  geo::Vector_t fEfieldOffsets = SCE->GetEfieldOffsets(geo::Point_t{TrkPos.X(),TrkPos.Y(),TrkPos.Z()});
521  Efield[beam_proton-1][planenum][iHit]=sqrt((efield + efield*fEfieldOffsets.X())*(efield + efield*fEfieldOffsets.X())+(efield*fEfieldOffsets.Y())*(efield*fEfieldOffsets.Y())+(efield*fEfieldOffsets.Z())*(efield*fEfieldOffsets.Z()));
522  } // loop over iHit..
523  } // loop over ical 2nd time...
524  } // loop over trks...
525  // stopping_muons->Fill(stopping_trks);
526  fEventTree->Fill();
527  } // end of analyze function
end
while True: pbar.update(maxval-len(onlies[E][S])) #print iS, "/", len(onlies[E][S]) found = False for...
intermediate_table::iterator iterator
int PdgCode() const
Definition: MCParticle.h:212
Double_t Total_energy[kMaxTracks]
constexpr std::uint32_t timeLow() const
Definition: Timestamp.h:29
double EndE() const
Definition: MCParticle.h:244
const TLorentzVector & EndPosition() const
Definition: MCParticle.h:225
std::vector< sim::TrackIDE > HitToTrackIDEs(detinfo::DetectorClocksData const &clockData, recob::Hit const &hit) const
Float_t trackthetayz[kMaxTracks]
Float_t EndPointx[kMaxTracks]
Float_t trkendcosxyz[kMaxTracks][3]
Handle< PROD > getHandle(SelectorBase const &) const
Definition: DataViewImpl.h:382
const simb::MCParticle * TrackIdToParticle_P(int id) const
geo::WireID WireID() const
Definition: Hit.h:233
Float_t trkendx[kMaxTracks]
Int_t TrackId_geant[kMaxTracks]
Float_t trkendy[kMaxTracks]
constexpr std::uint32_t timeHigh() const
Definition: Timestamp.h:34
Float_t StartPointy[kMaxTracks]
simb::Origin_t Origin() const
Definition: MCTruth.h:74
Float_t peakT_min[kMaxTracks]
Float_t trkstartcosxyz[kMaxTracks][3]
Vector_t VertexDirection() const
Definition: Track.h:132
Float_t trkhitz[kMaxTracks][3][3000]
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:224
Int_t pdg_all_beam[kMaxTracks]
Float_t StartPointz[kMaxTracks]
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Vector_t
Type for representation of momenta in 3D space.
Definition: geo_vectors.h:164
Float_t EndPointx_corrected[kMaxTracks]
std::string Process() const
Definition: MCParticle.h:215
int TrackId() const
Definition: MCParticle.h:210
Float_t peakT_max[kMaxTracks]
Float_t Efield[kMaxTracks][3][3000]
double dEdx(float dqdx, float Efield)
Definition: doAna.cpp:21
const double e
double Length(size_t p=0) const
Access to various track properties.
Definition: Track.h:167
Timestamp time() const
Double_t momentum[kMaxTracks]
Float_t trkstarty[kMaxTracks]
const art::Ptr< simb::MCTruth > & TrackIdToMCTruth_P(int id) const
Float_t trkdqdx[kMaxTracks][3][3000]
Point_t const & Vertex() const
Definition: Track.h:124
Float_t EndPointz[kMaxTracks]
Float_t trkpitch[kMaxTracks][3][3000]
SubRunNumber_t subRun() const
Definition: DataViewImpl.cc:78
Float_t StartPointx[kMaxTracks]
Float_t trkendz[kMaxTracks]
RunNumber_t run() const
Definition: DataViewImpl.cc:71
static int max(int a, int b)
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:493
Detector simulation of raw signals on wires.
Float_t EndPointy_corrected[kMaxTracks]
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:184
Float_t trkresrange[kMaxTracks][3][3000]
double Vx(const int i=0) const
Definition: MCParticle.h:221
Float_t trkstartx[kMaxTracks]
int ID() const
Definition: Track.h:198
Double_t EndE[kMaxTracks]
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
Vector_t EndDirection() const
Definition: Track.h:133
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:220
Float_t EndPointy[kMaxTracks]
double Vz(const int i=0) const
Definition: MCParticle.h:223
EventNumber_t event() const
Definition: EventID.h:116
Point_t const & End() const
Definition: Track.h:125
detail::Node< FrameID, bool > PlaneID
Definition: CRTID.h:125
Float_t trackthetaxz[kMaxTracks]
Float_t EndPointz_corrected[kMaxTracks]
std::vector< const sim::IDE * > HitToSimIDEs_Ps(detinfo::DetectorClocksData const &clockData, recob::Hit const &hit) const
Float_t trkhitx[kMaxTracks][3][3000]
Float_t trkdedx[kMaxTracks][3][3000]
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
Float_t trklen[kMaxTracks]
Float_t trkstartz[kMaxTracks]
Float_t trkhity[kMaxTracks][3][3000]
Int_t ntrkhits[kMaxTracks][3]
EventID id() const
Definition: Event.cc:34
double Vy(const int i=0) const
Definition: MCParticle.h:222
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track:
Definition: Track.h:49
QTextStream & endl(QTextStream &s)
void dune::protonanalysis::beginJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 180 of file protonanalysis_module.cc.

180  {
181  std::cout<<"job begin..."<<std::endl;
183  pdg_beamparticles=tfs->make<TH1F>("beam_pdg","beam particle pdgs;pdg;no of beam particles",6000,-2999.5,3000.5);
184  calibration_constants=tfs->make<TH1F>("calibration_constants","calibration constants;calibration constants;no of hits",6000,0.003000,0.009000);
185  electric_field_YZ=tfs->make<TH2F>("EField_YZ_x_0","Electric field for YZ plane at x= -300cm;Z coordinate;Y coordinate",139,0,695,120,0,600);
186  electric_field_YZ_3=tfs->make<TH2F>("EField_YZ_x_3","Electric field for YZ plane at x= -25cm;Z coordinate;Y coordinate",139,0,695,120,0,600);
187  EFieldX=tfs->make<TH1F>("EFieldX","Electric Field for Y=300cm && Z=350cm;X coordinate;Efield in kV/cm",120,-360,360);
188  EField3DX=tfs->make<TH3F>("EField3DX","Electric FieldX 3D",144,-360,360,120,0,600,139,0,695);
189 
190  pdg_beamparticles->Sumw2();
191  calibration_constants->Sumw2();
192  fEventTree = tfs->make<TTree>("Event", "Event Tree from Reco");
193  fEventTree->Branch("event", &event,"event/I");
194 
195 
196  fEventTree->Branch("evttime",&evttime,"evttime/D");
197  fEventTree->Branch("run", &run,"run/I");
198  fEventTree->Branch("subrun", &subrun,"surbrun/I");
199  fEventTree->Branch("year_month_date", &year_month_date,"year_month_date/I");
200  fEventTree->Branch("hour_min_sec", &hour_min_sec,"hour_min_sec/I");
201  fEventTree->Branch("beam_proton",&beam_proton,"beam_proton/I");//total tracks pasing broken tracks filter
202  fEventTree->Branch("tot_trks",&tot_trks,"tot_trks/I");
203 
204  fEventTree->Branch("total_trks",&total_trks,"total_trks/I");//total proton beam tracks
205  fEventTree->Branch("xprojectedlen",xprojectedlen,"xprojectedlen[beam_proton]/F");
206  fEventTree->Branch("trackthetaxz",trackthetaxz,"trackthetaxz[beam_proton]/F");
207  fEventTree->Branch("trackthetayz",trackthetayz,"trackthetayz[beam_proton]/F");
208  fEventTree->Branch("trkstartx",trkstartx,"trkstartx[beam_proton]/F");
209  fEventTree->Branch("trkstarty",trkstarty,"trkstarty[beam_proton]/F");
210  fEventTree->Branch("trkstartz",trkstartz,"trkstartz[beam_proton]/F");
211  fEventTree->Branch("trkendx",trkendx,"trkendx[beam_proton]/F");
212  fEventTree->Branch("trkendy",trkendy,"trkendy[beam_proton]/F");
213  fEventTree->Branch("trkendz",trkendz,"trkendz[beam_proton]/F");
214  fEventTree->Branch("trklen",trklen,"trklen[beam_proton]/F");
215  fEventTree->Branch("Total_energy",Total_energy,"Total_energy[beam_proton]/D");
216  fEventTree->Branch("EndE",EndE,"EndE[beam_proton]/D");
217  fEventTree->Branch("momentum",momentum,"momentum[beam_proton]/D");
218  fEventTree->Branch("peakT_max",peakT_max,"peakT_max[beam_proton]/F");
219  fEventTree->Branch("peakT_min",peakT_min,"peakT_min[beam_proton]/F");
220  fEventTree->Branch("pdg",pdg,"pdg[beam_proton]/I");
221  fEventTree->Branch("pdg_all_beam",pdg_all_beam,"pdg_all_beam[tot_trks]/I");
222  fEventTree->Branch("origin",origin,"origin[beam_proton]/I");
223  fEventTree->Branch("TrkID",TrkID,"TrkID[beam_proton]/I");
224  fEventTree->Branch("TrackId_geant",TrackId_geant,"TrackId_geant[beam_proton]/I");
225  fEventTree->Branch("EndPointx",EndPointx,"EndPointx[beam_proton]/F");
226  fEventTree->Branch("EndPointy",EndPointy,"EndPointy[beam_proton]/F");
227  fEventTree->Branch("EndPointz",EndPointz,"EndPointz[beam_proton]/F");
228  fEventTree->Branch("EndPointx_corrected",EndPointx_corrected,"EndPointx_corrected[beam_proton]/F");
229  fEventTree->Branch("EndPointy_corrected",EndPointy_corrected,"EndPointy_corrected[beam_proton]/F");
230  fEventTree->Branch("EndPointz_corrected",EndPointz_corrected,"EndPointz_corrected[beam_proton]/F");
231  fEventTree->Branch("StartPointx",StartPointx,"StartPointx[beam_proton]/F");
232  fEventTree->Branch("StartPointy",StartPointy,"StartPointy[beam_proton]/F");
233  fEventTree->Branch("StartPointz",StartPointz,"StartPointz[beam_proton]/F");
234  fEventTree->Branch("trkstartcosxyz",trkstartcosxyz,"trkstartcosxyz[beam_proton][3]/F");
235  fEventTree->Branch("trkendcosxyz",trkendcosxyz,"trkendcosxyz[beam_proton][3]/F");
236  fEventTree->Branch("ntrkhits",ntrkhits,"ntrkhits[beam_proton][3]/I");
237  fEventTree->Branch("trkdqdx",trkdqdx,"trkdqdx[beam_proton][3][3000]/F");
238  fEventTree->Branch("trkdedx",trkdedx,"trkdedx[beam_proton][3][3000]/F");
239  fEventTree->Branch("trkresrange",trkresrange,"trkresrange[beam_proton][3][3000]/F");
240  fEventTree->Branch("trkhitx",trkhitx,"trkhitx[beam_proton][3][3000]/F");
241  fEventTree->Branch("trkhity",trkhity,"trkhity[beam_proton][3][3000]/F");
242  fEventTree->Branch("trkhitz",trkhitz,"trkhitz[beam_proton][3][3000]/F");
243  fEventTree->Branch("trkpitch",trkpitch,"trkpitch[beam_proton][3][3000]/F");
244  fEventTree->Branch("Efield",Efield,"Efield[beam_proton][3][3000]/F");
245  }
Double_t Total_energy[kMaxTracks]
Float_t trackthetayz[kMaxTracks]
Float_t EndPointx[kMaxTracks]
Float_t trkendcosxyz[kMaxTracks][3]
Float_t trkendx[kMaxTracks]
Int_t TrackId_geant[kMaxTracks]
Float_t trkendy[kMaxTracks]
Float_t StartPointy[kMaxTracks]
Float_t peakT_min[kMaxTracks]
Float_t trkstartcosxyz[kMaxTracks][3]
Float_t trkhitz[kMaxTracks][3][3000]
Int_t pdg_all_beam[kMaxTracks]
Float_t StartPointz[kMaxTracks]
Float_t EndPointx_corrected[kMaxTracks]
Float_t peakT_max[kMaxTracks]
Float_t Efield[kMaxTracks][3][3000]
Double_t momentum[kMaxTracks]
Float_t trkstarty[kMaxTracks]
Float_t trkdqdx[kMaxTracks][3][3000]
Float_t EndPointz[kMaxTracks]
Float_t trkpitch[kMaxTracks][3][3000]
Float_t StartPointx[kMaxTracks]
Float_t trkendz[kMaxTracks]
Float_t EndPointy_corrected[kMaxTracks]
Float_t trkresrange[kMaxTracks][3][3000]
Float_t trkstartx[kMaxTracks]
Double_t EndE[kMaxTracks]
Float_t EndPointy[kMaxTracks]
Float_t trackthetaxz[kMaxTracks]
Float_t EndPointz_corrected[kMaxTracks]
Float_t trkhitx[kMaxTracks][3][3000]
Float_t trkdedx[kMaxTracks][3][3000]
Float_t trklen[kMaxTracks]
Float_t trkstartz[kMaxTracks]
Float_t trkhity[kMaxTracks][3][3000]
Int_t ntrkhits[kMaxTracks][3]
QTextStream & endl(QTextStream &s)
Event finding and building.
Float_t xprojectedlen[kMaxTracks]
void dune::protonanalysis::beginRun ( const art::Run run)

Definition at line 253 of file protonanalysis_module.cc.

253  {
254  mf::LogInfo("protonanalysis")<<"begin run..."<<std::endl;
255  }
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
QTextStream & endl(QTextStream &s)
void dune::protonanalysis::endJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 248 of file protonanalysis_module.cc.

248  {
249 
250  }
void dune::protonanalysis::reset ( )

Definition at line 530 of file protonanalysis_module.cc.

530  {
531  run = -99999;
532  subrun = -99999;
533  event = -99999;
534  evttime = -99999;
535  beam_proton = -99999;
536  total_trks=-99999;
537  year_month_date=-99999;
538  hour_min_sec=-99999;
539  for(int i=0; i<kMaxTracks; i++){
540  trackthetaxz[i]=-99999;
541  trackthetayz[i]=-99999;
542  trkstartx[i]=-99999;
543  trkstarty[i]=-99999;
544  trkstartz[i]=-99999;
545  trkendx[i]=-99999;
546  trkendy[i]=-99999;
547  trkendz[i]=-99999;
548  trklen[i]=-99999;
549  TrkID[i]=-99999;
550  pdg_all_beam[i]=-99999;
551  TrackId_geant[i]=-99999;
552  origin[i]=-99999;
553  EndE[i]=-99999;
554  pdg[i]=-99999;
555  Total_energy[i]=-999999;
556  EndPointx[i]=-99999;
557  EndPointy[i]=-99999;
558  EndPointz[i]=-99999;
559  EndPointx_corrected[i]=-99999;
560  EndPointy_corrected[i]=-99999;
561  EndPointz_corrected[i]=-99999;
562  StartPointx[i]=-99999;
563  StartPointy[i]=-99999;
564  StartPointz[i]=-99999;
565  pdg[i]=-99999;
566  EndE[i]=-99999;
567  momentum[i]=-99999;
568  peakT_max[i]=-99999;
569  peakT_min[i]=-99999;
570  xprojectedlen[i]=-99999;
571  trkstartcosxyz[i][0]=-99999;
572  trkstartcosxyz[i][1]=-99999;
573  trkstartcosxyz[i][2]=-99999;
574  trkendcosxyz[i][0]=-99999;
575  trkendcosxyz[i][1]=-99999;
576  trkendcosxyz[i][2]=-99999;
577  ntrkhits[i][0] = -99999;
578  ntrkhits[i][1] = -99999;
579  ntrkhits[i][2] = -99999;
580  for(int j=0; j<3; j++){
581  for(int k=0; k<3000; k++){
582  trkdqdx[i][j][k]=-99999;
583  trkdedx[i][j][k]=-99999;
584  trkresrange[i][j][k]=-99999;
585  trkhitx[i][j][k]=-99999;
586  trkhity[i][j][k]=-99999;
587  trkhitz[i][j][k]=-99999;
588  trkpitch[i][j][k]=-99999;
589  Efield[i][j][k]=-99999;
590  }
591  }
592  }
593  }
Double_t Total_energy[kMaxTracks]
Float_t trackthetayz[kMaxTracks]
Float_t EndPointx[kMaxTracks]
Float_t trkendcosxyz[kMaxTracks][3]
Float_t trkendx[kMaxTracks]
Int_t TrackId_geant[kMaxTracks]
Float_t trkendy[kMaxTracks]
Float_t StartPointy[kMaxTracks]
Float_t peakT_min[kMaxTracks]
Float_t trkstartcosxyz[kMaxTracks][3]
Float_t trkhitz[kMaxTracks][3][3000]
Int_t pdg_all_beam[kMaxTracks]
Float_t StartPointz[kMaxTracks]
Float_t EndPointx_corrected[kMaxTracks]
const int kMaxTracks
Float_t peakT_max[kMaxTracks]
Float_t Efield[kMaxTracks][3][3000]
Double_t momentum[kMaxTracks]
Float_t trkstarty[kMaxTracks]
Float_t trkdqdx[kMaxTracks][3][3000]
Float_t EndPointz[kMaxTracks]
Float_t trkpitch[kMaxTracks][3][3000]
Float_t StartPointx[kMaxTracks]
Float_t trkendz[kMaxTracks]
Float_t EndPointy_corrected[kMaxTracks]
Float_t trkresrange[kMaxTracks][3][3000]
Float_t trkstartx[kMaxTracks]
Double_t EndE[kMaxTracks]
Float_t EndPointy[kMaxTracks]
Float_t trackthetaxz[kMaxTracks]
Float_t EndPointz_corrected[kMaxTracks]
Float_t trkhitx[kMaxTracks][3][3000]
Float_t trkdedx[kMaxTracks][3][3000]
Float_t trklen[kMaxTracks]
Float_t trkstartz[kMaxTracks]
Float_t trkhity[kMaxTracks][3][3000]
Int_t ntrkhits[kMaxTracks][3]
Float_t xprojectedlen[kMaxTracks]

Member Data Documentation

Int_t dune::protonanalysis::beam_proton
private

Definition at line 114 of file protonanalysis_module.cc.

TH1F* dune::protonanalysis::calibration_constants
private

Definition at line 98 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::Efield[kMaxTracks][3][3000]
private

Definition at line 150 of file protonanalysis_module.cc.

TH3F* dune::protonanalysis::EField3DX
private

Definition at line 102 of file protonanalysis_module.cc.

TH1F* dune::protonanalysis::EFieldX
private

Definition at line 101 of file protonanalysis_module.cc.

TH2F* dune::protonanalysis::electric_field_YZ
private

Definition at line 99 of file protonanalysis_module.cc.

TH2F* dune::protonanalysis::electric_field_YZ_3
private

Definition at line 100 of file protonanalysis_module.cc.

Double_t dune::protonanalysis::EndE[kMaxTracks]
private

Definition at line 109 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::EndPointx[kMaxTracks]
private

Definition at line 131 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::EndPointx_corrected[kMaxTracks]
private

Definition at line 134 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::EndPointy[kMaxTracks]
private

Definition at line 132 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::EndPointy_corrected[kMaxTracks]
private

Definition at line 135 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::EndPointz[kMaxTracks]
private

Definition at line 133 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::EndPointz_corrected[kMaxTracks]
private

Definition at line 136 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::event
private

Definition at line 106 of file protonanalysis_module.cc.

Double_t dune::protonanalysis::evttime
private

Definition at line 111 of file protonanalysis_module.cc.

std::string dune::protonanalysis::fCalorimetryModuleLabel
private

Definition at line 157 of file protonanalysis_module.cc.

TTree* dune::protonanalysis::fEventTree
private

Definition at line 103 of file protonanalysis_module.cc.

std::string dune::protonanalysis::fHitsModuleLabel
private

Definition at line 155 of file protonanalysis_module.cc.

bool dune::protonanalysis::fSaveCaloInfo
private

Definition at line 158 of file protonanalysis_module.cc.

bool dune::protonanalysis::fSaveTrackInfo
private

Definition at line 159 of file protonanalysis_module.cc.

std::string dune::protonanalysis::fTrackModuleLabel
private

Definition at line 156 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::hour_min_sec
private

Definition at line 113 of file protonanalysis_module.cc.

Double_t dune::protonanalysis::momentum[kMaxTracks]
private

Definition at line 153 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::ntrkhits[kMaxTracks][3]
private

Definition at line 142 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::origin[kMaxTracks]
private

Definition at line 130 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::pdg[kMaxTracks]
private

Definition at line 110 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::pdg_all_beam[kMaxTracks]
private

Definition at line 117 of file protonanalysis_module.cc.

TH1F* dune::protonanalysis::pdg_beamparticles
private

Definition at line 97 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::peakT_max[kMaxTracks]
private

Definition at line 151 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::peakT_min[kMaxTracks]
private

Definition at line 152 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::run
private

Definition at line 104 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::StartPointx[kMaxTracks]
private

Definition at line 137 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::StartPointy[kMaxTracks]
private

Definition at line 138 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::StartPointz[kMaxTracks]
private

Definition at line 139 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::subrun
private

Definition at line 105 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::tot_trks
private

Definition at line 115 of file protonanalysis_module.cc.

Double_t dune::protonanalysis::Total_energy[kMaxTracks]
private

Definition at line 108 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::total_trks
private

Definition at line 116 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::TrackId_geant[kMaxTracks]
private

Definition at line 129 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trackthetaxz[kMaxTracks]
private

Definition at line 119 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trackthetayz[kMaxTracks]
private

Definition at line 120 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkdedx[kMaxTracks][3][3000]
private

Definition at line 144 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkdqdx[kMaxTracks][3][3000]
private

Definition at line 143 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkendcosxyz[kMaxTracks][3]
private

Definition at line 141 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkendx[kMaxTracks]
private

Definition at line 124 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkendy[kMaxTracks]
private

Definition at line 125 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkendz[kMaxTracks]
private

Definition at line 126 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkhitx[kMaxTracks][3][3000]
private

Definition at line 146 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkhity[kMaxTracks][3][3000]
private

Definition at line 147 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkhitz[kMaxTracks][3][3000]
private

Definition at line 148 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::TrkID[kMaxTracks]
private

Definition at line 128 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trklen[kMaxTracks]
private

Definition at line 127 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkpitch[kMaxTracks][3][3000]
private

Definition at line 149 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkresrange[kMaxTracks][3][3000]
private

Definition at line 145 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkstartcosxyz[kMaxTracks][3]
private

Definition at line 140 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkstartx[kMaxTracks]
private

Definition at line 121 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkstarty[kMaxTracks]
private

Definition at line 122 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::trkstartz[kMaxTracks]
private

Definition at line 123 of file protonanalysis_module.cc.

Float_t dune::protonanalysis::xprojectedlen[kMaxTracks]
private

Definition at line 118 of file protonanalysis_module.cc.

Int_t dune::protonanalysis::year_month_date
private

Definition at line 112 of file protonanalysis_module.cc.


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