Changeset 7a86a5c in git
- Timestamp:
- Nov 23, 2018, 2:37:05 AM (6 years ago)
- Branches:
- ImprovedOutputFile, Timing, llp, master
- Children:
- 7c7fe5e
- Parents:
- e5fa629
- Location:
- modules
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
modules/PdgCodeFilter.cc
re5fa629 r7a86a5c 76 76 fInvert = GetBool("Invert", false); 77 77 78 // no pileup 79 fRequireNotPileup = GetBool("RequireNotPileup", false); 80 78 81 fRequireStatus = GetBool("RequireStatus", false); 79 82 fStatus = GetInt("Status", 1); … … 127 130 if(fRequireStatus && (candidate->Status != fStatus)) continue; 128 131 if(fRequireCharge && (candidate->Charge != fCharge)) continue; 132 if(fRequireNotPileup && (candidate->IsPU >0 )) continue; 129 133 130 134 pass = kTRUE; -
modules/PdgCodeFilter.h
re5fa629 r7a86a5c 55 55 Bool_t fRequireCharge; //! 56 56 Int_t fCharge; //! 57 57 Bool_t fRequireNotPileup; //! 58 58 59 59 std::vector<Int_t> fPdgCodes; -
modules/StatusPidFilter.cc
re5fa629 r7a86a5c 153 153 { 154 154 // PT threshold 155 156 155 fPTMin = GetDouble("PTMin", 0.5); 156 157 // keep or remove pileup particles 158 fRequireNotPileup = GetBool("RequireNotPileup", false); 157 159 158 160 // import input array … … 227 229 if(!pass || (candidate->Momentum.Pt() < fPTMin && !(is_b_hadron || is_b_quark || is_tau_daughter || is_W_daughter)) ) continue; 228 230 231 // not pileup particles 232 if(fRequireNotPileup && (candidate->IsPU >0)) continue; 233 229 234 fOutputArray->Add(candidate); 230 235 } -
modules/StatusPidFilter.h
re5fa629 r7a86a5c 51 51 Double_t fPTMin; //! 52 52 53 Bool_t fRequireNotPileup; //! 54 53 55 TIterator *fItInputArray; //! 54 56
Note:
See TracChangeset
for help on using the changeset viewer.