Changeset f3c4047 in git for modules/PdgCodeFilter.cc
- Timestamp:
- Jun 26, 2015, 3:22:13 PM (9 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- 9458a020
- Parents:
- 28027d5 (diff), fa42514 (diff)
Note: this is a merge changeset, the changes displayed below correspond to the merge itself.
Use the(diff)
links above to see all the changes relative to each parent. - File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
modules/PdgCodeFilter.cc
r28027d5 rf3c4047 74 74 fPTMin = GetDouble("PTMin", 0.0); 75 75 76 fInvertPdg = GetBool("InvertPdg", false); 77 78 fRequireStatus = GetBool("RequireStatus", false); 79 fStatus = GetInt("Status", 1); 80 76 81 // import input array 77 82 fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles")); … … 109 114 Double_t pt; 110 115 116 const Bool_t requireStatus = fRequireStatus; 117 const Bool_t invertPdg = fInvertPdg; 118 const int reqStatus = fStatus; 119 111 120 fItInputArray->Reset(); 112 121 while((candidate = static_cast<Candidate*>(fItInputArray->Next()))) … … 116 125 pt = candidateMomentum.Pt(); 117 126 127 if(pt < fPTMin) continue; 128 if(requireStatus && (candidate->Status != reqStatus)) continue; 129 118 130 pass = kTRUE; 119 120 if(pt < fPTMin) pass = kFALSE;121 131 if(find(fPdgCodes.begin(), fPdgCodes.end(), pdgCode) != fPdgCodes.end()) pass = kFALSE; 122 132 133 if (invertPdg) pass = !pass; 123 134 if(pass) fOutputArray->Add(candidate); 124 135 }
Note:
See TracChangeset
for help on using the changeset viewer.