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

Public Member Functions

 ToFFilter (fhicl::ParameterSet const &p)
 
 ToFFilter (ToFFilter const &)=delete
 
 ToFFilter (ToFFilter &&)=delete
 
ToFFilteroperator= (ToFFilter const &)=delete
 
ToFFilteroperator= (ToFFilter &&)=delete
 
bool filter (art::Event &evt) override
 

Private Attributes

art::InputTag ftofLhitLabel
 
art::InputTag ftofChitLabel
 
art::InputTag ftofLflashLabel
 
art::InputTag ftofCflashLabel
 
art::InputTag ftofLflashhitLabel
 
art::InputTag ftofCflashhitLabel
 
bool fuse_Lhit
 
bool fuse_Chit
 
bool fuse_Lflsh
 
bool fuse_Cflsh
 
bool fuse_Lflsh_hit
 
bool fuse_Cflsh_hit
 
float ftof_Lhit_cut
 
float ftof_Chit_cut
 
float ftof_Lflsh_cut
 
float ftof_Cflsh_cut
 
float ftof_Lflshhit_cut
 
float ftof_Cflshhit_cut
 

Detailed Description

Definition at line 116 of file ToFFilter_module.cc.

Constructor & Destructor Documentation

sbnd::ToFFilter::ToFFilter ( fhicl::ParameterSet const &  p)
explicit

Definition at line 149 of file ToFFilter_module.cc.

149  :
150 EDFilter{pset},
151 ftofLhitLabel(pset.get<art::InputTag>("tofLhitLabel")),
152 ftofChitLabel(pset.get<art::InputTag>("tofChitLabel")),
153 ftofLflashLabel(pset.get<art::InputTag>("tofLflashLabel")),
154 ftofCflashLabel(pset.get<art::InputTag>("tofCflashLabel")),
155 ftofLflashhitLabel(pset.get<art::InputTag>("tofLflashhitLabel")),
156 ftofCflashhitLabel(pset.get<art::InputTag>("tofCflashhitLabel")),
157 fuse_Lhit(pset.get<bool>("use_Lhit")),
158 fuse_Chit(pset.get<bool>("use_Chit")),
159 fuse_Lflsh(pset.get<bool>("use_Lflsh")),
160 fuse_Cflsh(pset.get<bool>("use_Cflsh")),
161 fuse_Lflsh_hit(pset.get<bool>("use_Lflsh_hit")),
162 fuse_Cflsh_hit(pset.get<bool>("use_Cflsh_hit")),
163 ftof_Lhit_cut(pset.get<float>("tof_Lhit_cut")),
164 ftof_Chit_cut(pset.get<float>("tof_Chit_cut")),
165 ftof_Lflsh_cut(pset.get<float>("tof_Lflsh_cut")),
166 ftof_Cflsh_cut(pset.get<float>("tof_Cflsh_cut")),
167 ftof_Lflshhit_cut(pset.get<float>("tof_Lflshhit_cut")),
168 ftof_Cflshhit_cut(pset.get<float>("tof_Cflshhit_cut"))
169 {
170 }
art::InputTag ftofCflashLabel
art::InputTag ftofChitLabel
art::InputTag ftofCflashhitLabel
art::InputTag ftofLflashhitLabel
art::InputTag ftofLflashLabel
art::InputTag ftofLhitLabel
sbnd::ToFFilter::ToFFilter ( ToFFilter const &  )
delete
sbnd::ToFFilter::ToFFilter ( ToFFilter &&  )
delete

Member Function Documentation

bool sbnd::ToFFilter::filter ( art::Event &  evt)
override

Definition at line 172 of file ToFFilter_module.cc.

173 {
174  bool keep_event = true;
175 
176  // ================================== Calculatin ToF values using Largest optical hit method =========================
177 
178  if(fuse_Lhit){
179  art::Handle< std::vector<sbnd::ToF::ToF> > tofLhitListHandle;
180  std::vector< art::Ptr<sbnd::ToF::ToF> > tofLhitList;
181  if( evt.getByLabel(ftofLhitLabel,tofLhitListHandle))
182  art::fill_ptr_vector(tofLhitList, tofLhitListHandle);
183 
184  std::vector<bool> is_cos_like;
185 
186  for(auto const& tf_lhit : tofLhitList){
187  if(tf_lhit->tof >= ftof_Lhit_cut) is_cos_like.push_back(false);
188  else is_cos_like.push_back(true);
189  }
190 
191  if(is_cos_like.size()){
192  bool found_nu=false;
193  for(auto itr: is_cos_like){
194  if(!itr){
195  found_nu=true;
196  break;
197  }
198  }
199 
200  if(!found_nu) keep_event=false;
201  }
202  } // Using largest optical hit to calculate tof
203 
204  //======================================================================================================================
205 
206  // ================================== Calculatin ToF values using Closest optical hit method ===========================
207 
208  if(fuse_Chit){
209  art::Handle< std::vector<sbnd::ToF::ToF> > tofChitListHandle;
210  std::vector< art::Ptr<sbnd::ToF::ToF> > tofChitList;
211  if( evt.getByLabel(ftofChitLabel,tofChitListHandle))
212  art::fill_ptr_vector(tofChitList, tofChitListHandle);
213 
214  std::vector<bool> is_cos_like;
215 
216  for(auto const& tf_chit : tofChitList){
217  if(tf_chit->tof >= ftof_Chit_cut) is_cos_like.push_back(false);
218  else is_cos_like.push_back(true);
219  }
220 
221  if(is_cos_like.size()){
222  bool found_nu=false;
223  for(auto itr: is_cos_like){
224  if(!itr){
225  found_nu=true;
226  break;
227  }
228  }
229 
230  if(!found_nu) keep_event=false;
231  }
232  } // Using closest optical hit to calculate tof
233 
234  //======================================================================================================================
235 
236  //==================================== Calculation ToF values using Largest optical flash ==============================
237 
238  if(fuse_Lflsh){
239  art::Handle< std::vector<sbnd::ToF::ToF> > tofLflashListHandle;
240  std::vector< art::Ptr<sbnd::ToF::ToF> > tofLflashList;
241  if( evt.getByLabel(ftofLflashLabel,tofLflashListHandle))
242  art::fill_ptr_vector(tofLflashList, tofLflashListHandle);
243 
244  std::vector<bool> is_cos_like;
245 
246  for(auto const& tf_lflash : tofLflashList){
247  if(tf_lflash->tof >= ftof_Lflsh_cut) is_cos_like.push_back(false);
248  else is_cos_like.push_back(true);
249  }
250 
251  if(is_cos_like.size()){
252  bool found_nu=false;
253  for(auto itr: is_cos_like){
254  if(!itr){
255  found_nu=true;
256  break;
257  }
258  }
259 
260  if(!found_nu) keep_event=false;
261  }
262  } // Using largest optical flash to calculate tof
263 
264  //======================================================================================================================
265 
266  //==================================== Calculation ToF values using Largest optical flash ==============================
267 
268  if(fuse_Cflsh){
269  art::Handle< std::vector<sbnd::ToF::ToF> > tofCflashListHandle;
270  std::vector< art::Ptr<sbnd::ToF::ToF> > tofCflashList;
271  if( evt.getByLabel(ftofCflashLabel,tofCflashListHandle))
272  art::fill_ptr_vector(tofCflashList, tofCflashListHandle);
273 
274  std::vector<bool> is_cos_like;
275 
276  for(auto const& tf_cflash : tofCflashList){
277  if(tf_cflash->tof >= ftof_Cflsh_cut) is_cos_like.push_back(false);
278  else is_cos_like.push_back(true);
279  }
280 
281  if(is_cos_like.size()){
282  bool found_nu=false;
283  for(auto itr: is_cos_like){
284  if(!itr){
285  found_nu=true;
286  break;
287  }
288  }
289 
290  if(!found_nu) keep_event=false;
291  }
292  } // Using closest optical flash to calculate tof
293 
294  //======================================================================================================================
295 
296  //=========================Calculation ToF values using Earliest hit of the Largest flash ==============================
297 
298  if(fuse_Lflsh_hit){
299  art::Handle< std::vector<sbnd::ToF::ToF> > tofLflashhitListHandle;
300  std::vector< art::Ptr<sbnd::ToF::ToF> > tofLflashhitList;
301  if( evt.getByLabel(ftofLflashhitLabel,tofLflashhitListHandle))
302  art::fill_ptr_vector(tofLflashhitList, tofLflashhitListHandle);
303 
304  std::vector<bool> is_cos_like;
305 
306  for(auto const& tf_lflashhit : tofLflashhitList){
307  if(tf_lflashhit->tof >= ftof_Lflshhit_cut) is_cos_like.push_back(false);
308  else is_cos_like.push_back(true);
309  }
310 
311  if(is_cos_like.size()){
312  bool found_nu=false;
313  for(auto itr: is_cos_like){
314  if(!itr){
315  found_nu=true;
316  break;
317  }
318  }
319 
320  if(!found_nu) keep_event=false;
321  }
322  } // Using earliest hit of the largest optical flash to calculate tof
323 
324  //======================================================================================================================
325 
326  //=========================Calculation ToF values using Earliest hit of the Closest flash ==============================
327 
328  if(fuse_Cflsh_hit){
329  art::Handle< std::vector<sbnd::ToF::ToF> > tofCflashhitListHandle;
330  std::vector< art::Ptr<sbnd::ToF::ToF> > tofCflashhitList;
331  if( evt.getByLabel(ftofCflashhitLabel,tofCflashhitListHandle))
332  art::fill_ptr_vector(tofCflashhitList, tofCflashhitListHandle);
333 
334  std::vector<bool> is_cos_like;
335 
336  for(auto const& tf_cflashhit : tofCflashhitList){
337  if(tf_cflashhit->tof >= ftof_Cflshhit_cut) is_cos_like.push_back(false);
338  else is_cos_like.push_back(true);
339  }
340 
341  if(is_cos_like.size()){
342  bool found_nu=false;
343  for(auto itr: is_cos_like){
344  if(!itr){
345  found_nu=true;
346  break;
347  }
348  }
349 
350  if(!found_nu) keep_event=false;
351  }
352  } // Using earliest hit of the closest optical flash to calculate tof
353 
354  //======================================================================================================================
355 
356  return keep_event;
357 }
art::InputTag ftofCflashLabel
art::InputTag ftofChitLabel
art::InputTag ftofCflashhitLabel
art::InputTag ftofLflashhitLabel
TCEvent evt
Definition: DataStructs.cxx:8
art::InputTag ftofLflashLabel
art::InputTag ftofLhitLabel
ToFFilter& sbnd::ToFFilter::operator= ( ToFFilter const &  )
delete
ToFFilter& sbnd::ToFFilter::operator= ( ToFFilter &&  )
delete

Member Data Documentation

float sbnd::ToFFilter::ftof_Cflsh_cut
private

Definition at line 143 of file ToFFilter_module.cc.

float sbnd::ToFFilter::ftof_Cflshhit_cut
private

Definition at line 145 of file ToFFilter_module.cc.

float sbnd::ToFFilter::ftof_Chit_cut
private

Definition at line 141 of file ToFFilter_module.cc.

float sbnd::ToFFilter::ftof_Lflsh_cut
private

Definition at line 142 of file ToFFilter_module.cc.

float sbnd::ToFFilter::ftof_Lflshhit_cut
private

Definition at line 144 of file ToFFilter_module.cc.

float sbnd::ToFFilter::ftof_Lhit_cut
private

Definition at line 140 of file ToFFilter_module.cc.

art::InputTag sbnd::ToFFilter::ftofCflashhitLabel
private

Definition at line 133 of file ToFFilter_module.cc.

art::InputTag sbnd::ToFFilter::ftofCflashLabel
private

Definition at line 131 of file ToFFilter_module.cc.

art::InputTag sbnd::ToFFilter::ftofChitLabel
private

Definition at line 129 of file ToFFilter_module.cc.

art::InputTag sbnd::ToFFilter::ftofLflashhitLabel
private

Definition at line 132 of file ToFFilter_module.cc.

art::InputTag sbnd::ToFFilter::ftofLflashLabel
private

Definition at line 130 of file ToFFilter_module.cc.

art::InputTag sbnd::ToFFilter::ftofLhitLabel
private

Definition at line 128 of file ToFFilter_module.cc.

bool sbnd::ToFFilter::fuse_Cflsh
private

Definition at line 137 of file ToFFilter_module.cc.

bool sbnd::ToFFilter::fuse_Cflsh_hit
private

Definition at line 139 of file ToFFilter_module.cc.

bool sbnd::ToFFilter::fuse_Chit
private

Definition at line 135 of file ToFFilter_module.cc.

bool sbnd::ToFFilter::fuse_Lflsh
private

Definition at line 136 of file ToFFilter_module.cc.

bool sbnd::ToFFilter::fuse_Lflsh_hit
private

Definition at line 138 of file ToFFilter_module.cc.

bool sbnd::ToFFilter::fuse_Lhit
private

Definition at line 134 of file ToFFilter_module.cc.


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