CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/HLTrigger/HLTcore/plugins/TriggerSummaryProducerRAW.cc

Go to the documentation of this file.
00001 
00012 #include "HLTrigger/HLTcore/interface/TriggerSummaryProducerRAW.h"
00013 
00014 #include "DataFormats/Common/interface/Handle.h"
00015 #include "DataFormats/Common/interface/OrphanHandle.h"
00016 #include "DataFormats/HLTReco/interface/TriggerEventWithRefs.h"
00017 #include "DataFormats/Provenance/interface/Provenance.h"
00018 #include "FWCore/Framework/interface/ProcessMatch.h"
00019 #include "FWCore/Framework/interface/TriggerNamesService.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 #include "FWCore/ServiceRegistry/interface/Service.h"
00022 #include "FWCore/Utilities/interface/InputTag.h"
00023 
00024 #include <memory>
00025 #include<vector>
00026 
00027 //
00028 // constructors and destructor
00029 //
00030 TriggerSummaryProducerRAW::TriggerSummaryProducerRAW(const edm::ParameterSet& ps) : 
00031   pn_(ps.getParameter<std::string>("processName"))
00032 {
00033   if (pn_=="@") {
00034     
00035     edm::Service<edm::service::TriggerNamesService> tns;
00036     if (tns.isAvailable()) {
00037       pn_ = tns->getProcessName();
00038     } else {
00039       edm::LogError("TriggerSummaryProducerRaw") << "HLT Error: TriggerNamesService not available!";
00040       pn_="*";
00041     }
00042   }
00043 
00044   LogDebug("TriggerSummaryProducerRaw") << "Using process name: '" << pn_ <<"'";
00045   produces<trigger::TriggerEventWithRefs>();
00046 
00047   // Tell the getter what type of products to get and
00048   // also the process to get them from
00049   getterOfProducts_ = edm::GetterOfProducts<trigger::TriggerFilterObjectWithRefs>(edm::ProcessMatch(pn_), this);
00050   callWhenNewProductsRegistered(getterOfProducts_);
00051 }
00052 
00053 TriggerSummaryProducerRAW::~TriggerSummaryProducerRAW()
00054 {
00055 }
00056 
00057 //
00058 // member functions
00059 //
00060 
00061 // ------------ method called to produce the data  ------------
00062 void
00063 TriggerSummaryProducerRAW::produce(edm::Event& iEvent, const edm::EventSetup&)
00064 {
00065    using namespace std;
00066    using namespace edm;
00067    using namespace reco;
00068    using namespace trigger;
00069 
00070    std::vector<edm::Handle<trigger::TriggerFilterObjectWithRefs> > fobs;
00071    getterOfProducts_.fillHandles(iEvent, fobs);
00072 
00073    const unsigned int nfob(fobs.size());
00074    LogDebug("TriggerSummaryProducerRaw") << "Number of filter objects found: " << nfob;
00075 
00076    // construct single RAW product
00077    auto_ptr<TriggerEventWithRefs> product(new TriggerEventWithRefs(pn_,nfob));
00078    for (unsigned int ifob=0; ifob!=nfob; ++ifob) {
00079      const string& label    (fobs[ifob].provenance()->moduleLabel());
00080      const string& instance (fobs[ifob].provenance()->productInstanceName());
00081      const string& process  (fobs[ifob].provenance()->processName());
00082      const InputTag tag(label,instance,process);
00083      LogTrace("TriggerSummaryProducerRaw")
00084        << ifob << " " << tag << endl
00085        << " Sizes: "
00086        << " 1/" << fobs[ifob]->photonSize()
00087        << " 2/" << fobs[ifob]->electronSize()
00088        << " 3/" << fobs[ifob]->muonSize()
00089        << " 4/" << fobs[ifob]->jetSize()
00090        << " 5/" << fobs[ifob]->compositeSize()
00091        << " 6/" << fobs[ifob]->basemetSize()
00092        << " 7/" << fobs[ifob]->calometSize()
00093 
00094        << " 8/" << fobs[ifob]->pixtrackSize()
00095        << " 9/" << fobs[ifob]->l1emSize()
00096        << " A/" << fobs[ifob]->l1muonSize()
00097        << " B/" << fobs[ifob]->l1jetSize()
00098        << " C/" << fobs[ifob]->l1etmissSize()
00099        << " D/" << fobs[ifob]->l1hfringsSize()
00100        << " E/" << fobs[ifob]->pfjetSize()
00101        << " F/" << fobs[ifob]->pftauSize()
00102        << endl;
00103      product->addFilterObject(tag,*fobs[ifob]);
00104    }
00105 
00106    // place product in Event
00107    OrphanHandle<TriggerEventWithRefs> ref = iEvent.put(product);
00108    LogTrace("TriggerSummaryProducerRaw") << "Number of filter objects packed: " << ref->size();
00109 
00110    return;
00111 }