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

Public Member Functions

 NCDeltaRadiative (fhicl::ParameterSet const &p)
 
 NCDeltaRadiative (NCDeltaRadiative const &)=delete
 
 NCDeltaRadiative (NCDeltaRadiative &&)=delete
 
NCDeltaRadiativeoperator= (NCDeltaRadiative const &)=delete
 
NCDeltaRadiativeoperator= (NCDeltaRadiative &&)=delete
 
void cout_stuff (art::Event &e, bool passed)
 
void FillTree (art::Event &e, size_t const mct_index, size_t const parent_index, bool const is_nc_delta_radiative)
 
void Reset ()
 
bool filter (art::Event &e) override
 

Private Attributes

TTree * ftree
 
int frun
 
int fsubrun
 
int fevent
 
int fnu_pdg
 
int fccnc
 
int fmode
 
int finteraction_type
 
int fis_nc_delta_radiative
 
int fparent_status_code
 
int fparent_pdg
 

Detailed Description

Definition at line 30 of file NCDeltaRadiative_module.cc.

Constructor & Destructor Documentation

NCDeltaRadiative::NCDeltaRadiative ( fhicl::ParameterSet const &  p)
explicit

Definition at line 64 of file NCDeltaRadiative_module.cc.

64  :
65  art::EDFilter(p),
66  ftree(nullptr) {
67 
68  art::ServiceHandle<art::TFileService> tfs;
69  ftree = tfs->make<TTree>("NCDeltaRadiativeFilter", "");
70 
71  ftree->Branch("run", &frun, "run/I");
72  ftree->Branch("subrun", &fsubrun, "subrun/I");
73  ftree->Branch("event", &fevent, "event/I");
74  ftree->Branch("nu_pdg", &fnu_pdg, "nu_pdg/I");
75  ftree->Branch("ccnc", &fccnc, "ccnc/I");
76  ftree->Branch("mode", &fmode, "mode/I");
77  ftree->Branch("is_nc_delta_radiative", &fis_nc_delta_radiative, "is_nc_delta_radiative/I");
78  ftree->Branch("parent_status_code", &fparent_status_code, "parent_status_code/I");
79  ftree->Branch("parent_pdg", &fparent_pdg, "parent_pdg/I");
80 
81 }
pdgs p
Definition: selectors.fcl:22
art::ServiceHandle< art::TFileService > tfs
NCDeltaRadiative::NCDeltaRadiative ( NCDeltaRadiative const &  )
delete
NCDeltaRadiative::NCDeltaRadiative ( NCDeltaRadiative &&  )
delete

Member Function Documentation

void NCDeltaRadiative::cout_stuff ( art::Event &  e,
bool  passed = false 
)

Definition at line 84 of file NCDeltaRadiative_module.cc.

84  {
85 
86  art::ValidHandle<std::vector<simb::MCTruth>> const & ev_mct =
87  e.getValidHandle<std::vector<simb::MCTruth>>("generator");
88 
89  std::cout << passed << "\n"
90  << "==========================\n";
91  for(simb::MCTruth const & mct : *ev_mct) {
92  std::cout << "----------------------------\n";
93  for(int i = 0; i < mct.NParticles(); ++i) {
94  simb::MCParticle const & mcp = mct.GetParticle(i);
95  std::cout << mcp.TrackId() << " " << mcp.PdgCode() << " " << mcp.Mother() << " " << mcp.StatusCode() << "\n";
96  }
97  }
98 
99 }
do i e
BEGIN_PROLOG could also be cout
void NCDeltaRadiative::FillTree ( art::Event &  e,
size_t const  mct_index,
size_t const  parent_index,
bool const  is_nc_delta_radiative 
)

Definition at line 118 of file NCDeltaRadiative_module.cc.

121  {
122 
123  Reset();
124 
125  frun = e.id().run();
126  fsubrun = e.id().subRun();
127  fevent = e.id().event();
128 
129  art::ValidHandle<std::vector<simb::MCTruth>> const & ev_mct =
130  e.getValidHandle<std::vector<simb::MCTruth>>("generator");
131  simb::MCNeutrino const & mcn = ev_mct->at(mct_index).GetNeutrino();
132 
133  fnu_pdg = mcn.Nu().PdgCode();
134  fccnc = mcn.CCNC();
135  fmode = mcn.Mode();
136  finteraction_type = mcn.InteractionType();
137 
138  fis_nc_delta_radiative = is_nc_delta_radiative;
139  if(parent_index != SIZE_MAX) {
140  fparent_status_code = ev_mct->at(mct_index).GetParticle(parent_index).StatusCode();
141  fparent_pdg = ev_mct->at(mct_index).GetParticle(parent_index).PdgCode();
142  }
143 
144  ftree->Fill();
145 
146 }
do i e
bool NCDeltaRadiative::filter ( art::Event &  e)
override

Definition at line 149 of file NCDeltaRadiative_module.cc.

149  {
150 
151  art::ValidHandle<std::vector<simb::MCTruth>> const & ev_mct =
152  e.getValidHandle<std::vector<simb::MCTruth>>("generator");
153 
154  for(size_t i = 0; i < ev_mct->size(); ++i) {
155 
156  simb::MCTruth const & mct = ev_mct->at(i);
157  if(mct.GetNeutrino().CCNC() != 1) continue;
158 
159  std::vector<size_t> exiting_photon_parents;
160  for(int i = 0; i < mct.NParticles(); ++i) {
161  simb::MCParticle const & mcp = mct.GetParticle(i);
162  if(mcp.TrackId() != i) {
163  std::cout << "ERROR: " << __LINE__ << " " << __PRETTY_FUNCTION__ << "\nTrackId does not match index\n";
164  exit(1);
165  }
166  if(!(mcp.StatusCode() == 1 && mcp.PdgCode() == 22)) continue;
167  exiting_photon_parents.push_back(mcp.Mother());
168  }
169 
170  std::vector<size_t> in_nucleus_photons;
171  for(size_t const s : exiting_photon_parents) {
172  simb::MCParticle const & mcp = mct.GetParticle(s);
173  //CHEKC hardcode, TPC filter:
174  if(abs(mcp.Vx())>210 || abs(mcp.Vy())>210||mcp.Vz()>510 || mcp.Vz()<-1){
175  std::cout<<"OUTSIDE TPC x y z ="<<mcp.Vx()<<" "<<mcp.Vy()<<" "<<mcp.Vz()<<std::endl;
176  exit(0);
177  }
178  if(abs(mcp.PdgCode()) == 2114 || abs(mcp.PdgCode()) == 2214) {
179  if(ftree) FillTree(e, i, mcp.PdgCode(), true);
180  return true;
181  }
182  else if(mcp.PdgCode() == 22) {
183  in_nucleus_photons.push_back(mcp.Mother());
184  }
185  }
186 
187  for(size_t const s : in_nucleus_photons) {
188  simb::MCParticle const & mcp = mct.GetParticle(s);
189  if(abs(mcp.PdgCode()) == 2114 || abs(mcp.PdgCode()) == 2214) {
190  if(ftree) FillTree(e, i, s, true);
191  return true;
192  }
193  }
194 
195  }
196 
197  if(ftree) FillTree(e, 0, SIZE_MAX, false);
198  return false;
199 
200 }
void FillTree(art::Event &e, size_t const mct_index, size_t const parent_index, bool const is_nc_delta_radiative)
T abs(T value)
then echo File list $list not found else cat $list while read file do echo $file sed s
Definition: file_to_url.sh:60
do i e
BEGIN_PROLOG could also be cout
NCDeltaRadiative& NCDeltaRadiative::operator= ( NCDeltaRadiative const &  )
delete
NCDeltaRadiative& NCDeltaRadiative::operator= ( NCDeltaRadiative &&  )
delete
void NCDeltaRadiative::Reset ( )

Member Data Documentation

int NCDeltaRadiative::fccnc
private

Definition at line 38 of file NCDeltaRadiative_module.cc.

int NCDeltaRadiative::fevent
private

Definition at line 36 of file NCDeltaRadiative_module.cc.

int NCDeltaRadiative::finteraction_type
private

Definition at line 40 of file NCDeltaRadiative_module.cc.

int NCDeltaRadiative::fis_nc_delta_radiative
private

Definition at line 41 of file NCDeltaRadiative_module.cc.

int NCDeltaRadiative::fmode
private

Definition at line 39 of file NCDeltaRadiative_module.cc.

int NCDeltaRadiative::fnu_pdg
private

Definition at line 37 of file NCDeltaRadiative_module.cc.

int NCDeltaRadiative::fparent_pdg
private

Definition at line 43 of file NCDeltaRadiative_module.cc.

int NCDeltaRadiative::fparent_status_code
private

Definition at line 42 of file NCDeltaRadiative_module.cc.

int NCDeltaRadiative::frun
private

Definition at line 34 of file NCDeltaRadiative_module.cc.

int NCDeltaRadiative::fsubrun
private

Definition at line 35 of file NCDeltaRadiative_module.cc.

TTree* NCDeltaRadiative::ftree
private

Definition at line 32 of file NCDeltaRadiative_module.cc.


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