00001
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112 #include "FWCore/Framework/interface/EDLooper.h"
00113 #include "FWCore/Framework/interface/LooperFactory.h"
00114 #include "FWCore/Utilities/interface/InputTag.h"
00115 #include "FWCore/Framework/interface/Frameworkfwd.h"
00116 #include "FWCore/Framework/interface/Event.h"
00117 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00118
00119 #include <CLHEP/Vector/LorentzVector.h>
00120 #include <vector>
00121
00122 #include "FWCore/Framework/interface/EventSetup.h"
00123 #include "FWCore/Framework/interface/ESHandle.h"
00124 #include "FWCore/ServiceRegistry/interface/Service.h"
00125
00126 #include "MuScleFitBase.h"
00127 #include "Histograms.h"
00128 #include "MuScleFitPlotter.h"
00129 #include "MuonAnalysis/MomentumScaleCalibration/interface/Functions.h"
00130 #include "MuonAnalysis/MomentumScaleCalibration/interface/RootTreeHandler.h"
00131 #include "MuScleFitMuonSelector.h"
00132
00133 #include "DataFormats/TrackReco/interface/Track.h"
00134 #include "DataFormats/MuonReco/interface/Muon.h"
00135 #include "DataFormats/MuonReco/interface/CaloMuon.h"
00136 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00137
00138 #include "RecoMuon/TrackingTools/interface/MuonPatternRecoDumper.h"
00139 #include "RecoMuon/TrackingTools/interface/MuonServiceProxy.h"
00140
00141 #include "DataFormats/PatCandidates/interface/Muon.h"
00142 #include "DataFormats/PatCandidates/interface/CompositeCandidate.h"
00143 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00144 #include "DataFormats/Common/interface/TriggerResults.h"
00145 #include "FWCore/Common/interface/TriggerNames.h"
00146
00147 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
00148
00149 #include "HepPDT/defs.h"
00150 #include "HepPDT/TableBuilder.hh"
00151 #include "HepPDT/ParticleDataTable.hh"
00152
00153 #include "HepMC/GenParticle.h"
00154 #include "HepMC/GenEvent.h"
00155
00156 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00157 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00158 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
00159 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00160
00161 #include "TFile.h"
00162 #include "TTree.h"
00163 #include "TMinuit.h"
00164
00165
00166
00167
00168
00169 #ifdef USE_CALLGRIND
00170 #include "valgrind/callgrind.h"
00171 #endif
00172
00173
00174
00175
00176
00177
00178 namespace edm {
00179 class ParameterSet;
00180 class Event;
00181 class EventSetup;
00182 }
00183
00184
00185 class MuScleFit: public edm::EDLooper, MuScleFitBase
00186 {
00187 public:
00188
00189
00190 MuScleFit( const edm::ParameterSet& pset );
00191
00192
00193
00194 virtual ~MuScleFit();
00195
00196
00197
00198 void beginOfJobInConstructor();
00199
00200
00201 virtual void endOfJob();
00202
00203 virtual void startingNewLoop( unsigned int iLoop );
00204
00205 virtual edm::EDLooper::Status endOfLoop( const edm::EventSetup& eventSetup, unsigned int iLoop );
00206 virtual void endOfFastLoop( const unsigned int iLoop );
00207
00208 virtual edm::EDLooper::Status duringLoop( const edm::Event & event, const edm::EventSetup& eventSetup );
00213 virtual void duringFastLoop();
00214
00215 template<typename T>
00216 std::vector<reco::LeafCandidate> fillMuonCollection( const std::vector<T>& tracks );
00217 private:
00218
00219 protected:
00224 void selectMuons(const edm::Event & event);
00230 void selectMuons(const int maxEvents, const TString & treeFileName);
00231
00233 template<typename T>
00234 void takeSelectedMuonType(const T & muon, std::vector<reco::Track> & tracks);
00236 bool selGlobalMuon(const pat::Muon* aMuon);
00237 bool selTrackerMuon(const pat::Muon* aMuon);
00238
00240 bool checkDeltaR( reco::Particle::LorentzVector& genMu, reco::Particle::LorentzVector& recMu );
00242 void fillComparisonHistograms( const reco::Particle::LorentzVector & genMu, const reco::Particle::LorentzVector & recoMu, const std::string & inputName, const int charge );
00243
00245 void applySmearing( reco::Particle::LorentzVector & mu );
00247 void applyBias( reco::Particle::LorentzVector & mu, const int charge );
00248
00253 void checkParameters();
00254
00255 MuonServiceProxy *theService;
00256
00257
00258
00259 int numberOfSimTracks;
00260 int numberOfSimMuons;
00261 int numberOfSimVertices;
00262 int numberOfEwkZ;
00263
00264 bool ifHepMC;
00265 bool ifGenPart;
00266
00267
00268
00269 double minResMass_hwindow[6];
00270 double maxResMass_hwindow[6];
00271
00272
00273
00274 unsigned int maxLoopNumber;
00275 unsigned int loopCounter;
00276
00277 bool fastLoop;
00278
00279 MuScleFitPlotter *plotter;
00280
00281
00282
00283 reco::Particle::LorentzVector recMu1, recMu2;
00284 int iev;
00285 int totalEvents_;
00286
00287 bool compareToSimTracks_;
00288 edm::InputTag simTracksCollection_;
00289 bool PATmuons_;
00290 std::string genParticlesName_;
00291
00292
00293 std::string inputRootTreeFileName_;
00294
00295 std::string outputRootTreeFileName_;
00296
00297 int maxEventsFromRootTree_;
00298
00299 std::string triggerResultsLabel_;
00300 std::string triggerResultsProcess_;
00301 std::vector<std::string> triggerPath_;
00302 bool negateTrigger_;
00303 bool saveAllToTree_;
00304
00305 std::auto_ptr<MuScleFitMuonSelector> muonSelector_;
00306 };
00307
00308 template<typename T>
00309 std::vector<reco::LeafCandidate> MuScleFit::fillMuonCollection( const std::vector<T>& tracks )
00310 {
00311 std::vector<reco::LeafCandidate> muons;
00312 typename std::vector<T>::const_iterator track;
00313 for( track = tracks.begin(); track != tracks.end(); ++track ) {
00314 reco::Particle::LorentzVector mu;
00315 mu = reco::Particle::LorentzVector(track->px(),track->py(),track->pz(),
00316 sqrt(track->p()*track->p() + MuScleFitUtils::mMu2));
00317
00318
00319 MuScleFitUtils::goodmuon++;
00320 if (debug_>0)
00321 std::cout <<std::setprecision(9)<< "Muon #" << MuScleFitUtils::goodmuon
00322 << ": initial value Pt = " << mu.Pt() << std::endl;
00323
00324 applySmearing(mu);
00325 applyBias(mu, track->charge());
00326
00327 reco::LeafCandidate muon(track->charge(),mu);
00328
00329
00330 muons.push_back (muon);
00331 }
00332 return muons;
00333 }
00334
00335 template<typename T>
00336 void MuScleFit::takeSelectedMuonType(const T & muon, std::vector<reco::Track> & tracks)
00337 {
00338
00339
00340
00341 if(muon->isGlobalMuon() && theMuonType_==1)
00342 tracks.push_back(*(muon->globalTrack()));
00343 else if(muon->isStandAloneMuon() && theMuonType_==2)
00344 tracks.push_back(*(muon->outerTrack()));
00345 else if(muon->isTrackerMuon() && theMuonType_==3)
00346 tracks.push_back(*(muon->innerTrack()));
00347
00348 else if( theMuonType_ == 10 && !(muon->isStandAloneMuon()) )
00349 tracks.push_back(*(muon->innerTrack()));
00350 else if( theMuonType_ == 11 && muon->isGlobalMuon() )
00351 tracks.push_back(*(muon->innerTrack()));
00352 else if( theMuonType_ == 13 && muon->isTrackerMuon() )
00353 tracks.push_back(*(muon->innerTrack()));
00354 }
00355
00356
00357
00358
00359 MuScleFit::MuScleFit( const edm::ParameterSet& pset ) :
00360 MuScleFitBase( pset ),
00361 totalEvents_(0)
00362 {
00363 MuScleFitUtils::debug = debug_;
00364 if (debug_>0) std::cout << "[MuScleFit]: Constructor" << std::endl;
00365
00366 if ((theMuonType_<-4 || theMuonType_>5) && theMuonType_<10) {
00367 std::cout << "[MuScleFit]: Unknown muon type! Aborting." << std::endl;
00368 abort();
00369 }
00370
00371 loopCounter = 0;
00372
00373
00374
00375 minResMass_hwindow[0] = 71.1876;
00376 maxResMass_hwindow[0] = 111.188;
00377 minResMass_hwindow[1] = 10.15;
00378 maxResMass_hwindow[1] = 10.55;
00379 minResMass_hwindow[2] = 9.8;
00380 maxResMass_hwindow[2] = 10.2;
00381 minResMass_hwindow[3] = 9.25;
00382 maxResMass_hwindow[3] = 9.65;
00383 minResMass_hwindow[4] = 3.58;
00384 maxResMass_hwindow[4] = 3.78;
00385 minResMass_hwindow[5] = 3.0;
00386 maxResMass_hwindow[5] = 3.2;
00387
00388
00389
00390 maxLoopNumber = pset.getUntrackedParameter<int>("maxLoopNumber", 2);
00391 fastLoop = pset.getUntrackedParameter<bool>("FastLoop", true);
00392
00393
00394 MuScleFitUtils::doResolFit = pset.getParameter<std::vector<int> >("doResolFit");
00395 MuScleFitUtils::doScaleFit = pset.getParameter<std::vector<int> >("doScaleFit");
00396 MuScleFitUtils::doCrossSectionFit = pset.getParameter<std::vector<int> >("doCrossSectionFit");
00397 MuScleFitUtils::doBackgroundFit = pset.getParameter<std::vector<int> >("doBackgroundFit");
00398
00399
00400
00401 int biasType = pset.getParameter<int>("BiasType");
00402 MuScleFitUtils::BiasType = biasType;
00403
00404 MuScleFitUtils::biasFunction = scaleFunctionVecService( biasType );
00405 int smearType = pset.getParameter<int>("SmearType");
00406 MuScleFitUtils::SmearType = smearType;
00407 MuScleFitUtils::smearFunction = smearFunctionService( smearType );
00408
00409
00410
00411 int resolFitType = pset.getParameter<int>("ResolFitType");
00412 MuScleFitUtils::ResolFitType = resolFitType;
00413 MuScleFitUtils::resolutionFunction = resolutionFunctionService( resolFitType );
00414 MuScleFitUtils::resolutionFunctionForVec = resolutionFunctionVecService( resolFitType );
00415 int scaleType = pset.getParameter<int>("ScaleFitType");
00416 MuScleFitUtils::ScaleFitType = scaleType;
00417 MuScleFitUtils::scaleFunction = scaleFunctionService( scaleType );
00418 MuScleFitUtils::scaleFunctionForVec = scaleFunctionVecService( scaleType );
00419
00420
00421
00422 MuScleFitUtils::parBias = pset.getParameter<std::vector<double> >("parBias");
00423 MuScleFitUtils::parSmear = pset.getParameter<std::vector<double> >("parSmear");
00424 MuScleFitUtils::parResol = pset.getParameter<std::vector<double> >("parResol");
00425 MuScleFitUtils::parResolStep = pset.getUntrackedParameter<std::vector<double> >("parResolStep", std::vector<double>());
00426 MuScleFitUtils::parResolMin = pset.getUntrackedParameter<std::vector<double> >("parResolMin", std::vector<double>());
00427 MuScleFitUtils::parResolMax = pset.getUntrackedParameter<std::vector<double> >("parResolMax", std::vector<double>());
00428 MuScleFitUtils::parScale = pset.getParameter<std::vector<double> >("parScale");
00429 MuScleFitUtils::parScaleStep = pset.getUntrackedParameter<std::vector<double> >("parScaleStep", std::vector<double>());
00430 MuScleFitUtils::parScaleMin = pset.getUntrackedParameter<std::vector<double> >("parScaleMin", std::vector<double>());
00431 MuScleFitUtils::parScaleMax = pset.getUntrackedParameter<std::vector<double> >("parScaleMax", std::vector<double>());
00432 MuScleFitUtils::parCrossSection = pset.getParameter<std::vector<double> >("parCrossSection");
00433 MuScleFitUtils::parBgr = pset.getParameter<std::vector<double> >("parBgr");
00434 MuScleFitUtils::parResolFix = pset.getParameter<std::vector<int> >("parResolFix");
00435 MuScleFitUtils::parScaleFix = pset.getParameter<std::vector<int> >("parScaleFix");
00436 MuScleFitUtils::parCrossSectionFix = pset.getParameter<std::vector<int> >("parCrossSectionFix");
00437 MuScleFitUtils::parBgrFix = pset.getParameter<std::vector<int> >("parBgrFix");
00438 MuScleFitUtils::parResolOrder = pset.getParameter<std::vector<int> >("parResolOrder");
00439 MuScleFitUtils::parScaleOrder = pset.getParameter<std::vector<int> >("parScaleOrder");
00440 MuScleFitUtils::parCrossSectionOrder = pset.getParameter<std::vector<int> >("parCrossSectionOrder");
00441 MuScleFitUtils::parBgrOrder = pset.getParameter<std::vector<int> >("parBgrOrder");
00442
00443 MuScleFitUtils::resfind = pset.getParameter<std::vector<int> >("resfind");
00444 MuScleFitUtils::FitStrategy = pset.getParameter<int>("FitStrategy");
00445
00446
00447
00448 MuScleFitUtils::speedup = pset.getParameter<bool>("speedup");
00449
00450
00451 compareToSimTracks_ = pset.getParameter<bool>("compareToSimTracks");
00452 simTracksCollection_ = pset.getUntrackedParameter<edm::InputTag>("SimTracksCollection", edm::InputTag("g4SimHits"));
00453
00454 triggerResultsLabel_ = pset.getUntrackedParameter<std::string>("TriggerResultsLabel");
00455 triggerResultsProcess_ = pset.getUntrackedParameter<std::string>("TriggerResultsProcess");
00456 triggerPath_ = pset.getUntrackedParameter<std::vector<std::string> >("TriggerPath");
00457 negateTrigger_ = pset.getUntrackedParameter<bool>("NegateTrigger", false);
00458 saveAllToTree_ = pset.getUntrackedParameter<bool>("SaveAllToTree", false);
00459
00460 PATmuons_ = pset.getUntrackedParameter<bool>("PATmuons", false);
00461 genParticlesName_ = pset.getUntrackedParameter<std::string>("GenParticlesName", "genParticles");
00462
00463
00464
00465 MuScleFitUtils::useProbsFile_ = pset.getUntrackedParameter<bool>("UseProbsFile", true);
00466
00467
00468 MuScleFitUtils::sherpa_ = pset.getUntrackedParameter<bool>("Sherpa", false);
00469
00470 MuScleFitUtils::rapidityBinsForZ_ = pset.getUntrackedParameter<bool>("RapidityBinsForZ", true);
00471
00472
00473 MuScleFitUtils::separateRanges_ = pset.getUntrackedParameter<bool>("SeparateRanges", true);
00474 MuScleFitUtils::maxMuonPt_ = pset.getUntrackedParameter<double>("MaxMuonPt", 100000000.);
00475 MuScleFitUtils::minMuonPt_ = pset.getUntrackedParameter<double>("MinMuonPt", 0.);
00476 MuScleFitUtils::minMuonEtaFirstRange_ = pset.getUntrackedParameter<double>("MinMuonEtaFirstRange", -6.);
00477 MuScleFitUtils::maxMuonEtaFirstRange_ = pset.getUntrackedParameter<double>("MaxMuonEtaFirstRange", 6.);
00478 MuScleFitUtils::minMuonEtaSecondRange_ = pset.getUntrackedParameter<double>("MinMuonEtaSecondRange", -100.);
00479 MuScleFitUtils::maxMuonEtaSecondRange_ = pset.getUntrackedParameter<double>("MaxMuonEtaSecondRange", 100.);
00480 MuScleFitUtils::deltaPhiMinCut_ = pset.getUntrackedParameter<double>("DeltaPhiMinCut", -100.);
00481 MuScleFitUtils::deltaPhiMaxCut_ = pset.getUntrackedParameter<double>("DeltaPhiMaxCut", 100.);
00482
00483 MuScleFitUtils::debugMassResol_ = pset.getUntrackedParameter<bool>("DebugMassResol", false);
00484
00485
00486
00487
00488 checkParameters();
00489
00490
00491
00492 if (MuScleFitUtils::SmearType>0) {
00493 std::cout << "[MuScleFit-Constructor]: Generating random values for smearing" << std::endl;
00494 TF1 G("G", "[0]*exp(-0.5*pow(x,2))", -5., 5.);
00495 double norm = 1/sqrt(2*TMath::Pi());
00496 G.SetParameter (0,norm);
00497 for (int i=0; i<10000; i++) {
00498 for (int j=0; j<7; j++) {
00499 MuScleFitUtils::x[j][i] = G.GetRandom();
00500 }
00501 }
00502 }
00503 MuScleFitUtils::goodmuon = 0;
00504
00505 if(theMuonType_ > 0 && theMuonType_ < 4) {
00506 MuScleFitUtils::MuonTypeForCheckMassWindow = theMuonType_-1;
00507 MuScleFitUtils::MuonType = theMuonType_-1;
00508 }
00509 else if(theMuonType_ == 0 || theMuonType_ == 4 || theMuonType_ == 5 || theMuonType_ >= 10 || theMuonType_==-1 || theMuonType_==-2 || theMuonType_==-3 || theMuonType_==-4) {
00510 MuScleFitUtils::MuonTypeForCheckMassWindow = 2;
00511 MuScleFitUtils::MuonType = 2;
00512 }
00513 else{
00514 std::cout<<"Wrong muon type "<<theMuonType_<<std::endl;
00515 exit(1);
00516 }
00517
00518
00519 if( theMuonType_ == 2 ) {
00520 MuScleFitUtils::rapidityBinsForZ_ = false;
00521 }
00522
00523
00524
00525 MuScleFitUtils::massWindowHalfWidth[0][0] = 20.;
00526 MuScleFitUtils::massWindowHalfWidth[0][1] = 0.35;
00527 MuScleFitUtils::massWindowHalfWidth[0][2] = 0.35;
00528 MuScleFitUtils::massWindowHalfWidth[0][3] = 0.35;
00529 MuScleFitUtils::massWindowHalfWidth[0][4] = 0.2;
00530 MuScleFitUtils::massWindowHalfWidth[0][5] = 0.2;
00531 MuScleFitUtils::massWindowHalfWidth[1][0] = 50.;
00532 MuScleFitUtils::massWindowHalfWidth[1][1] = 2.5;
00533 MuScleFitUtils::massWindowHalfWidth[1][2] = 2.5;
00534 MuScleFitUtils::massWindowHalfWidth[1][3] = 2.5;
00535 MuScleFitUtils::massWindowHalfWidth[1][4] = 1.5;
00536 MuScleFitUtils::massWindowHalfWidth[1][5] = 1.5;
00537 MuScleFitUtils::massWindowHalfWidth[2][0] = 20.;
00538 MuScleFitUtils::massWindowHalfWidth[2][1] = 0.35;
00539 MuScleFitUtils::massWindowHalfWidth[2][2] = 0.35;
00540 MuScleFitUtils::massWindowHalfWidth[2][3] = 0.35;
00541 MuScleFitUtils::massWindowHalfWidth[2][4] = 0.2;
00542 MuScleFitUtils::massWindowHalfWidth[2][5] = 0.2;
00543
00544 muonSelector_.reset(new MuScleFitMuonSelector(theMuonLabel_, theMuonType_, PATmuons_,
00545 MuScleFitUtils::resfind,
00546 MuScleFitUtils::speedup, genParticlesName_,
00547 compareToSimTracks_, simTracksCollection_,
00548 MuScleFitUtils::sherpa_, debug_));
00549
00550 MuScleFitUtils::backgroundHandler = new BackgroundHandler( pset.getParameter<std::vector<int> >("BgrFitType"),
00551 pset.getParameter<std::vector<double> >("LeftWindowBorder"),
00552 pset.getParameter<std::vector<double> >("RightWindowBorder"),
00553 MuScleFitUtils::ResMass,
00554 MuScleFitUtils::massWindowHalfWidth[MuScleFitUtils::MuonTypeForCheckMassWindow] );
00555
00556 MuScleFitUtils::crossSectionHandler = new CrossSectionHandler( MuScleFitUtils::parCrossSection, MuScleFitUtils::resfind );
00557
00558
00559
00560
00561 MuScleFitUtils::normalizeLikelihoodByEventNumber_ = pset.getUntrackedParameter<bool>("NormalizeLikelihoodByEventNumber", true);
00562 if(debug_>0) std::cout << "End of MuScleFit constructor" << std::endl;
00563
00564 inputRootTreeFileName_ = pset.getParameter<std::string>("InputRootTreeFileName");
00565 outputRootTreeFileName_ = pset.getParameter<std::string>("OutputRootTreeFileName");
00566 maxEventsFromRootTree_ = pset.getParameter<int>("MaxEventsFromRootTree");
00567
00568 MuScleFitUtils::startWithSimplex_ = pset.getParameter<bool>("StartWithSimplex");
00569 MuScleFitUtils::computeMinosErrors_ = pset.getParameter<bool>("ComputeMinosErrors");
00570 MuScleFitUtils::minimumShapePlots_ = pset.getParameter<bool>("MinimumShapePlots");
00571
00572 beginOfJobInConstructor();
00573 }
00574
00575
00576
00577 MuScleFit::~MuScleFit () {
00578 if (debug_>0) std::cout << "[MuScleFit]: Destructor" << std::endl;
00579 std::cout << "Total number of analyzed events = " << totalEvents_ << std::endl;
00580
00581 if( !(outputRootTreeFileName_.empty()) ) {
00582
00583 if( !(inputRootTreeFileName_.empty() && (int(MuScleFitUtils::SavedPair.size()) != totalEvents_)) ) {
00584 std::cout << "Saving muon pairs to root tree" << std::endl;
00585 RootTreeHandler rootTreeHandler;
00586 if( MuScleFitUtils::speedup ) {
00587
00588 rootTreeHandler.writeTree(outputRootTreeFileName_, &(muonPairs_), theMuonType_, 0, saveAllToTree_);
00589 }
00590 else {
00591
00592 rootTreeHandler.writeTree(outputRootTreeFileName_, &(muonPairs_), theMuonType_, &(genMuonPairs_), saveAllToTree_ );
00593 }
00594 }
00595 else {
00596 std::cout << "ERROR: events in the vector = " << MuScleFitUtils::SavedPair.size() << " != totalEvents = " << totalEvents_ << std::endl;
00597 }
00598 }
00599 }
00600
00601
00602
00603 void MuScleFit::beginOfJobInConstructor()
00604
00605
00606 {
00607 if (debug_>0) std::cout << "[MuScleFit]: beginOfJob" << std::endl;
00608
00609 if( MuScleFitUtils::useProbsFile_ ) {
00610 readProbabilityDistributionsFromFile();
00611 }
00612
00613 if (debug_>0) std::cout << "[MuScleFit]: beginOfJob" << std::endl;
00614
00615
00616
00617 for (unsigned int i=0; i<(maxLoopNumber); i++) {
00618 std::stringstream ss;
00619 ss << i;
00620 std::string rootFileName = ss.str() + "_" + theRootFileName_;
00621 theFiles_.push_back (new TFile(rootFileName.c_str(), "RECREATE"));
00622 }
00623 if (debug_>0) std::cout << "[MuScleFit]: Root file created" << std::endl;
00624
00625 std::cout << "creating plotter" << std::endl;
00626 plotter = new MuScleFitPlotter(theGenInfoRootFileName_);
00627 plotter->debug = debug_;
00628 }
00629
00630
00631
00632 void MuScleFit::endOfJob () {
00633 if (debug_>0) std::cout << "[MuScleFit]: endOfJob" << std::endl;
00634 }
00635
00636
00637
00638 void MuScleFit::startingNewLoop( unsigned int iLoop )
00639 {
00640 if (debug_>0) std::cout << "[MuScleFit]: Starting loop # " << iLoop << std::endl;
00641
00642
00643
00644 MuScleFitUtils::goodmuon = 0;
00645
00646
00647
00648 MuScleFitUtils::counter_resprob = 0;
00649
00650
00651
00652 fillHistoMap(theFiles_[iLoop], iLoop);
00653
00654 loopCounter = iLoop;
00655 MuScleFitUtils::loopCounter = loopCounter;
00656
00657 iev = 0;
00658 MuScleFitUtils::iev_ = 0;
00659
00660 MuScleFitUtils::oldNormalization_ = 0;
00661 }
00662
00663
00664
00665 edm::EDLooper::Status MuScleFit::endOfLoop( const edm::EventSetup& eventSetup, unsigned int iLoop )
00666 {
00667 unsigned int iFastLoop = 1;
00668
00669
00670 if( !(inputRootTreeFileName_.empty()) ) {
00671 selectMuons(maxEventsFromRootTree_, inputRootTreeFileName_);
00672
00673 totalEvents_ = MuScleFitUtils::SavedPair.size();
00674 iFastLoop = 0;
00675 }
00676 else {
00677 endOfFastLoop(iLoop);
00678 }
00679
00680
00681 if( fastLoop == true ) {
00682 for( ; iFastLoop<maxLoopNumber; ++iFastLoop ) {
00683
00684 std::cout << "Starting fast loop number " << iFastLoop << std::endl;
00685
00686
00687
00688 startingNewLoop(iFastLoop);
00689
00690
00691
00692
00693 while( iev<totalEvents_ ) {
00694 if( iev%50000 == 0 ) {
00695 std::cout << "Fast looping on event number " << iev << std::endl;
00696 }
00697
00698 duringFastLoop();
00699 }
00700 std::cout << "End of fast loop number " << iFastLoop << ". Ran on " << iev << " events" << std::endl;
00701 endOfFastLoop(iFastLoop);
00702 }
00703 }
00704
00705 if (iFastLoop>=maxLoopNumber-1) {
00706 return kStop;
00707 } else {
00708 return kContinue;
00709 }
00710 }
00711
00712 void MuScleFit::endOfFastLoop( const unsigned int iLoop )
00713 {
00714
00715
00716 if( loopCounter == 0 ) {
00717
00718
00719 delete plotter;
00720 }
00721
00722 std::cout << "Ending loop # " << iLoop << std::endl;
00723
00724
00725
00726
00727 writeHistoMap(iLoop);
00728
00729
00730
00731
00732 TDirectory * likelihoodDir = theFiles_[iLoop]->mkdir("likelihood");
00733 likelihoodDir->cd();
00734 MuScleFitUtils::minimizeLikelihood();
00735
00736
00737 theFiles_[iLoop]->Close();
00738
00739 delete theFiles_[iLoop];
00740
00741
00742
00743 clearHistoMap();
00744 }
00745
00746
00747
00748 edm::EDLooper::Status MuScleFit::duringLoop( const edm::Event & event, const edm::EventSetup& eventSetup )
00749 {
00750 edm::Handle<edm::TriggerResults> triggerResults;
00751 event.getByLabel(edm::InputTag(triggerResultsLabel_.c_str(), "", triggerResultsProcess_.c_str()), triggerResults);
00752
00753 bool isFired = false;
00754
00755 if(triggerPath_[0] == "")
00756 isFired = true;
00757 else if(triggerPath_[0] == "All"){
00758 isFired =triggerResults->accept();
00759 if(debug_>0)
00760 std::cout<<"Trigger "<<isFired<<std::endl;
00761 }
00762 else{
00763 bool changed;
00764 HLTConfigProvider hltConfig;
00765 hltConfig.init(event.getRun(), eventSetup, triggerResultsProcess_, changed);
00766
00767
00768 const edm::TriggerNames triggerNames = event.triggerNames(*triggerResults);
00769
00770 for (unsigned i=0; i<triggerNames.size(); i++) {
00771 std::string hltName = triggerNames.triggerName(i);
00772
00773
00774 for ( unsigned int ipath=0; ipath<triggerPath_.size(); ipath++ ) {
00775 if ( hltName.find(triggerPath_[ipath]) != std::string::npos ) {
00776 unsigned int triggerIndex( hltConfig.triggerIndex(hltName) );
00777
00778
00779 if (triggerIndex < triggerResults->size()) {
00780 isFired = triggerResults->accept(triggerIndex);
00781 if(debug_>0)
00782 std::cout << triggerPath_[ipath] <<" "<< hltName << " " << isFired<<std::endl;
00783 }
00784 }
00785 }
00786 }
00787
00788 }
00789
00790 if( negateTrigger_ && isFired ) return kContinue;
00791 else if( !(negateTrigger_) && !isFired ) return kContinue;
00792
00793 #ifdef USE_CALLGRIND
00794 CALLGRIND_START_INSTRUMENTATION;
00795 #endif
00796
00797 if (debug_>0) {
00798 std::cout << "[MuScleFit-duringLoop]: loopCounter = " << loopCounter
00799 << " Run: " << event.id().run() << " Event: " << event.id().event() << std::endl;
00800 }
00801
00802
00803
00804
00805
00806
00807 if( loopCounter == 0 ) {
00808
00809 if( !fastLoop || inputRootTreeFileName_.empty() ) {
00810 if( debug_ > 0 ) std::cout << "Reading from edm event" << std::endl;
00811 selectMuons(event);
00812 duringFastLoop();
00813 ++totalEvents_;
00814 }
00815 }
00816
00817 return kContinue;
00818
00819 #ifdef USE_CALLGRIND
00820 CALLGRIND_STOP_INSTRUMENTATION;
00821 CALLGRIND_DUMP_STATS;
00822 #endif
00823 }
00824
00825 void MuScleFit::selectMuons(const edm::Event & event)
00826 {
00827 recMu1 = reco::Particle::LorentzVector(0,0,0,0);
00828 recMu2 = reco::Particle::LorentzVector(0,0,0,0);
00829
00830 std::vector<reco::LeafCandidate> muons;
00831 muonSelector_->selectMuons(event, muons, genMuonPairs_, MuScleFitUtils::simPair, plotter);
00832
00833
00834
00835
00836 std::pair<reco::Particle::LorentzVector, reco::Particle::LorentzVector> recMuFromBestRes =
00837 MuScleFitUtils::findBestRecoRes(muons);
00838
00839 if (MuScleFitUtils::ResFound) {
00840 if (debug_>0) {
00841 std::cout <<std::setprecision(9)<< "Pt after findbestrecores: " << (recMuFromBestRes.first).Pt() << " "
00842 << (recMuFromBestRes.second).Pt() << std::endl;
00843 std::cout << "recMu1 = " << recMu1 << std::endl;
00844 std::cout << "recMu2 = " << recMu2 << std::endl;
00845 }
00846 recMu1 = recMuFromBestRes.first;
00847 recMu2 = recMuFromBestRes.second;
00848 if (debug_>0) {
00849 std::cout << "after recMu1 = " << recMu1 << std::endl;
00850 std::cout << "after recMu2 = " << recMu2 << std::endl;
00851 std::cout << "mu1.pt = " << recMu1.Pt() << std::endl;
00852 std::cout << "mu2.pt = " << recMu2.Pt() << std::endl;
00853 }
00854 MuScleFitUtils::SavedPair.push_back( std::make_pair( recMu1, recMu2 ) );
00855 } else {
00856 MuScleFitUtils::SavedPair.push_back( std::make_pair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.) ) );
00857 }
00858
00859
00860
00861 muonPairs_.push_back(MuonPair(MuScleFitUtils::SavedPair.back().first,
00862 MuScleFitUtils::SavedPair.back().second,
00863 event.run(), event.id().event()));
00864
00865 if( MuScleFitUtils::speedup == false ) {
00866 MuScleFitUtils::genPair.push_back(std::make_pair( genMuonPairs_.back().mu1, genMuonPairs_.back().mu2 ));
00867 }
00868 }
00869
00870 void MuScleFit::selectMuons(const int maxEvents, const TString & treeFileName)
00871 {
00872 std::cout << "Reading muon pairs from Root Tree in " << treeFileName << std::endl;
00873 RootTreeHandler rootTreeHandler;
00874 std::vector<std::pair<int, int> > evtRun;
00875 if( MuScleFitUtils::speedup ) {
00876 rootTreeHandler.readTree(maxEvents, inputRootTreeFileName_, &(MuScleFitUtils::SavedPair), theMuonType_, &evtRun);
00877 }
00878 else {
00879 rootTreeHandler.readTree(maxEvents, inputRootTreeFileName_, &(MuScleFitUtils::SavedPair), theMuonType_, &evtRun, &(MuScleFitUtils::genPair));
00880 }
00881
00882 std::vector<std::pair<int, int> >::iterator evtRunIt = evtRun.begin();
00883 std::vector<std::pair<lorentzVector,lorentzVector> >::iterator it = MuScleFitUtils::SavedPair.begin();
00884 std::vector<std::pair<lorentzVector,lorentzVector> >::iterator genIt;
00885 if(MuScleFitUtils::speedup == false) genIt = MuScleFitUtils::genPair.begin();
00886 for( ; it != MuScleFitUtils::SavedPair.end(); ++it, ++evtRunIt ) {
00887
00888
00889
00890
00891 double pt1 = it->first.pt();
00892
00893 double pt2 = it->second.pt();
00894
00895 double eta1 = it->first.eta();
00896
00897 double eta2 = it->second.eta();
00898
00899
00900 bool dontPass = false;
00901 bool eta1InFirstRange;
00902 bool eta2InFirstRange;
00903 bool eta1InSecondRange;
00904 bool eta2InSecondRange;
00905
00906 if( MuScleFitUtils::separateRanges_ ) {
00907 eta1InFirstRange = eta1 >= MuScleFitUtils::minMuonEtaFirstRange_ && eta1 < MuScleFitUtils::maxMuonEtaFirstRange_;
00908 eta2InFirstRange = eta2 >= MuScleFitUtils::minMuonEtaFirstRange_ && eta2 < MuScleFitUtils::maxMuonEtaFirstRange_;
00909 eta1InSecondRange = eta1 >= MuScleFitUtils::minMuonEtaSecondRange_ && eta1 < MuScleFitUtils::maxMuonEtaSecondRange_;
00910 eta2InSecondRange = eta2 >= MuScleFitUtils::minMuonEtaSecondRange_ && eta2 < MuScleFitUtils::maxMuonEtaSecondRange_;
00911
00912
00913 if( !(pt1 >= MuScleFitUtils::minMuonPt_ && pt1 < MuScleFitUtils::maxMuonPt_ &&
00914 pt2 >= MuScleFitUtils::minMuonPt_ && pt2 < MuScleFitUtils::maxMuonPt_ &&
00915 eta1InFirstRange && eta2InSecondRange ) ) {
00916 dontPass = true;
00917 }
00918 }
00919 else {
00920 eta1 = fabs(eta1);
00921 eta2 = fabs(eta2);
00922 eta1InFirstRange = eta1 >= MuScleFitUtils::minMuonEtaFirstRange_ && eta1 < MuScleFitUtils::maxMuonEtaFirstRange_;
00923 eta2InFirstRange = eta2 >= MuScleFitUtils::minMuonEtaFirstRange_ && eta2 < MuScleFitUtils::maxMuonEtaFirstRange_;
00924 eta1InSecondRange = eta1 >= MuScleFitUtils::minMuonEtaSecondRange_ && eta1 < MuScleFitUtils::maxMuonEtaSecondRange_;
00925 eta2InSecondRange = eta2 >= MuScleFitUtils::minMuonEtaSecondRange_ && eta2 < MuScleFitUtils::maxMuonEtaSecondRange_;
00926 if( !(pt1 >= MuScleFitUtils::minMuonPt_ && pt1 < MuScleFitUtils::maxMuonPt_ &&
00927 pt2 >= MuScleFitUtils::minMuonPt_ && pt2 < MuScleFitUtils::maxMuonPt_ &&
00928 ( ((eta1InFirstRange && !eta2InFirstRange) && (eta2InSecondRange && !eta1InSecondRange)) ||
00929 ((eta2InFirstRange && !eta1InFirstRange) && (eta1InSecondRange && !eta2InSecondRange)) )) ) {
00930 dontPass = true;
00931 }
00932 }
00933
00934
00935 double deltaPhi = MuScleFitUtils::deltaPhi(it->first.phi(), it->second.phi());
00936 if( (deltaPhi <= MuScleFitUtils::deltaPhiMinCut_) || (deltaPhi >= MuScleFitUtils::deltaPhiMaxCut_) ) dontPass = true;
00937
00938 if( dontPass ) {
00939
00940 it->first = reco::Particle::LorentzVector(0,0,0,0);
00941 it->second = reco::Particle::LorentzVector(0,0,0,0);
00942 }
00943
00944
00945 if( (MuScleFitUtils::SmearType != 0) || (MuScleFitUtils::BiasType != 0) ) {
00946 applySmearing(it->first);
00947 applyBias(it->first, -1);
00948 applySmearing(it->second);
00949 applyBias(it->second, 1);
00950 }
00951 muonPairs_.push_back(MuonPair(it->first, it->second,
00952 evtRunIt->second, evtRunIt->first));
00953
00954
00955 if( MuScleFitUtils::speedup == false ) {
00956 genMuonPairs_.push_back(GenMuonPair(genIt->first, genIt->second, 0));
00957 ++genIt;
00958 }
00959 }
00960 plotter->fillTreeRec(MuScleFitUtils::SavedPair);
00961 if( !(MuScleFitUtils::speedup) ) {
00962 plotter->fillTreeGen(MuScleFitUtils::genPair);
00963 }
00964 }
00965
00966 void MuScleFit::duringFastLoop()
00967 {
00968
00969
00970 MuScleFitUtils::ResFound = false;
00971 recMu1 = (MuScleFitUtils::SavedPair[iev].first);
00972 recMu2 = (MuScleFitUtils::SavedPair[iev].second);
00973 if (recMu1.Pt()>0 && recMu2.Pt()>0) {
00974 MuScleFitUtils::ResFound = true;
00975 if (debug_>0) std::cout << "Ev = " << iev << ": found muons in tree with Pt = "
00976 << recMu1.Pt() << " " << recMu2.Pt() << std::endl;
00977 }
00978
00979 if( debug_>0 ) std::cout << "About to start lik par correction and histo filling; ResFound is "
00980 << MuScleFitUtils::ResFound << std::endl;
00981
00982
00983 if( MuScleFitUtils::ResFound ) {
00984
00985
00986
00987
00988
00989 double weight = MuScleFitUtils::computeWeight( (recMu1+recMu2).mass(), iev, true );
00990 if (debug_>0) {
00991 std::cout << "Loop #" << loopCounter << "Event #" << iev << ": before correction Pt1 = "
00992 << recMu1.Pt() << " Pt2 = " << recMu2.Pt() << std::endl;
00993 }
00994
00995
00996 if ( loopCounter>0 ) {
00997 if ( MuScleFitUtils::doScaleFit[loopCounter-1] ) {
00998 recMu1 = (MuScleFitUtils::applyScale(recMu1, MuScleFitUtils::parvalue[loopCounter-1], -1));
00999 recMu2 = (MuScleFitUtils::applyScale(recMu2, MuScleFitUtils::parvalue[loopCounter-1], 1));
01000 }
01001 }
01002 if (debug_>0) {
01003 std::cout << "Loop #" << loopCounter << "Event #" << iev << ": after correction Pt1 = "
01004 << recMu1.Pt() << " Pt2 = " << recMu2.Pt() << std::endl;
01005 }
01006
01007 reco::Particle::LorentzVector bestRecRes( recMu1+recMu2 );
01008
01009
01010
01011
01012 mapHisto_["hRecBestMu"]->Fill(recMu1, -1,weight);
01013 mapHisto_["hRecBestMuVSEta"]->Fill(recMu1);
01014 mapHisto_["hRecBestMu"]->Fill(recMu2, +1,weight);
01015 mapHisto_["hRecBestMuVSEta"]->Fill(recMu2);
01016 mapHisto_["hDeltaRecBestMu"]->Fill(recMu1, recMu2);
01017
01018 mapHisto_["hRecBestRes"]->Fill(bestRecRes,+1, weight);
01019 mapHisto_["hRecBestResAllEvents"]->Fill(bestRecRes,+1, 1.);
01020
01021
01022
01023
01024
01025
01026 mapHisto_["hRecBestResVSMu"]->Fill (recMu1, bestRecRes, -1, weight);
01027 mapHisto_["hRecBestResVSMu"]->Fill (recMu2, bestRecRes, +1, weight);
01028
01029 mapHisto_["hRecBestResVSMu"]->Fill(recMu1, recMu2, bestRecRes, weight);
01030
01031
01032
01033
01034
01035
01036
01037
01038
01039
01040
01041 mapHisto_["hRecBestResVSRes"]->Fill (bestRecRes, bestRecRes, +1, weight);
01042
01043
01044
01045
01046
01047
01048 std::vector<double> * parval;
01049 std::vector<double> initpar;
01050
01051
01052 if (loopCounter==0) {
01053 initpar = MuScleFitUtils::parResol;
01054 initpar.insert( initpar.end(), MuScleFitUtils::parScale.begin(), MuScleFitUtils::parScale.end() );
01055 initpar.insert( initpar.end(), MuScleFitUtils::parCrossSection.begin(), MuScleFitUtils::parCrossSection.end() );
01056 initpar.insert( initpar.end(), MuScleFitUtils::parBgr.begin(), MuScleFitUtils::parBgr.end() );
01057 parval = &initpar;
01058 } else {
01059 parval = &(MuScleFitUtils::parvalue[loopCounter-1]);
01060 }
01061
01062
01063
01064 if( !MuScleFitUtils::speedup ) {
01065
01066
01067 if(checkDeltaR(MuScleFitUtils::genPair[iev].first,recMu1)) {
01068 fillComparisonHistograms( MuScleFitUtils::genPair[iev].first, recMu1, "Gen", -1 );
01069 }
01070 if(checkDeltaR(MuScleFitUtils::genPair[iev].second,recMu2)){
01071 fillComparisonHistograms( MuScleFitUtils::genPair[iev].second, recMu2, "Gen", +1 );
01072 }
01073 if( compareToSimTracks_ ) {
01074
01075 if(checkDeltaR(MuScleFitUtils::simPair[iev].first,recMu1)){
01076 fillComparisonHistograms( MuScleFitUtils::simPair[iev].first, recMu1, "Sim", -1 );
01077 }
01078 if(checkDeltaR(MuScleFitUtils::simPair[iev].second,recMu2)){
01079 fillComparisonHistograms( MuScleFitUtils::simPair[iev].second, recMu2, "Sim", +1 );
01080 }
01081 }
01082 }
01083
01084
01085
01086
01087
01088 mapHisto_["hFunctionResolPt"]->Fill( recMu1, MuScleFitUtils::resolutionFunctionForVec->sigmaPt(recMu1.Pt(), recMu1.Eta(), *parval ), -1 );
01089 mapHisto_["hFunctionResolCotgTheta"]->Fill( recMu1, MuScleFitUtils::resolutionFunctionForVec->sigmaCotgTh(recMu1.Pt(), recMu1.Eta(), *parval ), -1 );
01090 mapHisto_["hFunctionResolPhi"]->Fill( recMu1, MuScleFitUtils::resolutionFunctionForVec->sigmaPhi(recMu1.Pt(), recMu1.Eta(), *parval ), -1 );
01091 mapHisto_["hFunctionResolPt"]->Fill( recMu2, MuScleFitUtils::resolutionFunctionForVec->sigmaPt(recMu2.Pt(), recMu2.Eta(), *parval ), +1 );
01092 mapHisto_["hFunctionResolCotgTheta"]->Fill( recMu2, MuScleFitUtils::resolutionFunctionForVec->sigmaCotgTh(recMu2.Pt(), recMu2.Eta(), *parval ), +1 );
01093 mapHisto_["hFunctionResolPhi"]->Fill( recMu2, MuScleFitUtils::resolutionFunctionForVec->sigmaPhi(recMu2.Pt(), recMu2.Eta(), *parval ), +1 );
01094
01095
01096
01097 if( debug_ > 0 ) std::cout << "mass = " << bestRecRes.mass() << std::endl;
01098 if (weight!=0.) {
01099 double massResol;
01100 double prob;
01101 double deltalike;
01102 if (loopCounter==0) {
01103 std::vector<double> initpar;
01104 for (int i=0; i<(int)(MuScleFitUtils::parResol.size()); i++) {
01105 initpar.push_back(MuScleFitUtils::parResol[i]);
01106 }
01107 for (int i=0; i<(int)(MuScleFitUtils::parScale.size()); i++) {
01108 initpar.push_back(MuScleFitUtils::parScale[i]);
01109 }
01110
01111
01112
01113 MuScleFitUtils::crossSectionHandler->addParameters(initpar);
01114
01115 for (int i=0; i<(int)(MuScleFitUtils::parBgr.size()); i++) {
01116 initpar.push_back(MuScleFitUtils::parBgr[i]);
01117 }
01118 massResol = MuScleFitUtils::massResolution( recMu1, recMu2, initpar );
01119
01120 prob = MuScleFitUtils::massProb( bestRecRes.mass(), bestRecRes.Eta(), bestRecRes.Rapidity(), massResol,
01121 initpar, true, recMu1.eta(), recMu2.eta() );
01122 } else {
01123 massResol = MuScleFitUtils::massResolution( recMu1, recMu2,
01124 MuScleFitUtils::parvalue[loopCounter-1] );
01125
01126
01127 prob = MuScleFitUtils::massProb( bestRecRes.mass(), bestRecRes.Eta(), bestRecRes.Rapidity(),
01128 massResol, MuScleFitUtils::parvalue[loopCounter-1], true,
01129 recMu1.eta(), recMu2.eta() );
01130 }
01131 if( debug_ > 0 ) std::cout << "inside weight: mass = " << bestRecRes.mass() << ", prob = " << prob << std::endl;
01132 if (prob>0) {
01133 if( debug_ > 0 ) std::cout << "inside prob: mass = " << bestRecRes.mass() << ", prob = " << prob << std::endl;
01134
01135 deltalike = log(prob)*weight;
01136 mapHisto_["hLikeVSMu"]->Fill(recMu1, deltalike);
01137 mapHisto_["hLikeVSMu"]->Fill(recMu2, deltalike);
01138 mapHisto_["hLikeVSMuMinus"]->Fill(recMu1, deltalike);
01139 mapHisto_["hLikeVSMuPlus"]->Fill(recMu2, deltalike);
01140
01141 double recoMass = (recMu1+recMu2).mass();
01142 if( recoMass != 0 ) {
01143
01144 mapHisto_["hResolMassVSMu"]->Fill(recMu1, massResol, -1);
01145 mapHisto_["hResolMassVSMu"]->Fill(recMu2, massResol, +1);
01146 mapHisto_["hFunctionResolMassVSMu"]->Fill(recMu1, massResol/recoMass, -1);
01147 mapHisto_["hFunctionResolMassVSMu"]->Fill(recMu2, massResol/recoMass, +1);
01148 }
01149
01150 if( MuScleFitUtils::debugMassResol_ ) {
01151 mapHisto_["hdMdPt1"]->Fill(recMu1, MuScleFitUtils::massResolComponents.dmdpt1, -1);
01152 mapHisto_["hdMdPt2"]->Fill(recMu2, MuScleFitUtils::massResolComponents.dmdpt2, +1);
01153 mapHisto_["hdMdPhi1"]->Fill(recMu1, MuScleFitUtils::massResolComponents.dmdphi1, -1);
01154 mapHisto_["hdMdPhi2"]->Fill(recMu2, MuScleFitUtils::massResolComponents.dmdphi2, +1);
01155 mapHisto_["hdMdCotgTh1"]->Fill(recMu1, MuScleFitUtils::massResolComponents.dmdcotgth1, -1);
01156 mapHisto_["hdMdCotgTh2"]->Fill(recMu2, MuScleFitUtils::massResolComponents.dmdcotgth2, +1);
01157 }
01158
01159 if( !MuScleFitUtils::speedup ) {
01160 double genMass = (MuScleFitUtils::genPair[iev].first + MuScleFitUtils::genPair[iev].second).mass();
01161
01162 if( genMass != 0 ) {
01163 mapHisto_["hGenResVSMu"]->Fill((MuScleFitUtils::genPair[iev].first), (MuScleFitUtils::genPair[iev].first + MuScleFitUtils::genPair[iev].second), -1);
01164 mapHisto_["hGenResVSMu"]->Fill((MuScleFitUtils::genPair[iev].second), (MuScleFitUtils::genPair[iev].first + MuScleFitUtils::genPair[iev].second), +1);
01165 double diffMass = (recoMass - genMass)/genMass;
01166
01167
01168 double pt1 = recMu1.pt();
01169 double eta1 = recMu1.eta();
01170 double pt2 = recMu2.pt();
01171 double eta2 = recMu2.eta();
01172
01173 if( diffMass == diffMass ) {
01174
01175 mapHisto_["hDeltaMassOverGenMassVsPt"]->Fill(pt1, diffMass);
01176 mapHisto_["hDeltaMassOverGenMassVsPt"]->Fill(pt2, diffMass);
01177 mapHisto_["hDeltaMassOverGenMassVsEta"]->Fill(eta1, diffMass);
01178 mapHisto_["hDeltaMassOverGenMassVsEta"]->Fill(eta2, diffMass);
01179
01180 mapHisto_["hMassResolutionVsPtEta"]->Fill(pt1, eta1, diffMass, diffMass);
01181 mapHisto_["hMassResolutionVsPtEta"]->Fill(pt2, eta2, diffMass, diffMass);
01182 }
01183 else {
01184 std::cout << "Error, there is a nan: recoMass = " << recoMass << ", genMass = " << genMass << std::endl;
01185 }
01186 }
01187
01188 double massRes = MuScleFitUtils::massResolution(recMu1, recMu2, MuScleFitUtils::parResol);
01189 mapHisto_["hFunctionResolMass"]->Fill( recMu1, std::pow(massRes,2), -1 );
01190 mapHisto_["hFunctionResolMass"]->Fill( recMu2, std::pow(massRes,2), +1 );
01191 }
01192
01193 mapHisto_["hMass_P"]->Fill(bestRecRes.mass(), prob);
01194 if( debug_ > 0 ) std::cout << "mass = " << bestRecRes.mass() << ", prob = " << prob << std::endl;
01195 mapHisto_["hMass_fine_P"]->Fill(bestRecRes.mass(), prob);
01196
01197 mapHisto_["hMassProbVsRes"]->Fill(bestRecRes, bestRecRes, +1, prob);
01198 mapHisto_["hMassProbVsMu"]->Fill(recMu1, bestRecRes, -1, prob);
01199 mapHisto_["hMassProbVsMu"]->Fill(recMu2, bestRecRes, +1, prob);
01200 mapHisto_["hMassProbVsRes_fine"]->Fill(bestRecRes, bestRecRes, +1, prob);
01201 mapHisto_["hMassProbVsMu_fine"]->Fill(recMu1, bestRecRes, -1, prob);
01202 mapHisto_["hMassProbVsMu_fine"]->Fill(recMu2, bestRecRes, +1, prob);
01203 }
01204 }
01205 }
01206
01207
01208
01209 if (loopCounter>0) {
01210 if (debug_>0) std::cout << "[MuScleFit]: filling the pair" << std::endl;
01211 MuScleFitUtils::SavedPair[iev] = std::make_pair( recMu1, recMu2 );
01212 }
01213
01214 iev++;
01215 MuScleFitUtils::iev_++;
01216
01217
01218 }
01219
01220 bool MuScleFit::checkDeltaR(reco::Particle::LorentzVector& genMu, reco::Particle::LorentzVector& recMu){
01221
01222 double deltaR = sqrt(MuScleFitUtils::deltaPhi(recMu.Phi(),genMu.Phi()) * MuScleFitUtils::deltaPhi(recMu.Phi(),genMu.Phi()) +
01223 ((recMu.Eta()-genMu.Eta()) * (recMu.Eta()-genMu.Eta())));
01224 if(deltaR<0.01)
01225 return true;
01226 else if( debug_ > 0 ) {
01227 std::cout<<"Reco muon "<<recMu<<" with eta "<<recMu.Eta()<<" and phi "<<recMu.Phi()<<std::endl
01228 <<" DOES NOT MATCH with generated muon from resonance: "<<std::endl
01229 <<genMu<<" with eta "<<genMu.Eta()<<" and phi "<<genMu.Phi()<<std::endl;
01230 }
01231 return false;
01232 }
01233
01234 void MuScleFit::fillComparisonHistograms( const reco::Particle::LorentzVector & genMu, const reco::Particle::LorentzVector & recMu,
01235 const std::string & inputName, const int charge )
01236 {
01237 std::string name(inputName + "VSMu");
01238 mapHisto_["hResolPt"+name]->Fill(recMu, (-genMu.Pt()+recMu.Pt())/genMu.Pt(), charge);
01239 mapHisto_["hResolTheta"+name]->Fill(recMu, (-genMu.Theta()+recMu.Theta()), charge);
01240 mapHisto_["hResolCotgTheta"+name]->Fill(recMu,(-cos(genMu.Theta())/sin(genMu.Theta())
01241 +cos(recMu.Theta())/sin(recMu.Theta())), charge);
01242 mapHisto_["hResolEta"+name]->Fill(recMu, (-genMu.Eta()+recMu.Eta()),charge);
01243 mapHisto_["hResolPhi"+name]->Fill(recMu, MuScleFitUtils::deltaPhiNoFabs(recMu.Phi(), genMu.Phi()), charge);
01244
01245
01246 if( (genMu.Pt() != 0) && (recMu.Pt() != 0) ) {
01247 mapHisto_["hPtRecoVsPt"+inputName]->Fill(genMu.Pt(), recMu.Pt());
01248 }
01249 }
01250
01251 void MuScleFit::applySmearing( reco::Particle::LorentzVector & mu )
01252 {
01253 if( MuScleFitUtils::SmearType>0 ) {
01254 mu = MuScleFitUtils::applySmearing( mu );
01255 if (debug_>0) std::cout << "Muon #" << MuScleFitUtils::goodmuon
01256 << ": after smearing Pt = " << mu.Pt() << std::endl;
01257 }
01258 }
01259
01260 void MuScleFit::applyBias( reco::Particle::LorentzVector & mu, const int charge )
01261 {
01262 if( MuScleFitUtils::BiasType>0 ) {
01263 mu = MuScleFitUtils::applyBias( mu, charge );
01264 if (debug_>0) std::cout << "Muon #" << MuScleFitUtils::goodmuon
01265 << ": after bias Pt = " << mu.Pt() << std::endl;
01266 }
01267 }
01268
01269
01270
01271 void MuScleFit::checkParameters() {
01272
01273
01274 if( MuScleFitUtils::doResolFit.size() != maxLoopNumber ) {
01275 std::cout << "[MuScleFit-Constructor]: wrong size of resolution fits selector = " << MuScleFitUtils::doResolFit.size() << std::endl;
01276 std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
01277 abort();
01278 }
01279 if( MuScleFitUtils::doScaleFit.size() != maxLoopNumber ) {
01280 std::cout << "[MuScleFit-Constructor]: wrong size of scale fits selector = " << MuScleFitUtils::doScaleFit.size() << std::endl;
01281 std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
01282 abort();
01283 }
01284 if( MuScleFitUtils::doCrossSectionFit.size() != maxLoopNumber ) {
01285 std::cout << "[MuScleFit-Constructor]: wrong size of cross section fits selector = " << MuScleFitUtils::doCrossSectionFit.size() << std::endl;
01286 std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
01287 abort();
01288 }
01289 if( MuScleFitUtils::doBackgroundFit.size() != maxLoopNumber ) {
01290 std::cout << "[MuScleFit-Constructor]: wrong size of background fits selector = " << MuScleFitUtils::doBackgroundFit.size() << std::endl;
01291 std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
01292 abort();
01293 }
01294
01295
01296
01297 if ((MuScleFitUtils::BiasType==1 && MuScleFitUtils::parBias.size()!=2) ||
01298 (MuScleFitUtils::BiasType==2 && MuScleFitUtils::parBias.size()!=2) ||
01299 (MuScleFitUtils::BiasType==3 && MuScleFitUtils::parBias.size()!=4) ||
01300 (MuScleFitUtils::BiasType==4 && MuScleFitUtils::parBias.size()!=3) ||
01301 (MuScleFitUtils::BiasType==5 && MuScleFitUtils::parBias.size()!=3) ||
01302 (MuScleFitUtils::BiasType==6 && MuScleFitUtils::parBias.size()!=3) ||
01303 (MuScleFitUtils::BiasType==7 && MuScleFitUtils::parBias.size()!=4) ||
01304 (MuScleFitUtils::BiasType==8 && MuScleFitUtils::parBias.size()!=4) ||
01305 (MuScleFitUtils::BiasType==9 && MuScleFitUtils::parBias.size()!=2) ||
01306 (MuScleFitUtils::BiasType==10 && MuScleFitUtils::parBias.size()!=3) ||
01307 (MuScleFitUtils::BiasType==11 && MuScleFitUtils::parBias.size()!=4) ||
01308 (MuScleFitUtils::BiasType==12 && MuScleFitUtils::parBias.size()!=6) ||
01309 (MuScleFitUtils::BiasType==13 && MuScleFitUtils::parBias.size()!=8) ||
01310 MuScleFitUtils::BiasType<0 || MuScleFitUtils::BiasType>13) {
01311 std::cout << "[MuScleFit-Constructor]: Wrong bias type or number of parameters: aborting!" << std::endl;
01312 abort();
01313 }
01314
01315
01316 if ((MuScleFitUtils::SmearType==1 && MuScleFitUtils::parSmear.size()!=3) ||
01317 (MuScleFitUtils::SmearType==2 && MuScleFitUtils::parSmear.size()!=4) ||
01318 (MuScleFitUtils::SmearType==3 && MuScleFitUtils::parSmear.size()!=5) ||
01319 (MuScleFitUtils::SmearType==4 && MuScleFitUtils::parSmear.size()!=6) ||
01320 (MuScleFitUtils::SmearType==5 && MuScleFitUtils::parSmear.size()!=7) ||
01321 (MuScleFitUtils::SmearType==6 && MuScleFitUtils::parSmear.size()!=16) ||
01322 (MuScleFitUtils::SmearType==7 && MuScleFitUtils::parSmear.size()!=0) ||
01323 MuScleFitUtils::SmearType<0 || MuScleFitUtils::SmearType>7) {
01324 std::cout << "[MuScleFit-Constructor]: Wrong smear type or number of parameters: aborting!" << std::endl;
01325 abort();
01326 }
01327
01328
01329 if (MuScleFitUtils::parResol.size()!=MuScleFitUtils::parResolFix.size() ||
01330 MuScleFitUtils::parResol.size()!=MuScleFitUtils::parResolOrder.size()) {
01331 std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Resol: aborting!" << std::endl;
01332 abort();
01333 }
01334 if (MuScleFitUtils::parScale.size()!=MuScleFitUtils::parScaleFix.size() ||
01335 MuScleFitUtils::parScale.size()!=MuScleFitUtils::parScaleOrder.size()) {
01336 std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Scale: aborting!" << std::endl;
01337 abort();
01338 }
01339 if (MuScleFitUtils::parCrossSection.size()!=MuScleFitUtils::parCrossSectionFix.size() ||
01340 MuScleFitUtils::parCrossSection.size()!=MuScleFitUtils::parCrossSectionOrder.size()) {
01341 std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Bgr: aborting!" << std::endl;
01342 abort();
01343 }
01344 if (MuScleFitUtils::parBgr.size()!=MuScleFitUtils::parBgrFix.size() ||
01345 MuScleFitUtils::parBgr.size()!=MuScleFitUtils::parBgrOrder.size()) {
01346 std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Bgr: aborting!" << std::endl;
01347 abort();
01348 }
01349
01350
01351
01352 if (MuScleFitUtils::resfind.size()!=6) {
01353 std::cout << "[MuScleFit-Constructor]: resfind must have 6 elements (1 Z, 3 Y, 2 Psi): aborting!" << std::endl;
01354 abort();
01355 }
01356 }
01357
01358 bool MuScleFit::selGlobalMuon(const pat::Muon* aMuon) {
01359
01360 reco::TrackRef iTrack = aMuon->innerTrack();
01361 const reco::HitPattern& p = iTrack->hitPattern();
01362
01363 reco::TrackRef gTrack = aMuon->globalTrack();
01364 const reco::HitPattern& q = gTrack->hitPattern();
01365
01366 return (
01367 iTrack->found() > 11 &&
01368 gTrack->chi2()/gTrack->ndof() < 20.0 &&
01369 q.numberOfValidMuonHits() > 0 &&
01370 iTrack->chi2()/iTrack->ndof() < 4.0 &&
01371 aMuon->muonID("TrackerMuonArbitrated") &&
01372 aMuon->muonID("TMLastStationAngTight") &&
01373 p.pixelLayersWithMeasurement() > 1 &&
01374 fabs(iTrack->dxy()) < 3.0 &&
01375 fabs(iTrack->dz()) < 15.0 );
01376 }
01377
01378
01379 bool MuScleFit::selTrackerMuon(const pat::Muon* aMuon) {
01380
01381 reco::TrackRef iTrack = aMuon->innerTrack();
01382 const reco::HitPattern& p = iTrack->hitPattern();
01383
01384 return (
01385 iTrack->found() > 11 &&
01386 iTrack->chi2()/iTrack->ndof() < 4.0 &&
01387 aMuon->muonID("TrackerMuonArbitrated") &&
01388 aMuon->muonID("TMLastStationAngTight") &&
01389 p.pixelLayersWithMeasurement() > 1 &&
01390 fabs(iTrack->dxy()) < 3.0 &&
01391 fabs(iTrack->dz()) < 15.0 );
01392
01393 }
01394
01395
01396 #include "FWCore/Framework/interface/MakerMacros.h"
01397 DEFINE_FWK_LOOPER(MuScleFit);