CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Functions
FWLiteWithPythonConfig.cc File Reference
#include <TH1F.h>
#include <TROOT.h>
#include <TFile.h>
#include <TSystem.h>
#include "DataFormats/FWLite/interface/Event.h"
#include "DataFormats/Common/interface/Handle.h"
#include "FWCore/FWLite/interface/FWLiteEnabler.h"
#include "DataFormats/FWLite/interface/InputSource.h"
#include "DataFormats/FWLite/interface/OutputFiles.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSetReader/interface/ParameterSetReader.h"
#include "DataFormats/MuonReco/interface/Muon.h"
#include "DataFormats/PatCandidates/interface/Muon.h"
#include "PhysicsTools/FWLite/interface/TFileService.h"

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Function Documentation

int main ( int  argc,
char *  argv[] 
)

===============================================================================================================================================================================================


variant2: for each run define phi-averaged A for normalization channel (Dref,16) and then, divide Rijk on it, i.e. get RRijk

eta=27

eta=25

eta=23

eta=22

eta=21

eta=26

eta=24

eta=19

eta=17

eta=25

eta=23

eta=22

eta=21

eta=26

eta=24

eta=20

eta=19

eta=18

eta=27 L1=1

eta=25 L1=1

eta=23 L1=1

eta=22 L1=1

eta=21 L1=1

eta=29 L1=1

eta=26 L1=1

eta=24 L1=1

eta=20 L1=1

eta=19 L1=1

eta=18 L1=1

eta=17 L1=1

eta=28 L7=1

eta=27 L7=1

eta=25 L7=1

eta=23 L7=1

eta=22 L7=1

eta=21 L7=1

eta=26 L7=1

eta=24 L7=1

eta=20 L7=1

eta=19 L7=1

eta=18 L7=1

eta=17 L7=1

eta=27

eta=25

eta=23

eta=22

eta=21

eta=26

eta=24

eta=19

eta=17

eta=25

eta=23

eta=22

eta=21

eta=26

eta=24

eta=20

eta=19

eta=18

eta=27 L1=1

eta=25 L1=1

eta=23 L1=1

eta=22 L1=1

eta=21 L1=1

eta=26 L1=1

eta=24 L1=1

eta=20 L1=1

eta=19 L1=1

eta=18 L1=1

eta=17 L1=1

eta=28 L7=1

eta=27 L7=1

eta=25 L7=1

eta=23 L7=1

eta=22 L7=1

eta=21 L7=1

eta=26 L7=1

eta=24 L7=1

eta=20 L7=1

eta=19 L7=1

eta=18 L7=1

eta=17 L7=1

eta=27

eta=28

errA with average Amplitudes

errA with average Amplitudes

errA with average Amplitudes

errA with average Amplitudes

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

RBX:

errA with average Amplitudes

errA with average Amplitudes

errA with average Amplitudes

errA with average Amplitudes

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

RBX:

errA with average Amplitudes

errA with average Amplitudes

errA with average Amplitudes

errA with average Amplitudes

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

Summed Amplitude Plots:

RBX:

Prepare maps of good/bad channels:

Prepare maps of good/bad channels:

Prepare maps of good/bad channels:

Prepare maps of good/bad channels:

Definition at line 19 of file FWLiteWithPythonConfig.cc.

References fwlite::Event::atEnd(), gather_cfg::cout, DeadROC_duringRun::dir, FWLiteEnabler::enable(), makeMEIFBenchmarkPlots::ev, beamvalidation::exit(), fwlite::OutputFiles::file(), compareTotals::fs, edm::ParameterSet::getParameter(), TFileDirectory::make(), ResonanceBuilder::mass, TFileDirectory::mkdir(), dumpRecoGeometry_cfg::Muon, patZpeak::muons, LaserDQM_cfg::process, edm::readPSetsFrom(), TrackerOfflineValidation_Standalone_cff::TFileService, and fwlite::Event::toBegin().

19  {
20  // define what muon you are using; this is necessary as FWLite is not
21  // capable of reading edm::Views
22  using reco::Muon;
23 
24  // ----------------------------------------------------------------------
25  // First Part:
26  //
27  // * enable FWLite
28  // * book the histograms of interest
29  // * open the input file
30  // ----------------------------------------------------------------------
31 
32  // load framework libraries
33  gSystem->Load("libFWCoreFWLite");
35 
36  // parse arguments
37  if (argc < 2) {
38  std::cout << "Usage : " << argv[0] << " [parameters.py]" << std::endl;
39  return 0;
40  }
41 
42  if (!edm::readPSetsFrom(argv[1])->existsAs<edm::ParameterSet>("process")) {
43  std::cout << " ERROR: ParametersSet 'process' is missing in your configuration file" << std::endl;
44  exit(0);
45  }
46  // get the python configuration
47  const edm::ParameterSet& process = edm::readPSetsFrom(argv[1])->getParameter<edm::ParameterSet>("process");
48  fwlite::InputSource inputHandler_(process);
49  fwlite::OutputFiles outputHandler_(process);
50 
51  // now get each parameter
52  const edm::ParameterSet& ana = process.getParameter<edm::ParameterSet>("muonAnalyzer");
53  edm::InputTag muons_(ana.getParameter<edm::InputTag>("muons"));
54 
55  // book a set of histograms
56  fwlite::TFileService fs = fwlite::TFileService(outputHandler_.file());
57  TFileDirectory dir = fs.mkdir("analyzeBasicPat");
58  TH1F* muonPt_ = dir.make<TH1F>("muonPt", "pt", 100, 0., 300.);
59  TH1F* muonEta_ = dir.make<TH1F>("muonEta", "eta", 100, -3., 3.);
60  TH1F* muonPhi_ = dir.make<TH1F>("muonPhi", "phi", 100, -5., 5.);
61  TH1F* mumuMass_ = dir.make<TH1F>("mumuMass", "mass", 90, 30., 120.);
62 
63  // loop the events
64  int ievt = 0;
65  int maxEvents_(inputHandler_.maxEvents());
66  for (unsigned int iFile = 0; iFile < inputHandler_.files().size(); ++iFile) {
67  // open input file (can be located on castor)
68  TFile* inFile = TFile::Open(inputHandler_.files()[iFile].c_str());
69  if (inFile) {
70  // ----------------------------------------------------------------------
71  // Second Part:
72  //
73  // * loop the events in the input file
74  // * receive the collections of interest via fwlite::Handle
75  // * fill the histograms
76  // * after the loop close the input file
77  // ----------------------------------------------------------------------
78  fwlite::Event ev(inFile);
79  for (ev.toBegin(); !ev.atEnd(); ++ev, ++ievt) {
80  edm::EventBase const& event = ev;
81  // break loop if maximal number of events is reached
82  if (maxEvents_ > 0 ? ievt + 1 > maxEvents_ : false)
83  break;
84  // simple event counter
85  if (inputHandler_.reportAfter() != 0 ? (ievt > 0 && ievt % inputHandler_.reportAfter() == 0) : false)
86  std::cout << " processing event: " << ievt << std::endl;
87 
88  // Handle to the muon collection
90  event.getByLabel(muons_, muons);
91 
92  // loop muon collection and fill histograms
93  for (std::vector<Muon>::const_iterator mu1 = muons->begin(); mu1 != muons->end(); ++mu1) {
94  muonPt_->Fill(mu1->pt());
95  muonEta_->Fill(mu1->eta());
96  muonPhi_->Fill(mu1->phi());
97  if (mu1->pt() > 20 && fabs(mu1->eta()) < 2.1) {
98  for (std::vector<Muon>::const_iterator mu2 = muons->begin(); mu2 != muons->end(); ++mu2) {
99  if (mu2 > mu1) { // prevent double conting
100  if (mu1->charge() * mu2->charge() < 0) { // check only muon pairs of unequal charge
101  if (mu2->pt() > 20 && fabs(mu2->eta()) < 2.1) {
102  mumuMass_->Fill((mu1->p4() + mu2->p4()).mass());
103  }
104  }
105  }
106  }
107  }
108  }
109  }
110  // close input file
111  inFile->Close();
112  }
113  // break loop if maximal number of events is reached:
114  // this has to be done twice to stop the file loop as well
115  if (maxEvents_ > 0 ? ievt + 1 > maxEvents_ : false)
116  break;
117  }
118  return 0;
119 }
static void enable()
enable automatic library loading
std::unique_ptr< edm::ParameterSet > readPSetsFrom(std::string const &fileOrString)
T * make(const Args &...args) const
make new ROOT object
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
tuple argc
Definition: dir2webdir.py:39
tuple muons
Definition: patZpeak.py:41
tuple cout
Definition: gather_cfg.py:144
tuple process
Definition: LaserDQM_cfg.py:3