CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/DQMOffline/Trigger/src/GeneralHLTOffline.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    GeneralHLTOffline
00004 // Class:      GeneralHLTOffline
00005 // 
00013 //
00014 // Original Author:  Jason Michael Slaunwhite,512 1-008,`+41227670494,
00015 //         Created:  Fri Aug  5 10:34:47 CEST 2011
00016 // $Id: GeneralHLTOffline.cc,v 1.3 2012/05/07 15:23:27 jcarson Exp $
00017 //
00018 //
00019 
00020 // system include files
00021 #include <memory>
00022 
00023 // user include files
00024 #include "DQMServices/Core/interface/DQMStore.h"
00025 #include "DQMServices/Core/interface/MonitorElement.h"
00026 
00027 #include "DataFormats/Common/interface/TriggerResults.h"
00028 #include "DataFormats/Common/interface/TriggerResults.h"
00029 #include "DataFormats/HLTReco/interface/TriggerObject.h"
00030 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
00031 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
00032 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
00033 
00034 #include "FWCore/Framework/interface/Frameworkfwd.h"
00035 #include "FWCore/Framework/interface/EDAnalyzer.h"
00036 #include "FWCore/Framework/interface/Event.h"
00037 #include "FWCore/Framework/interface/MakerMacros.h"
00038 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00039 #include "FWCore/ServiceRegistry/interface/Service.h"
00040 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00041 
00042 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
00043 
00044 #include "TMath.h"
00045 
00046 #include "TStyle.h"
00047 //
00048 // class declaration
00049 //
00050 
00051 using namespace edm;
00052 using namespace trigger;
00053 using std::vector;
00054 using std::string;
00055 
00056 
00057 class GeneralHLTOffline : public edm::EDAnalyzer {
00058    public:
00059       explicit GeneralHLTOffline(const edm::ParameterSet&);
00060       ~GeneralHLTOffline();
00061 
00062   //static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
00063 
00064 
00065    private:
00066       virtual void beginJob() ;
00067       virtual void analyze(const edm::Event&, const edm::EventSetup&);
00068       virtual void endJob() ;
00069 
00070       virtual void beginRun(edm::Run const&, edm::EventSetup const&);
00071       virtual void endRun(edm::Run const&, edm::EventSetup const&);
00072       virtual void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&);
00073       virtual void endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&);
00074       virtual void setupHltMatrix(std::string, int);
00075   virtual void fillHltMatrix(std::string, std::string, double, double, bool);
00076 
00077       // ----------member data ---------------------------
00078 
00079 
00080   bool debugPrint;
00081   bool outputPrint;
00082 
00083   std::string plotDirectoryName;
00084 
00085   DQMStore * dbe;
00086 
00087   HLTConfigProvider hltConfig_;
00088 
00089   vector< vector<string> > PDsVectorPathsVector;
00090   
00091 };
00092 
00093 //
00094 // constants, enums and typedefs
00095 //
00096 
00097 //
00098 // static data member definitions
00099 //
00100 
00101 //
00102 // constructors and destructor
00103 //
00104 GeneralHLTOffline::GeneralHLTOffline(const edm::ParameterSet& iConfig)
00105 
00106 {
00107    //now do what ever initialization is needed
00108 
00109   debugPrint = false;
00110   outputPrint = false;
00111 
00112   if (debugPrint) std::cout << "Inside Constructor" << std::endl;
00113 
00114   plotDirectoryName = iConfig.getUntrackedParameter<std::string>("dirname", "HLT/General");
00115 
00116   if (debugPrint) std::cout << "Got plot dirname = " << plotDirectoryName << std::endl;
00117 
00118 
00119 }
00120 
00121 
00122 GeneralHLTOffline::~GeneralHLTOffline()
00123 {
00124  
00125    // do anything here that needs to be done at desctruction time
00126    // (e.g. close files, deallocate resources etc.)
00127 
00128 }
00129 
00130 
00131 //
00132 // member functions
00133 //
00134 
00135 // ------------ method called for each event  ------------
00136 void
00137 GeneralHLTOffline::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00138 {
00139    using namespace edm;
00140    using std::string;
00141 
00142    if (debugPrint) std::cout << "Inside analyze" << std::endl;
00143 
00144     // Access Trigger Results
00145    edm::Handle<edm::TriggerResults> triggerResults;
00146    iEvent.getByLabel(InputTag("TriggerResults","", "HLT"), triggerResults);
00147    
00148    if (!triggerResults.isValid()) {
00149      if (debugPrint) std::cout << "Trigger results not valid" << std::endl;
00150      return; }
00151 
00152     if (debugPrint) std::cout << "Found triggerResults" << std::endl;
00153 
00154    edm::Handle<trigger::TriggerEvent>         aodTriggerEvent;   
00155    iEvent.getByLabel(InputTag("hltTriggerSummaryAOD", "", "HLT"), aodTriggerEvent);
00156    
00157    if ( !aodTriggerEvent.isValid() ) { 
00158      if (debugPrint) std::cout << "No AOD trigger summary found! Returning..."; 
00159      return; 
00160    }
00161 
00162   std::vector<std::string> nameStreams = hltConfig_.streamNames();
00163 
00164   const TriggerObjectCollection objects = aodTriggerEvent->getObjects();
00165 
00166   bool streamAfound = false;
00167   int i = 0;
00168 
00169   for (vector<string>::iterator streamName = nameStreams.begin(); 
00170        streamName != nameStreams.end(); ++streamName) {
00171     
00172     if  (hltConfig_.streamName(i) == "A") {
00173       if (debugPrint) std::cout << " Stream A not found " << std::endl;
00174       streamAfound = true;
00175     }
00176     else
00177       if (debugPrint) std::cout << " Stream A found " << std::endl;
00178     i++;
00179   }
00180 
00181    if (streamAfound) {
00182      vector<string> datasetNames =  hltConfig_.streamContent("A");
00183      // Loop over PDs
00184      for (unsigned int iPD = 0; iPD < datasetNames.size(); iPD++) { 
00185        
00186        //     if (datasetNames[iPD] != "SingleMu" && datasetNames[iPD] != "SingleElectron" && datasetNames[iPD] != "Jet") continue;  
00187        
00188        unsigned int keyTracker[1000]; // Array to eliminate double counts by tracking what Keys have already been fired
00189        for(unsigned int irreproduceableIterator = 0; irreproduceableIterator < 1000; irreproduceableIterator++) {
00190          keyTracker[irreproduceableIterator] = 1001;
00191        }
00192        // Loop over Paths in each PD
00193        for (unsigned int iPath = 0; iPath < PDsVectorPathsVector[iPD].size(); iPath++) { //Andrew - where does PDsVectorPathsVector get defined?
00194          
00195          std::string pathName = PDsVectorPathsVector[iPD][iPath];
00196          
00197          if (debugPrint) std::cout << "Looking at path " << pathName << std::endl;
00198          
00199          unsigned int index = hltConfig_.triggerIndex(pathName);
00200          
00201          if (debugPrint) std::cout << "Index = " << index << " triggerResults->size() = " << triggerResults->size() << std::endl;
00202          
00203          if (index < triggerResults->size()) {
00204            if(triggerResults->accept(index)) {
00205              
00206              if (debugPrint) std::cout << "We fired path " << pathName << std::endl;
00207              
00208              // look up module labels for this path
00209              
00210              vector<std::string> modulesThisPath = hltConfig_.moduleLabels(pathName);
00211              
00212              if (debugPrint) std::cout << "Looping over module labels " << std::endl;
00213              
00214              // Loop backward through module names
00215              for ( int iModule = (modulesThisPath.size()-1); iModule >= 0; iModule--) {
00216                
00217                if (debugPrint) std::cout << "Module name is " << modulesThisPath[iModule] << std::endl;
00218                
00219                // check to see if you have savetags information
00220                if (hltConfig_.saveTags(modulesThisPath[iModule])) {
00221                  
00222                  if (debugPrint) std::cout << "For path " << pathName << " this module " << modulesThisPath[iModule] <<" is a saveTags module of type " << hltConfig_.moduleType(modulesThisPath[iModule]) << std::endl;
00223                  
00224                  if (hltConfig_.moduleType(modulesThisPath[iModule]) == "HLTLevel1GTSeed") break;
00225                  
00226                  InputTag moduleWhoseResultsWeWant(modulesThisPath[iModule], "", "HLT");
00227                  
00228                  unsigned int indexOfModuleInAodTriggerEvent = aodTriggerEvent->filterIndex(moduleWhoseResultsWeWant);
00229                  
00230                  if ( indexOfModuleInAodTriggerEvent < aodTriggerEvent->sizeFilters() ) {
00231                    const Keys &keys = aodTriggerEvent->filterKeys( indexOfModuleInAodTriggerEvent );
00232                    if (debugPrint) std::cout << "Got Keys for index " << indexOfModuleInAodTriggerEvent <<", size of keys is " << keys.size() << std::endl;
00233                    
00234                    for ( size_t iKey = 0; iKey < keys.size(); iKey++ ) {
00235                      TriggerObject foundObject = objects[keys[iKey]];
00236                      bool first_count = false;
00237                      
00238                      if(keyTracker[iKey] != iKey) first_count = true;
00239                      
00240                      if (debugPrint || outputPrint) std::cout << "This object has id (pt, eta, phi) = "
00241                                                               << " " << foundObject.id() << " " 
00242                                                               << std::setw(10) << foundObject.pt()
00243                                                               << ", " << std::setw(10) << foundObject.eta() 
00244                                                               << ", " << std::setw(10) << foundObject.phi()
00245                                                               << "    for path = " << std::setw(20) << pathName
00246                                                               << " module " << std::setw(40) << modulesThisPath[iModule]
00247                                                               << " iKey " << iKey << std::endl;
00248                      
00249                      fillHltMatrix(datasetNames[iPD],pathName,foundObject.eta(),foundObject.phi(),first_count);
00250                      
00251                      keyTracker[iKey] = iKey;
00252                      
00253                    }// end for each key               
00254                  }// end if filter in aodTriggerEvent
00255                  
00256                  
00257                  // OK, we found the last module. No need to look at the others.
00258                  // get out of the loop
00259                  
00260                  break;
00261                }// end if saveTags
00262              }//end Loop backward through module names   
00263            }// end if(triggerResults->accept(index))
00264          }// end if (index < triggerResults->size())
00265        }// end Loop over Paths in each PD
00266      }//end Loop over PDs
00267    }
00268 
00269    
00270 
00271 }
00272 
00273 
00274 // ------------ method called once each job just before starting event loop  ------------
00275 void 
00276 GeneralHLTOffline::beginJob()
00277 {using namespace edm;
00278  using std::string;
00279 
00280 
00281   
00282   if (debugPrint) std::cout << "Inside begin job" << std::endl; 
00283 
00284   dbe = Service<DQMStore>().operator->();
00285  
00286  if (dbe) {
00287 
00288     dbe->setCurrentFolder(plotDirectoryName);
00289     
00290   }
00291 
00292     
00293 }
00294 
00295 // ------------ method called once each job just after ending the event loop  ------------
00296 void 
00297 GeneralHLTOffline::endJob() 
00298 {
00299 }
00300 
00301 // ------------ method called when starting to processes a run  ------------
00302 void 
00303 GeneralHLTOffline::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup)
00304 {
00305 
00306   if (debugPrint) std::cout << "Inside beginRun" << std::endl;
00307 
00308 
00309   bool changed = true;
00310   if (hltConfig_.init(iRun, iSetup, "HLT", changed)) {
00311     if(debugPrint)
00312       if(debugPrint) std::cout << "HLT config with process name " 
00313                 << "HLT" << " successfully extracted" << std::endl;
00314   } else {
00315     if (debugPrint)
00316       if (debugPrint) std::cout << "Warning, didn't find process HLT" << std::endl;
00317   }
00318 
00319   if (debugPrint) std::cout << " About to access stream A content " << std::endl;
00320 
00321   std::vector<std::string> nameStreams = hltConfig_.streamNames();
00322 
00323   bool streamAfound = false;
00324   int i = 0;
00325   for (vector<string>::iterator streamName = nameStreams.begin(); 
00326        streamName != nameStreams.end(); ++streamName) {
00327     
00328     if  (hltConfig_.streamName(i) == "A") {
00329       if (debugPrint) std::cout << " Stream A found " << std::endl;
00330       streamAfound = true;
00331     }
00332     else
00333       if (debugPrint) std::cout << " Stream A not found " << std::endl;
00334     i++;
00335   }
00336 
00337 
00338   if (streamAfound) {
00339     vector<string> datasetNames =  hltConfig_.streamContent("A");
00340   
00341     if (debugPrint) std::cout << " Size of Stream A dataset " << datasetNames.size() << std::endl;
00342 
00343     for (unsigned int i=0;i<datasetNames.size();i++) {
00344       
00345       if (debugPrint) std::cout << "This is dataset " << datasetNames[i] <<std::endl;
00346       
00347       vector<string> datasetPaths = hltConfig_.datasetContent(datasetNames[i]);
00348 
00349       if (debugPrint) std::cout << "datasetPaths.size() = " << datasetPaths.size() << std::endl;
00350       
00351       PDsVectorPathsVector.push_back(datasetPaths);
00352       
00353       if (debugPrint) std::cout <<"Found PD: " << datasetNames[i]  << std::endl;     
00354       setupHltMatrix(datasetNames[i],i);   
00355       
00356     }// end of loop over dataset names
00357   }
00358 
00359 }// end of beginRun
00360 
00361 // ------------ method called when ending the processing of a run  ------------
00362 void GeneralHLTOffline::endRun(edm::Run const&, edm::EventSetup const&)
00363 {
00364 }
00365 
00366 void GeneralHLTOffline::setupHltMatrix(std::string label, int iPD) {
00367 
00368 std::string h_name;
00369 std::string h_title;
00370 std::string h_name_1dEta;
00371 std::string h_name_1dPhi;
00372 std::string h_title_1dEta;
00373 std::string h_title_1dPhi;
00374 std::string h_name_1dEtaPath;
00375 std::string h_name_1dPhiPath;
00376 std::string h_title_1dEtaPath;
00377 std::string h_title_1dPhiPath;
00378 std::string pathName;
00379 std::string PD_Folder;
00380 std::string Path_Folder;
00381 std::string HLTMenu;
00382 
00383 PD_Folder = TString("HLT/GeneralHLTOffline");
00384 if (label != "SingleMu" && label != "SingleElectron" && label != "Jet")  PD_Folder = TString("HLT/GeneralHLTOffline/"+label); 
00385 
00386 dbe->setCurrentFolder(PD_Folder.c_str());
00387 HLTMenu = hltConfig_.tableName();
00388 dbe->bookString("hltMenuName",HLTMenu.c_str());
00389 
00390 
00391 
00392 
00393 h_name = "HLT_"+label+"_EtaVsPhi";
00394 h_title = "HLT_"+label+"_EtaVsPhi";
00395 h_name_1dEta = "HLT_"+label+"_1dEta";
00396 h_name_1dPhi = "HLT_"+label+"_1dPhi";
00397 h_title_1dEta = label+" Occupancy Vs Eta";
00398 h_title_1dPhi = label+" Occupancy Vs Phi";
00399 //Int_t numBinsEta = 12;
00400 //Int_t numBinsPhi = 8;
00401 Int_t numBinsEta = 30;
00402 Int_t numBinsPhi = 34;
00403 Int_t numBinsEtaFine = 60;
00404 Int_t numBinsPhiFine = 66;
00405 Double_t EtaMax = 2.610;
00406 Double_t PhiMax = 17.0*TMath::Pi()/16.0;
00407 Double_t PhiMaxFine = 33.0*TMath::Pi()/32.0;
00408 
00409   //Double_t eta_bins[] = {-2.610,-2.349,-2.088,-1.827,-1.740,-1.653,-1.566,-1.479,-1.392,-1.305,-1.044,-0.783,-0.522,-0.261,0,0.261,0.522,0.783,1.044,1.305,1.392,1.479,1.566,1.653,1.740,1.827,2.088,2.349,2.610}; //Has narrower bins in Barrel/Endcap border region
00410 //Double_t eta_bins[] = {-2.610,-2.175,-1.740,-1.305,-0.870,-0.435,0,0.435,0.870,1.305,1.740,2.175,2.610};
00411 //Double_t phi_bins[] = {-TMath::Pi(),-3*TMath::Pi()/4,-TMath::Pi()/2,-TMath::Pi()/4,0,TMath::Pi()/4,TMath::Pi()/2,3*TMath::Pi()/4,TMath::Pi()};
00412 
00413  TH2F * hist_EtaVsPhi = new TH2F(h_name.c_str(),h_title.c_str(),numBinsEta,-EtaMax,EtaMax,numBinsPhi,-PhiMax,PhiMax);
00414  TH1F * hist_1dEta = new TH1F(h_name_1dEta.c_str(),h_title_1dEta.c_str(),numBinsEtaFine,-EtaMax,EtaMax);
00415  TH1F * hist_1dPhi = new TH1F(h_name_1dPhi.c_str(),h_title_1dPhi.c_str(),numBinsPhiFine,-PhiMaxFine,PhiMaxFine);
00416 
00417  hist_EtaVsPhi->SetMinimum(0);
00418  hist_1dEta->SetMinimum(0);
00419  hist_1dPhi->SetMinimum(0);
00420 
00421  // Do not comment out these pointers since it is the booking that is the work here. 
00422  // MonitorElement * ME_EtaVsPhi = dbe->book2D(h_name.c_str(),hist_EtaVsPhi);
00423  dbe->book2D(h_name.c_str(),hist_EtaVsPhi);
00424  // MonitorElement * ME_1dEta = dbe->book1D(h_name_1dEta.c_str(),hist_1dEta);
00425  dbe->book1D(h_name_1dEta.c_str(),hist_1dEta);
00426  // MonitorElement * ME_1dPhi = dbe->book1D(h_name_1dPhi.c_str(),hist_1dPhi);
00427  dbe->book1D(h_name_1dPhi.c_str(),hist_1dPhi);
00428 
00429   for (unsigned int iPath = 0; iPath < PDsVectorPathsVector[iPD].size(); iPath++) { 
00430     pathName = PDsVectorPathsVector[iPD][iPath];
00431     h_name_1dEtaPath = "HLT_"+pathName+"_1dEta";
00432     h_name_1dPhiPath = "HLT_"+pathName+"_1dPhi";
00433     h_title_1dEtaPath = pathName+" Occupancy Vs Eta";
00434     h_title_1dPhiPath = pathName+"Occupancy Vs Phi";
00435     Path_Folder = TString("HLT/GeneralHLTOffline/"+label+"/Paths");
00436     dbe->setCurrentFolder(Path_Folder.c_str());
00437 
00438     // Do not comment out these pointers since it is the booking that is the work here. 
00439     //     MonitorElement * ME_1dEta = dbe->book1D(h_name_1dEtaPath.c_str(),h_title_1dEtaPath.c_str(),numBinsEtaFine,-EtaMax,EtaMax);
00440     dbe->book1D(h_name_1dEtaPath.c_str(),h_title_1dEtaPath.c_str(),numBinsEtaFine,-EtaMax,EtaMax);
00441      //     MonitorElement * ME_1dPhi = dbe->book1D(h_name_1dPhiPath.c_str(),h_title_1dPhiPath.c_str(),numBinsPhiFine,-PhiMaxFine,PhiMaxFine);
00442     dbe->book1D(h_name_1dPhiPath.c_str(),h_title_1dPhiPath.c_str(),numBinsPhiFine,-PhiMaxFine,PhiMaxFine);
00443   
00444     if (debugPrint) std::cout << "book1D for " << pathName << std::endl;
00445   }
00446    
00447  if (debugPrint) std::cout << "Success setupHltMatrix( " << label << " , " << iPD << " )" << std::cout;
00448 } //End setupHltMatrix
00449 
00450 void GeneralHLTOffline::fillHltMatrix(std::string label, std::string path,double Eta, double Phi, bool first_count) {
00451 
00452   if (debugPrint) std::cout << "Inside fillHltMatrix( " << label << " , " << path << " ) " << std::endl;
00453 
00454   std::string fullPathToME;
00455   std::string fullPathToME1dEta;
00456   std::string fullPathToME1dPhi;
00457   std::string fullPathToME1dEtaPath;
00458   std::string fullPathToME1dPhiPath;
00459 
00460  fullPathToME = "HLT/GeneralHLTOffline/HLT_"+label+"_EtaVsPhi"; 
00461  fullPathToME1dEta = "HLT/GeneralHLTOffline/HLT_"+label+"_1dEta";
00462  fullPathToME1dPhi = "HLT/GeneralHLTOffline/HLT_"+label+"_1dPhi";
00463 
00464 if (label != "SingleMu" && label != "SingleElectron" && label != "Jet") {
00465  fullPathToME = "HLT/GeneralHLTOffline/"+label+"/HLT_"+label+"_EtaVsPhi"; 
00466  fullPathToME1dEta = "HLT/GeneralHLTOffline/"+label+"/HLT_"+label+"_1dEta";
00467  fullPathToME1dPhi = "HLT/GeneralHLTOffline/"+label+"/HLT_"+label+"_1dPhi";
00468  }
00469  
00470  fullPathToME1dEtaPath = "HLT/GeneralHLTOffline/"+label+"/Paths/HLT_"+path+"_1dEta";
00471  fullPathToME1dPhiPath = "HLT/GeneralHLTOffline/"+label+"/Paths/HLT_"+path+"_1dPhi";
00472 
00473   if (debugPrint) std::cout << "fullPathToME = " << std::endl;
00474 
00475   MonitorElement * ME_2d = dbe->get(fullPathToME);
00476   MonitorElement * ME_1dEta = dbe->get(fullPathToME1dEta);
00477   MonitorElement * ME_1dPhi = dbe->get(fullPathToME1dPhi);  
00478   MonitorElement * ME_1dEtaPath = dbe->get(fullPathToME1dEtaPath);
00479   MonitorElement * ME_1dPhiPath = dbe->get(fullPathToME1dPhiPath);
00480 
00481   if (debugPrint) std::cout << "MonitorElement * " << std::endl;
00482 
00483   TH2F * hist_2d = ME_2d->getTH2F();
00484   TH1F * hist_1dEta = ME_1dEta->getTH1F();
00485   TH1F * hist_1dPhi = ME_1dPhi->getTH1F();
00486   TH1F * hist_1dEtaPath = ME_1dEtaPath->getTH1F();
00487   TH1F * hist_1dPhiPath = ME_1dPhiPath->getTH1F();
00488 
00489   if (debugPrint) std::cout << "TH2F *" << std::endl;
00490 
00491   //int i=2;
00492   //if (Eta>1.305 && Eta<1.872) i=0;
00493   //if (Eta<-1.305 && Eta>-1.872) i=0;
00494   //for (int ii=i; ii<3; ++ii) hist_2d->Fill(Eta,Phi); //Scales narrow bins in Barrel/Endcap border region
00495 
00496   if(first_count) {
00497     hist_1dEta->Fill(Eta);
00498     hist_1dPhi->Fill(Phi); 
00499     hist_2d->Fill(Eta,Phi); }
00500     hist_1dEtaPath->Fill(Eta); 
00501     hist_1dPhiPath->Fill(Phi);
00502 
00503  if (debugPrint) std::cout << "hist->Fill" << std::endl;
00504 
00505 } //End fillHltMatrix
00506 
00507 // ------------ method called when starting to processes a luminosity block  ------------
00508 void GeneralHLTOffline::beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&)
00509 {
00510 }
00511 
00512 // ------------ method called when ending the processing of a luminosity block  ------------
00513 void 
00514 GeneralHLTOffline::endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&)
00515 {
00516 }
00517 
00518 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
00519 // void
00520 // GeneralHLTOffline::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00521 //   //The following says we do not know what parameters are allowed so do no validation
00522 //   // Please change this to state exactly what you do use, even if it is no parameters
00523 //   edm::ParameterSetDescription desc;
00524 //   desc.setUnknown();
00525 //   descriptions.addDefault(desc);
00526 // }
00527 
00528 //define this as a plug-in for online (not offline ?)
00529 DEFINE_FWK_MODULE(GeneralHLTOffline);