![]() |
![]() |
#include <EventWithHistoryFilter.h>
Public Member Functions | |
EventWithHistoryFilter () | |
EventWithHistoryFilter (const edm::ParameterSet &iConfig) | |
const bool | selected (const EventWithHistory &he, const edm::Event &iEvent, const edm::EventSetup &iSetup) const |
const bool | selected (const EventWithHistory &he, const edm::EventSetup &iSetup) const |
const bool | selected (const edm::Event &event, const edm::EventSetup &iSetup) const |
void | set (const edm::ParameterSet &iConfig) |
Private Member Functions | |
const int | getAPVLatency (const edm::EventSetup &iSetup) const |
const int | getAPVMode (const edm::EventSetup &iSetup) const |
const std::vector< int > | getAPVPhase (const edm::Event &iEvent) const |
const bool | is_selected (const EventWithHistory &he, const edm::EventSetup &iSetup, const std::vector< int > apvphases) const |
const bool | isAPVLatencyNotNeeded () const |
const bool | isAPVModeNotNeeded () const |
const bool | isAPVPhaseNotNeeded () const |
const bool | isCutInactive (const std::vector< int > &range) const |
const bool | isInRange (const long long bx, const std::vector< int > &range, const bool extra) const |
void | printConfig () const |
Private Attributes | |
std::vector< int > | _apvmodes |
std::string | _APVPhaseLabel |
std::vector< int > | _bxcyclerange |
std::vector< int > | _bxcyclerangelat |
std::vector< int > | _bxrange |
std::vector< int > | _bxrangelat |
std::vector< int > | _dbxcyclerange |
std::vector< int > | _dbxcyclerangelat |
std::vector< int > | _dbxrange |
std::vector< int > | _dbxrangelat |
std::vector< int > | _dbxtrpltrange |
edm::InputTag | _historyProduct |
bool | _noAPVPhase |
std::string | _partition |
Definition at line 14 of file EventWithHistoryFilter.h.
EventWithHistoryFilter::EventWithHistoryFilter | ( | ) |
Definition at line 17 of file EventWithHistoryFilter.cc.
References printConfig().
: _historyProduct(), _partition(), _APVPhaseLabel(), _apvmodes(), _dbxrange(), _dbxrangelat(), _bxrange(), _bxrangelat(), _bxcyclerange(), _bxcyclerangelat(), _dbxcyclerange(), _dbxcyclerangelat(), _dbxtrpltrange(), _noAPVPhase(true) { printConfig(); }
EventWithHistoryFilter::EventWithHistoryFilter | ( | const edm::ParameterSet & | iConfig | ) |
Definition at line 32 of file EventWithHistoryFilter.cc.
References _noAPVPhase, isAPVPhaseNotNeeded(), and printConfig().
: _historyProduct(iConfig.getUntrackedParameter<edm::InputTag>("historyProduct",edm::InputTag("consecutiveHEs"))), _partition(iConfig.getUntrackedParameter<std::string>("partitionName","Any")), _APVPhaseLabel(iConfig.getUntrackedParameter<std::string>("APVPhaseLabel","APVPhases")), _apvmodes(iConfig.getUntrackedParameter<std::vector<int> >("apvModes",std::vector<int>())), _dbxrange(iConfig.getUntrackedParameter<std::vector<int> >("dbxRange",std::vector<int>())), _dbxrangelat(iConfig.getUntrackedParameter<std::vector<int> >("dbxRangeLtcyAware",std::vector<int>())), _bxrange(iConfig.getUntrackedParameter<std::vector<int> >("absBXRange",std::vector<int>())), _bxrangelat(iConfig.getUntrackedParameter<std::vector<int> >("absBXRangeLtcyAware",std::vector<int>())), _bxcyclerange(iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRange",std::vector<int>())), _bxcyclerangelat(iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRangeLtcyAware",std::vector<int>())), _dbxcyclerange(iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRange",std::vector<int>())), _dbxcyclerangelat(iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRangeLtcyAware",std::vector<int>())), _dbxtrpltrange(iConfig.getUntrackedParameter<std::vector<int> >("dbxTripletRange",std::vector<int>())) { _noAPVPhase = isAPVPhaseNotNeeded(); printConfig(); }
const int EventWithHistoryFilter::getAPVLatency | ( | const edm::EventSetup & | iSetup | ) | const [private] |
Definition at line 165 of file EventWithHistoryFilter.cc.
References edm::EventSetup::get(), and isAPVLatencyNotNeeded().
Referenced by is_selected().
{ if(isAPVLatencyNotNeeded()) return -1; edm::ESHandle<SiStripLatency> apvlat; iSetup.get<SiStripLatencyRcd>().get(apvlat); const int latency = apvlat->singleLatency()!=255 ? apvlat->singleLatency(): -1; // thrown an exception if latency value is invalid /* if(latency < 0 && !isAPVLatencyNotNeeded()) throw cms::Exception("InvalidAPVLatency") << " invalid APVLatency found "; */ return latency; }
const int EventWithHistoryFilter::getAPVMode | ( | const edm::EventSetup & | iSetup | ) | const [private] |
Definition at line 183 of file EventWithHistoryFilter.cc.
References edm::EventSetup::get(), isAPVModeNotNeeded(), and mode.
Referenced by is_selected().
{ if(isAPVModeNotNeeded()) return -1; edm::ESHandle<SiStripLatency> apvlat; iSetup.get<SiStripLatencyRcd>().get(apvlat); int mode = -1; if(apvlat->singleReadOutMode()==1) mode = 47; if(apvlat->singleReadOutMode()==0) mode = 37; // thrown an exception if mode value is invalid /* if(mode < 0 && !isAPVModeNotNeeded()) throw cms::Exception("InvalidAPVMode") << " invalid APVMode found "; */ return mode; }
const std::vector< int > EventWithHistoryFilter::getAPVPhase | ( | const edm::Event & | iEvent | ) | const [private] |
Definition at line 203 of file EventWithHistoryFilter.cc.
References _APVPhaseLabel, _noAPVPhase, _partition, and edm::Event::getByLabel().
Referenced by selected().
{ if(_noAPVPhase) { const std::vector<int> dummy; return dummy; } edm::Handle<APVCyclePhaseCollection> apvPhases; iEvent.getByLabel(_APVPhaseLabel,apvPhases); const std::vector<int> phases = apvPhases->getPhases(_partition.c_str()); /* if(!_noAPVPhase) { if(phases.size()==0) throw cms::Exception("NoPartitionAPVPhase") << " No APV phase for partition " << _partition.c_str() << " : check if a proper partition has been chosen "; } */ return phases; }
const bool EventWithHistoryFilter::is_selected | ( | const EventWithHistory & | he, |
const edm::EventSetup & | iSetup, | ||
const std::vector< int > | apvphases | ||
) | const [private] |
Definition at line 98 of file EventWithHistoryFilter.cc.
References _apvmodes, _bxcyclerange, _bxcyclerangelat, _bxrange, _bxrangelat, _dbxcyclerange, _dbxcyclerangelat, _dbxrange, _dbxrangelat, _dbxtrpltrange, EventWithHistory::absoluteBX(), EventWithHistory::absoluteBXinCycle(), EventWithHistory::deltaBX(), EventWithHistory::deltaBXinCycle(), EventWithHistory::depth(), getAPVLatency(), getAPVMode(), isAPVModeNotNeeded(), isCutInactive(), isInRange(), and selected().
Referenced by selected().
{ const int latency = getAPVLatency(iSetup); bool selected = true; if(!isAPVModeNotNeeded()) { const int apvmode = getAPVMode(iSetup); bool modeok = false; for(std::vector<int>::const_iterator wantedmode =_apvmodes.begin();wantedmode!=_apvmodes.end();++wantedmode) { modeok = modeok || (apvmode == *wantedmode); } if(!modeok) return false; } selected = selected && (isCutInactive(_dbxrange) || isInRange(he.deltaBX(),_dbxrange,he.depth()!=0)); selected = selected && (isCutInactive(_dbxrangelat) || isInRange(he.deltaBX()-latency,_dbxrangelat,he.depth()!=0 && latency>=0)); selected = selected && (isCutInactive(_bxrange) || isInRange(he.absoluteBX()%70,_bxrange,1)); selected = selected && (isCutInactive(_bxrangelat) || isInRange((he.absoluteBX()-latency)%70,_bxrangelat,latency>=0)); // loop on all the phases and require that the cut is fulfilled for at least one of them bool phaseselected; phaseselected = isCutInactive(_bxcyclerange); for(std::vector<int>::const_iterator phase=apvphases.begin();phase!=apvphases.end();++phase) { phaseselected = phaseselected || isInRange(he.absoluteBXinCycle(*phase)%70,_bxcyclerange,*phase>0); } selected = selected && phaseselected; phaseselected = isCutInactive(_bxcyclerangelat); for(std::vector<int>::const_iterator phase=apvphases.begin();phase!=apvphases.end();++phase) { phaseselected = phaseselected || isInRange((he.absoluteBXinCycle(*phase)-latency)%70,_bxcyclerangelat, *phase>=0 && latency>=0); } selected = selected && phaseselected; phaseselected = isCutInactive(_dbxcyclerange); for(std::vector<int>::const_iterator phase=apvphases.begin();phase!=apvphases.end();++phase) { phaseselected = phaseselected || isInRange(he.deltaBXinCycle(*phase),_dbxcyclerange,he.depth()!=0 && *phase>=0); } selected = selected && phaseselected; phaseselected = isCutInactive(_dbxcyclerangelat); for(std::vector<int>::const_iterator phase=apvphases.begin();phase!=apvphases.end();++phase) { phaseselected = phaseselected || isInRange(he.deltaBXinCycle(*phase)-latency,_dbxcyclerangelat, he.depth()!=0 && *phase>=0 && latency>=0); } selected = selected && phaseselected; // end of phase-dependent cuts selected = selected && (isCutInactive(_dbxtrpltrange) || isInRange(he.deltaBX(1,2),_dbxtrpltrange,he.depth()>1)); return selected; }
const bool EventWithHistoryFilter::isAPVLatencyNotNeeded | ( | ) | const [private] |
Definition at line 225 of file EventWithHistoryFilter.cc.
References _bxcyclerangelat, _bxrangelat, _dbxcyclerangelat, _dbxrangelat, and isCutInactive().
Referenced by getAPVLatency().
{ return isCutInactive(_bxrangelat) && isCutInactive(_dbxrangelat) && isCutInactive(_bxcyclerangelat) && isCutInactive(_dbxcyclerangelat); }
const bool EventWithHistoryFilter::isAPVModeNotNeeded | ( | ) | const [private] |
Definition at line 245 of file EventWithHistoryFilter.cc.
References _apvmodes.
Referenced by getAPVMode(), and is_selected().
{ return (_apvmodes.size()==0) ; }
const bool EventWithHistoryFilter::isAPVPhaseNotNeeded | ( | ) | const [private] |
Definition at line 235 of file EventWithHistoryFilter.cc.
References _bxcyclerange, _bxcyclerangelat, _dbxcyclerange, _dbxcyclerangelat, and isCutInactive().
Referenced by EventWithHistoryFilter(), and set().
{ return isCutInactive(_bxcyclerange) && isCutInactive(_dbxcyclerange) && isCutInactive(_bxcyclerangelat) && isCutInactive(_dbxcyclerangelat); }
const bool EventWithHistoryFilter::isCutInactive | ( | const std::vector< int > & | range | ) | const [private] |
Definition at line 251 of file EventWithHistoryFilter.cc.
Referenced by is_selected(), isAPVLatencyNotNeeded(), isAPVPhaseNotNeeded(), and printConfig().
{
return ((range.size()==0 ||
(range.size()==1 && range[0]<0) ||
(range.size()==2 && range[0]<0 && range[1]<0)));
}
const bool EventWithHistoryFilter::isInRange | ( | const long long | bx, |
const std::vector< int > & | range, | ||
const bool | extra | ||
) | const [private] |
Definition at line 259 of file EventWithHistoryFilter.cc.
Referenced by is_selected().
{ bool cut1 = range.size()<1 || range[0]<0 || (extra && bx >= range[0]); bool cut2 = range.size()<2 || range[1]<0 || (extra && bx <= range[1]); if(range.size()>=2 && range[0]>=0 && range[1]>=0 && (range[0] > range[1])) { return cut1 || cut2; } else { return cut1 && cut2; } }
void EventWithHistoryFilter::printConfig | ( | ) | const [private] |
Definition at line 273 of file EventWithHistoryFilter.cc.
References _APVPhaseLabel, _bxcyclerange, _bxcyclerangelat, _bxrange, _bxrangelat, _dbxcyclerange, _dbxcyclerangelat, _dbxrange, _dbxrangelat, _dbxtrpltrange, _historyProduct, _partition, edm::InputTag::instance(), isCutInactive(), edm::InputTag::label(), and edm::InputTag::process().
Referenced by EventWithHistoryFilter(), and set().
{ std::string msgcategory = "EventWithHistoryFilterConfiguration"; if(!( isCutInactive(_bxrange) && isCutInactive(_bxrangelat) && isCutInactive(_bxcyclerange) && isCutInactive(_bxcyclerangelat) && isCutInactive(_dbxrange) && isCutInactive(_dbxrangelat) && isCutInactive(_dbxcyclerange) && isCutInactive(_dbxcyclerangelat) && isCutInactive(_dbxtrpltrange) )) { edm::LogInfo(msgcategory.c_str()) << "historyProduct: " << _historyProduct.label() << " " << _historyProduct.instance() << " " << _historyProduct.process() << " " << " APVCyclePhase: " << _APVPhaseLabel; edm::LogVerbatim(msgcategory.c_str()) << "-----------------------"; edm::LogVerbatim(msgcategory.c_str()) << "List of active cuts:"; if(!isCutInactive(_bxrange)) { edm::LogVerbatim(msgcategory.c_str()) << "......................"; if(_bxrange.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBX lower limit " << _bxrange[0]; if(_bxrange.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBX upper limit " << _bxrange[1]; edm::LogVerbatim(msgcategory.c_str()) << "......................"; } if(!isCutInactive(_bxrangelat)) { edm::LogVerbatim(msgcategory.c_str()) << "......................"; if(_bxrangelat.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXLtcyAware lower limit " << _bxrangelat[0]; if(_bxrangelat.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXLtcyAware upper limit " << _bxrangelat[1]; edm::LogVerbatim(msgcategory.c_str()) << "......................"; } if(!isCutInactive(_bxcyclerange)) { edm::LogVerbatim(msgcategory.c_str()) << "......................"; edm::LogVerbatim(msgcategory.c_str()) <<"absoluteBXinCycle partition: " << _partition; if(_bxcyclerange.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXinCycle lower limit " << _bxcyclerange[0]; if(_bxcyclerange.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXinCycle upper limit " << _bxcyclerange[1]; edm::LogVerbatim(msgcategory.c_str()) << "......................"; } if(!isCutInactive(_bxcyclerangelat)) { edm::LogVerbatim(msgcategory.c_str()) << "......................"; edm::LogVerbatim(msgcategory.c_str()) <<"absoluteBXinCycleLtcyAware partition: " << _partition; if(_bxcyclerangelat.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXinCycleLtcyAware lower limit " << _bxcyclerangelat[0]; if(_bxcyclerangelat.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXinCycleLtcyAware upper limit " << _bxcyclerangelat[1]; edm::LogVerbatim(msgcategory.c_str()) << "......................"; } if(!isCutInactive(_dbxrange)) { edm::LogVerbatim(msgcategory.c_str()) << "......................"; if(_dbxrange.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "deltaBX lower limit " << _dbxrange[0]; if(_dbxrange.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "deltaBX upper limit " << _dbxrange[1]; edm::LogVerbatim(msgcategory.c_str()) << "......................"; } if(!isCutInactive(_dbxrangelat)) { edm::LogVerbatim(msgcategory.c_str()) << "......................"; if(_dbxrangelat.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXLtcyAware lower limit " << _dbxrangelat[0]; if(_dbxrangelat.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXLtcyAware upper limit " << _dbxrangelat[1]; edm::LogVerbatim(msgcategory.c_str()) << "......................"; } if(!isCutInactive(_dbxcyclerange)) { edm::LogVerbatim(msgcategory.c_str()) << "......................"; edm::LogVerbatim(msgcategory.c_str()) <<"deltaBXinCycle partition: " << _partition; if(_dbxcyclerange.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXinCycle lower limit " << _dbxcyclerange[0]; if(_dbxcyclerange.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXinCycle upper limit " << _dbxcyclerange[1]; edm::LogVerbatim(msgcategory.c_str()) << "......................"; } if(!isCutInactive(_dbxcyclerangelat)) { edm::LogVerbatim(msgcategory.c_str()) << "......................"; edm::LogVerbatim(msgcategory.c_str()) <<"deltaBXinCycleLtcyAware partition: " << _partition; if(_dbxcyclerangelat.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXinCycleLtcyAware lower limit " << _dbxcyclerangelat[0]; if(_dbxcyclerangelat.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXinCycleLtcyAware upper limit " << _dbxcyclerangelat[1]; edm::LogVerbatim(msgcategory.c_str()) << "......................"; } if(!isCutInactive(_dbxtrpltrange)) { edm::LogVerbatim(msgcategory.c_str()) << "......................"; if(_dbxtrpltrange.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "TripletIsolation lower limit " << _dbxtrpltrange[0]; if(_dbxtrpltrange.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "TripletIsolation upper limit " << _dbxtrpltrange[1]; edm::LogVerbatim(msgcategory.c_str()) << "......................"; } edm::LogVerbatim(msgcategory.c_str()) << "-----------------------"; } }
const bool EventWithHistoryFilter::selected | ( | const edm::Event & | event, |
const edm::EventSetup & | iSetup | ||
) | const |
Definition at line 86 of file EventWithHistoryFilter.cc.
References _historyProduct, getAPVPhase(), and is_selected().
{ const std::vector<int> apvphases = getAPVPhase(event); edm::Handle<EventWithHistory> hEvent; event.getByLabel(_historyProduct,hEvent); return is_selected(*hEvent,iSetup,apvphases); }
const bool EventWithHistoryFilter::selected | ( | const EventWithHistory & | he, |
const edm::EventSetup & | iSetup | ||
) | const |
Definition at line 72 of file EventWithHistoryFilter.cc.
References is_selected().
Referenced by is_selected().
{ const std::vector<int> dummy; return is_selected(he,iSetup,dummy); }
const bool EventWithHistoryFilter::selected | ( | const EventWithHistory & | he, |
const edm::Event & | iEvent, | ||
const edm::EventSetup & | iSetup | ||
) | const |
Definition at line 79 of file EventWithHistoryFilter.cc.
References getAPVPhase(), and is_selected().
{ const std::vector<int> apvphases = getAPVPhase(iEvent); return is_selected(he,iSetup,apvphases); }
void EventWithHistoryFilter::set | ( | const edm::ParameterSet & | iConfig | ) |
Definition at line 51 of file EventWithHistoryFilter.cc.
References _APVPhaseLabel, _bxcyclerange, _bxcyclerangelat, _bxrange, _bxrangelat, _dbxcyclerange, _dbxcyclerangelat, _dbxrange, _dbxrangelat, _dbxtrpltrange, _historyProduct, _noAPVPhase, _partition, edm::ParameterSet::getUntrackedParameter(), isAPVPhaseNotNeeded(), and printConfig().
{ _historyProduct = iConfig.getUntrackedParameter<edm::InputTag>("historyProduct",edm::InputTag("consecutiveHEs")); _partition = iConfig.getUntrackedParameter<std::string>("partitionName","Any"); _APVPhaseLabel = iConfig.getUntrackedParameter<std::string>("APVPhaseLabel","APVPhases"); _dbxrange = iConfig.getUntrackedParameter<std::vector<int> >("dbxRange",std::vector<int>()); _dbxrangelat = iConfig.getUntrackedParameter<std::vector<int> >("dbxRangeLtcyAware",std::vector<int>()); _bxrange = iConfig.getUntrackedParameter<std::vector<int> >("absBXRange",std::vector<int>()); _bxrangelat = iConfig.getUntrackedParameter<std::vector<int> >("absBXRangeLtcyAware",std::vector<int>()); _bxcyclerange = iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRange",std::vector<int>()); _bxcyclerangelat = iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRangeLtcyAware",std::vector<int>()); _dbxcyclerange = iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRange",std::vector<int>()); _dbxcyclerangelat = iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRangeLtcyAware",std::vector<int>()); _dbxtrpltrange = iConfig.getUntrackedParameter<std::vector<int> >("dbxTripletRange",std::vector<int>()); _noAPVPhase = isAPVPhaseNotNeeded(); printConfig(); }
std::vector<int> EventWithHistoryFilter::_apvmodes [private] |
Definition at line 41 of file EventWithHistoryFilter.h.
Referenced by is_selected(), and isAPVModeNotNeeded().
std::string EventWithHistoryFilter::_APVPhaseLabel [private] |
Definition at line 40 of file EventWithHistoryFilter.h.
Referenced by getAPVPhase(), printConfig(), and set().
std::vector<int> EventWithHistoryFilter::_bxcyclerange [private] |
Definition at line 46 of file EventWithHistoryFilter.h.
Referenced by is_selected(), isAPVPhaseNotNeeded(), printConfig(), and set().
std::vector<int> EventWithHistoryFilter::_bxcyclerangelat [private] |
Definition at line 47 of file EventWithHistoryFilter.h.
Referenced by is_selected(), isAPVLatencyNotNeeded(), isAPVPhaseNotNeeded(), printConfig(), and set().
std::vector<int> EventWithHistoryFilter::_bxrange [private] |
Definition at line 44 of file EventWithHistoryFilter.h.
Referenced by is_selected(), printConfig(), and set().
std::vector<int> EventWithHistoryFilter::_bxrangelat [private] |
Definition at line 45 of file EventWithHistoryFilter.h.
Referenced by is_selected(), isAPVLatencyNotNeeded(), printConfig(), and set().
std::vector<int> EventWithHistoryFilter::_dbxcyclerange [private] |
Definition at line 48 of file EventWithHistoryFilter.h.
Referenced by is_selected(), isAPVPhaseNotNeeded(), printConfig(), and set().
std::vector<int> EventWithHistoryFilter::_dbxcyclerangelat [private] |
Definition at line 49 of file EventWithHistoryFilter.h.
Referenced by is_selected(), isAPVLatencyNotNeeded(), isAPVPhaseNotNeeded(), printConfig(), and set().
std::vector<int> EventWithHistoryFilter::_dbxrange [private] |
Definition at line 42 of file EventWithHistoryFilter.h.
Referenced by is_selected(), printConfig(), and set().
std::vector<int> EventWithHistoryFilter::_dbxrangelat [private] |
Definition at line 43 of file EventWithHistoryFilter.h.
Referenced by is_selected(), isAPVLatencyNotNeeded(), printConfig(), and set().
std::vector<int> EventWithHistoryFilter::_dbxtrpltrange [private] |
Definition at line 50 of file EventWithHistoryFilter.h.
Referenced by is_selected(), printConfig(), and set().
Definition at line 38 of file EventWithHistoryFilter.h.
Referenced by printConfig(), selected(), and set().
bool EventWithHistoryFilter::_noAPVPhase [private] |
Definition at line 51 of file EventWithHistoryFilter.h.
Referenced by EventWithHistoryFilter(), getAPVPhase(), and set().
std::string EventWithHistoryFilter::_partition [private] |
Definition at line 39 of file EventWithHistoryFilter.h.
Referenced by getAPVPhase(), printConfig(), and set().