CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/DQM/BeamMonitor/plugins/AlcaBeamMonitor.cc

Go to the documentation of this file.
00001 /*
00002  * \file AlcaBeamMonitor.cc
00003  * \author Lorenzo Uplegger/FNAL
00004  * $Date: 2011/03/22 03:07:39 $
00005  * $Revision: 1.9 $
00006  *
00007  */
00008 
00009 #include "FWCore/Framework/interface/MakerMacros.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/Framework/interface/ESHandle.h"
00012 #include "FWCore/ServiceRegistry/interface/Service.h"
00013 #include "CondFormats/DataRecord/interface/BeamSpotObjectsRcd.h"
00014 #include "CondFormats/BeamSpotObjects/interface/BeamSpotObjects.h"
00015 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00016 //#include "DataFormats/Scalers/interface/BeamSpotOnline.h"
00017 #include "DataFormats/Common/interface/View.h"
00018 #include "DataFormats/Common/interface/Handle.h"
00019 #include "RecoVertex/BeamSpotProducer/interface/BeamFitter.h"
00020 #include "RecoVertex/BeamSpotProducer/interface/PVFitter.h"
00021 #include "DQM/BeamMonitor/plugins/AlcaBeamMonitor.h"
00022 #include "DQMServices/Core/interface/QReport.h"
00023 #include "FWCore/Framework/interface/Run.h"
00024 #include "FWCore/Framework/interface/LuminosityBlock.h"
00025 #include "FWCore/Framework/interface/EventSetup.h"
00026 #include <numeric>
00027 //#include <iostream>
00028 
00029 using namespace std;
00030 using namespace edm;
00031 using namespace reco;
00032 
00033 //----------------------------------------------------------------------------------------------------------------------
00034 AlcaBeamMonitor::AlcaBeamMonitor( const ParameterSet& ps ) : 
00035   parameters_           (ps),
00036   monitorName_          (parameters_.getUntrackedParameter<string>("MonitorName","YourSubsystemName")),
00037   primaryVertexLabel_   (parameters_.getUntrackedParameter<InputTag>("PrimaryVertexLabel")),
00038   beamSpotLabel_        (parameters_.getUntrackedParameter<InputTag>("BeamSpotLabel")),
00039   trackLabel_           (parameters_.getUntrackedParameter<InputTag>("TrackLabel")),
00040   scalerLabel_          (parameters_.getUntrackedParameter<InputTag>("ScalerLabel")),
00041   numberOfValuesToSave_ (0)
00042 {
00043   dbe_ = Service<DQMStore>().operator->();
00044   
00045   if (monitorName_ != "" ) monitorName_ = monitorName_+"/" ;
00046   
00047   theBeamFitter_ = new BeamFitter(parameters_);
00048   theBeamFitter_->resetTrkVector();
00049   theBeamFitter_->resetLSRange();
00050   theBeamFitter_->resetRefTime();
00051   theBeamFitter_->resetPVFitter();
00052 
00053   thePVFitter_ = new PVFitter(parameters_);
00054 
00055 
00056   varNamesV_.push_back("x");
00057   varNamesV_.push_back("y");
00058   varNamesV_.push_back("z");
00059   varNamesV_.push_back("sigmaX");
00060   varNamesV_.push_back("sigmaY");
00061   varNamesV_.push_back("sigmaZ");
00062 
00063   histoByCategoryNames_.insert( pair<string,string>("run",        "Coordinate"));
00064   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex fit-DataBase"));
00065   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex fit-BeamFit"));
00066   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex fit-Scalers"));
00067   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex-DataBase"));
00068   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex-BeamFit"));
00069   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex-Scalers"));
00070 
00071   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased BeamSpotFit"));  
00072   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased PrimaryVertex"));
00073   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased DataBase"));     
00074   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased Scalers"));      
00075   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased PrimaryVertex-DataBase fit"));
00076   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased PrimaryVertex-Scalers fit"));
00077   histoByCategoryNames_.insert( pair<string,string>("validation", "Lumibased Scalers-DataBase fit"));
00078   histoByCategoryNames_.insert( pair<string,string>("validation", "Lumibased PrimaryVertex-DataBase"));
00079   histoByCategoryNames_.insert( pair<string,string>("validation", "Lumibased PrimaryVertex-Scalers"));
00080 
00081 
00082   for(vector<string>::iterator itV=varNamesV_.begin(); itV!=varNamesV_.end(); itV++){
00083     for(multimap<string,string>::iterator itM=histoByCategoryNames_.begin(); itM!=histoByCategoryNames_.end(); itM++){
00084       if(itM->first=="run"){
00085         histosMap_[*itV][itM->first][itM->second] = 0;
00086       }
00087       else{
00088         positionsMap_[*itV][itM->first][itM->second] = 3*numberOfValuesToSave_;//value, error, ok 
00089         ++numberOfValuesToSave_;
00090       }
00091     }
00092   }
00093   
00094 //  beamSpotsMap_["BF"] = map<LuminosityBlockNumber_t,BeamSpot>();//For each lumi the beamfitter will have a result
00095 //  beamSpotsMap_["PV"] = map<LuminosityBlockNumber_t,BeamSpot>();//For each lumi the PVfitter will have a result
00096 //  beamSpotsMap_["DB"] = map<LuminosityBlockNumber_t,BeamSpot>();//For each lumi we take the values that are stored in the database, already collapsed then
00097 //  beamSpotsMap_["SC"] = map<LuminosityBlockNumber_t,BeamSpot>();//For each lumi we take the beamspot value in the file that is the same as the scaler for the alca reco stream
00098 //  beamSpotsMap_["BF"] = 0;//For each lumi the beamfitter will have a result
00099 //  beamSpotsMap_["PV"] = 0;//For each lumi the PVfitter will have a result
00100 //  beamSpotsMap_["DB"] = 0;//For each lumi we take the values that are stored in the database, already collapsed then
00101 //  beamSpotsMap_["SC"] = 0;//For each lumi we take the beamspot value in the file that is the same as the scaler for the alca reco stream
00102 }
00103 
00104 
00105 AlcaBeamMonitor::~AlcaBeamMonitor() {
00106   if(theBeamFitter_ != 0){
00107     delete theBeamFitter_;
00108   }
00109   
00110   if(thePVFitter_ != 0){
00111     delete thePVFitter_;
00112   }
00113 }
00114 
00115 
00116 //----------------------------------------------------------------------------------------------------------------------
00117 void AlcaBeamMonitor::beginJob() {
00118   string name;
00119   string title;
00120   dbe_->setCurrentFolder(monitorName_+"Debug");
00121   for(HistosContainer::iterator itM=histosMap_.begin(); itM!=histosMap_.end(); itM++){
00122     for(map<string,MonitorElement*>::iterator itMM=itM->second["run"].begin(); itMM!=itM->second["run"].end(); itMM++){
00123       name = string("h") + itM->first + itMM->first;
00124       title = itM->first + "_{0} " + itMM->first;
00125       if(itM->first == "x" || itM->first == "y"){
00126         if(itMM->first == "Coordinate"){
00127           itMM->second = dbe_->book1D(name,title,1001,-0.2525,0.2525);
00128         }
00129         else if(itMM->first == "PrimaryVertex fit-DataBase" || itMM->first == "PrimaryVertex fit-BeamFit" || itMM->first == "PrimaryVertex fit-Scalers"
00130              || itMM->first == "PrimaryVertex-DataBase" || itMM->first == "PrimaryVertex-BeamFit" || itMM->first == "PrimaryVertex-Scalers"){
00131           itMM->second = dbe_->book1D(name,title,1001,-0.02525,0.02525);
00132         }
00133         else{
00134           //assert(0);
00135         }
00136       }
00137       else if(itM->first == "z"){
00138         if(itMM->first == "Coordinate"){
00139           itMM->second = dbe_->book1D(name,title,101,-5.05,5.05);
00140         }
00141         else if(itMM->first == "PrimaryVertex fit-DataBase" || itMM->first == "PrimaryVertex fit-BeamFit" || itMM->first == "PrimaryVertex fit-Scalers"){
00142           itMM->second = dbe_->book1D(name,title,101,-0.505,0.505);
00143         }
00144         else if(itMM->first == "PrimaryVertex-DataBase" || itMM->first == "PrimaryVertex-BeamFit" || itMM->first == "PrimaryVertex-Scalers"){
00145           itMM->second = dbe_->book1D(name,title,1001,-5.005,5.005);
00146         }
00147         else{
00148           //assert(0);
00149         }
00150       }
00151       else if(itM->first == "sigmaX" || itM->first == "sigmaY"){
00152         if(itMM->first == "Coordinate"){
00153           itMM->second = dbe_->book1D(name,title,100,0,0.015);
00154         }
00155         else if(itMM->first == "PrimaryVertex fit-DataBase" || itMM->first == "PrimaryVertex fit-BeamFit" || itMM->first == "PrimaryVertex fit-Scalers"
00156              || itMM->first == "PrimaryVertex-DataBase" || itMM->first == "PrimaryVertex-BeamFit" || itMM->first == "PrimaryVertex-Scalers"){
00157           itMM->second = 0;
00158         }
00159         else{
00160           //assert(0);
00161         }
00162       }
00163       else if(itM->first == "sigmaZ"){
00164         if(itMM->first == "Coordinate"){
00165           itMM->second = dbe_->book1D(name,title,110,0,11);
00166         }
00167         else if(itMM->first == "PrimaryVertex fit-DataBase" || itMM->first == "PrimaryVertex fit-BeamFit" || itMM->first == "PrimaryVertex fit-Scalers"
00168              || itMM->first == "PrimaryVertex-DataBase" || itMM->first == "PrimaryVertex-BeamFit" || itMM->first == "PrimaryVertex-Scalers"){
00169           itMM->second = dbe_->book1D(name,title,101,-5.05,5.05);
00170         }
00171         else{
00172           //assert(0);
00173         }
00174       }
00175       else{
00176         //assert(0);
00177       }
00178       if(itMM->second != 0){
00179         if(itMM->first == "Coordinate"){                                
00180           itMM->second->setAxisTitle(itM->first + "_{0} (cm)",1);  
00181         }
00182         else if(itMM->first == "PrimaryVertex fit-DataBase" || itMM->first == "PrimaryVertex fit-BeamFit" || itMM->first == "PrimaryVertex fit-Scalers"
00183              || itMM->first == "PrimaryVertex-DataBase" || itMM->first == "PrimaryVertex-BeamFit" || itMM->first == "PrimaryVertex-Scalers"){
00184           itMM->second->setAxisTitle(itMM->first + " " + itM->first + "_{0} (cm)",1);  
00185         }
00186         itMM->second->setAxisTitle("Entries",2);
00187       }         
00188     }
00189   }
00190   dbe_->setCurrentFolder(monitorName_+"Service");
00191   theValuesContainer_ = dbe_->bookProfile("hHistoLumiValues","Histo Lumi Values", 3*numberOfValuesToSave_, 0., 3*numberOfValuesToSave_, 100., -100., 9000., " ");
00192   theValuesContainer_->setLumiFlag();
00193 
00194 }
00195 
00196 //----------------------------------------------------------------------------------------------------------------------
00197 void AlcaBeamMonitor::beginRun(const edm::Run& r, const EventSetup& context) {
00198   // create and cd into new folder
00199   dbe_->setCurrentFolder(monitorName_+"Validation");
00200   //Book histograms
00201   hD0Phi0_ = dbe_->bookProfile("hD0Phi0","d_{0} vs. #phi_{0} (All Tracks)",63,-3.15,3.15,100,-0.1,0.1,"");
00202   hD0Phi0_->setAxisTitle("#phi_{0} (rad)",1);
00203   hD0Phi0_->setAxisTitle("d_{0} (cm)",2);
00204 
00205   dbe_->setCurrentFolder(monitorName_+"Debug");
00206   hDxyBS_ = dbe_->book1D("hDxyBS","dxy_{0} w.r.t. Beam spot (All Tracks)",100,-0.1,0.1);
00207   hDxyBS_->setAxisTitle("dxy_{0} w.r.t. Beam spot (cm)",1);
00208 }
00209 
00210 //----------------------------------------------------------------------------------------------------------------------
00211 void AlcaBeamMonitor::beginLuminosityBlock(const LuminosityBlock& iLumi, const EventSetup& iSetup) {
00212   // Always create a beamspot group for each lumi weather we have results or not! Each Beamspot will be of unknown type!
00213   
00214   vertices_.clear();
00215   theValuesContainer_->Reset();
00216   beamSpotsMap_.clear();
00217   
00218   //Read BeamSpot from DB
00219   ESHandle<BeamSpotObjects> bsDBHandle;
00220   try{
00221     iSetup.get<BeamSpotObjectsRcd>().get(bsDBHandle);
00222   }
00223   catch( cms::Exception& exception ){                                 
00224     LogInfo("AlcaBeamMonitor") 
00225       << exception.what(); 
00226     return;           
00227   }                                   
00228   if(bsDBHandle.isValid()) { // check the product
00229     const BeamSpotObjects *spotDB = bsDBHandle.product();
00230 
00231     // translate from BeamSpotObjects to reco::BeamSpot
00232     BeamSpot::Point apoint( spotDB->GetX(), spotDB->GetY(), spotDB->GetZ() );
00233   
00234     BeamSpot::CovarianceMatrix matrix;
00235     for ( int i=0; i<7; ++i ) {
00236       for ( int j=0; j<7; ++j ) {
00237         matrix(i,j) = spotDB->GetCovariance(i,j);
00238       }
00239     }
00240   
00241     beamSpotsMap_["DB"] = BeamSpot( apoint,
00242                                     spotDB->GetSigmaZ(),
00243                                     spotDB->Getdxdz(),
00244                                     spotDB->Getdydz(),
00245                                     spotDB->GetBeamWidthX(),
00246                                     matrix );
00247 
00248     BeamSpot* aSpot = &(beamSpotsMap_["DB"]);
00249 
00250     aSpot->setBeamWidthY( spotDB->GetBeamWidthY() );
00251     aSpot->setEmittanceX( spotDB->GetEmittanceX() );
00252     aSpot->setEmittanceY( spotDB->GetEmittanceY() );
00253     aSpot->setbetaStar( spotDB->GetBetaStar() );
00254 
00255     if ( spotDB->GetBeamType() == 2 ) {
00256       aSpot->setType( reco::BeamSpot::Tracker );
00257     } else{
00258       aSpot->setType( reco::BeamSpot::Fake );
00259     }
00260     //LogInfo("AlcaBeamMonitor")
00261     //  << *aSpot << std::endl;
00262   }
00263   else {
00264     LogInfo("AlcaBeamMonitor") 
00265       << "Database BeamSpot is not valid at lumi: " << iLumi.id().luminosityBlock(); 
00266   }
00267 }
00268 
00269 //----------------------------------------------------------------------------------------------------------------------
00270 void AlcaBeamMonitor::analyze(const Event& iEvent, const EventSetup& iSetup ){
00271   
00272   //------ BeamFitter 
00273   theBeamFitter_->readEvent(iEvent);
00274   //------ PVFitter 
00275   thePVFitter_->readEvent(iEvent);
00276   
00277   if(beamSpotsMap_.find("DB") != beamSpotsMap_.end()){
00278     //------ Tracks
00279     Handle<reco::TrackCollection> TrackCollection;
00280     iEvent.getByLabel(trackLabel_, TrackCollection);
00281     const reco::TrackCollection *tracks = TrackCollection.product();
00282     for ( reco::TrackCollection::const_iterator track = tracks->begin(); track != tracks->end(); ++track ) {    
00283       hD0Phi0_->Fill(track->phi(), -1*track->dxy(beamSpotsMap_["DB"].position()));
00284       hDxyBS_->Fill(-1*track->dxy(beamSpotsMap_["DB"].position()));
00285     }
00286   }
00287   
00288   //------ Primary Vertices
00289   Handle<VertexCollection > PVCollection;
00290   if (iEvent.getByLabel(primaryVertexLabel_, PVCollection )) {
00291     vertices_.push_back(*PVCollection.product());
00292   }
00293 
00294   if(beamSpotsMap_.find("SC") == beamSpotsMap_.end()){
00295     //BeamSpot from file for this stream is = to the scalar BeamSpot
00296     Handle<BeamSpot> recoBeamSpotHandle;
00297     try{
00298       iEvent.getByLabel(scalerLabel_,recoBeamSpotHandle);
00299     }
00300     catch( cms::Exception& exception ){                               
00301       LogInfo("AlcaBeamMonitor") 
00302       << exception.what(); 
00303       return;         
00304     }                                 
00305     beamSpotsMap_["SC"] = *recoBeamSpotHandle;
00306     if ( beamSpotsMap_["SC"].BeamWidthX() != 0 ) {
00307       beamSpotsMap_["SC"].setType( reco::BeamSpot::Tracker );
00308     } else{
00309       beamSpotsMap_["SC"].setType( reco::BeamSpot::Fake );
00310     }
00311   }
00312 }
00313 
00314 
00315 //----------------------------------------------------------------------------------------------------------------------
00316 void AlcaBeamMonitor::endLuminosityBlock(const LuminosityBlock& iLumi, const EventSetup& iSetup) {
00317   if (theBeamFitter_->runPVandTrkFitter()) {
00318     beamSpotsMap_["BF"] = theBeamFitter_->getBeamSpot();
00319   }
00320   theBeamFitter_->resetTrkVector();
00321   theBeamFitter_->resetLSRange();
00322   theBeamFitter_->resetRefTime();
00323   theBeamFitter_->resetPVFitter();
00324 
00325   if ( thePVFitter_->runFitter() ) {
00326     beamSpotsMap_["PV"] = thePVFitter_->getBeamSpot();
00327   }
00328   thePVFitter_->resetAll();
00329 
00330   //    "PV,BF..."      Value,Error
00331   map<std::string,pair<double,double> >   resultsMap;
00332   vector<pair<double,double> >  vertexResults;
00333   MonitorElement* histo=0;
00334   int position = 0;
00335   for(vector<string>::iterator itV=varNamesV_.begin(); itV!=varNamesV_.end(); itV++){
00336     resultsMap.clear();
00337     for(BeamSpotContainer::iterator itBS = beamSpotsMap_.begin(); itBS != beamSpotsMap_.end(); itBS++){
00338       if(itBS->second.type() == BeamSpot::Tracker){
00339         if(*itV == "x"){
00340           resultsMap[itBS->first] = pair<double,double>(itBS->second.x0(),itBS->second.x0Error());
00341         }
00342         else if(*itV == "y"){
00343           resultsMap[itBS->first] = pair<double,double>(itBS->second.y0(),itBS->second.y0Error());
00344         }
00345         else if(*itV == "z"){
00346           resultsMap[itBS->first] = pair<double,double>(itBS->second.z0(),itBS->second.z0Error());
00347         }
00348         else if(*itV == "sigmaX"){
00349           resultsMap[itBS->first] = pair<double,double>(itBS->second.BeamWidthX(),itBS->second.BeamWidthXError());
00350         }
00351         else if(*itV == "sigmaY"){
00352           resultsMap[itBS->first] = pair<double,double>(itBS->second.BeamWidthY(),itBS->second.BeamWidthYError());
00353         }
00354         else if(*itV == "sigmaZ"){
00355           resultsMap[itBS->first] = pair<double,double>(itBS->second.sigmaZ(),itBS->second.sigmaZ0Error());
00356         }
00357         else{
00358           LogInfo("AlcaBeamMonitor")
00359             << "The histosMap_ has been built with the name " << *itV << " that I can't recognize!";
00360           //assert(0);
00361         }
00362       }
00363     }
00364     vertexResults.clear();
00365     for(vector<VertexCollection>::iterator itPV = vertices_.begin(); itPV != vertices_.end(); itPV++){
00366       if(itPV->size() != 0){
00367         for(VertexCollection::const_iterator pv = itPV->begin(); pv != itPV->end(); pv++) {
00368           if (pv->isFake() || pv->tracksSize()<10)  continue;
00369           if(*itV == "x"){                                                                                    
00370             vertexResults.push_back(pair<double,double>(pv->x(),pv->xError()));       
00371           }                                                                                                           
00372           else if(*itV == "y"){                                                                               
00373             vertexResults.push_back(pair<double,double>(pv->y(),pv->yError()));       
00374           }                                                                                                           
00375           else if(*itV == "z"){                                                                               
00376             vertexResults.push_back(pair<double,double>(pv->z(),pv->zError()));       
00377           }                                                                                                           
00378           else if(*itV != "sigmaX" && *itV != "sigmaY" && *itV != "sigmaZ"){                  
00379             LogInfo("AlcaBeamMonitor")                                                                                
00380               << "The histosMap_ has been built with the name " << *itV << " that I can't recognize!";
00381             //assert(0);                                                                                              
00382           }                                                                                                           
00383         }
00384       }
00385     }
00386 /*
00387   histoByCategoryNames_.insert( pair<string,string>("run",        "Coordinate"));
00388   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex fit-DataBase"));
00389   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex fit-BeamFit"));
00390   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex fit-Scalers"));
00391   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex-DataBase"));
00392   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex-BeamFit"));
00393   histoByCategoryNames_.insert( pair<string,string>("run",        "PrimaryVertex-Scalers"));
00394 
00395   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased BeamSpotFit"));  
00396   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased PrimaryVertex"));
00397   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased DataBase"));     
00398   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased Scalers"));      
00399   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased PrimaryVertex-DataBase fit"));
00400   histoByCategoryNames_.insert( pair<string,string>("lumi",       "Lumibased PrimaryVertex-Scalers fit"));
00401   histoByCategoryNames_.insert( pair<string,string>("validation", "Lumibased Scalers-DataBase fit"));
00402   histoByCategoryNames_.insert( pair<string,string>("validation", "Lumibased PrimaryVertex-DataBase"));
00403   histoByCategoryNames_.insert( pair<string,string>("validation", "Lumibased PrimaryVertex-Scalers"));
00404 */
00405     for(multimap<string,string>::iterator itM=histoByCategoryNames_.begin(); itM!=histoByCategoryNames_.end(); itM++){
00406       if(itM->first == "run" && (histo = histosMap_[*itV][itM->first][itM->second]) == 0){
00407         continue;
00408       }
00409       else if(itM->first != "run"){
00410         position = positionsMap_[*itV][itM->first][itM->second];
00411       }
00412       if(itM->second == "Coordinate"){
00413         if(beamSpotsMap_.find("DB") != beamSpotsMap_.end()){
00414           histo->Fill(resultsMap["DB"].first);
00415         }
00416       }
00417       else if(itM->second == "PrimaryVertex fit-DataBase"){
00418         if(resultsMap.find("PV") != resultsMap.end() && resultsMap.find("DB") != resultsMap.end()){
00419           histo->Fill(resultsMap["PV"].first-resultsMap["DB"].first);
00420         }
00421       }
00422       else if(itM->second == "PrimaryVertex fit-BeamFit"){
00423         if(resultsMap.find("PV") != resultsMap.end() && resultsMap.find("BF") != resultsMap.end()){
00424           histo->Fill(resultsMap["PV"].first-resultsMap["BF"].first);
00425         }
00426       }
00427       else if(itM->second == "PrimaryVertex fit-Scalers"){
00428         if(resultsMap.find("PV") != resultsMap.end() && resultsMap.find("SC") != resultsMap.end()){
00429           histo->Fill(resultsMap["PV"].first-resultsMap["SC"].first);
00430         }
00431       }
00432       else if(itM->second == "PrimaryVertex-DataBase"){
00433         if(resultsMap.find("PV") != resultsMap.end() && resultsMap.find("DB") != resultsMap.end()){
00434           for(vector<pair<double,double> >::iterator itPV=vertexResults.begin(); itPV!=vertexResults.end(); itPV++){
00435             histo->Fill(itPV->first-resultsMap["DB"].first);
00436           }
00437         }
00438       }
00439       else if(itM->second == "PrimaryVertex-BeamFit"){
00440         if(resultsMap.find("PV") != resultsMap.end() && resultsMap.find("BF") != resultsMap.end()){
00441           for(vector<pair<double,double> >::iterator itPV=vertexResults.begin(); itPV!=vertexResults.end(); itPV++){
00442             histo->Fill(itPV->first-resultsMap["BF"].first);
00443           }
00444         }
00445       }
00446       else if(itM->second == "PrimaryVertex-Scalers"){
00447         if(resultsMap.find("PV") != resultsMap.end() && resultsMap.find("SC") != resultsMap.end()){
00448           for(vector<pair<double,double> >::iterator itPV=vertexResults.begin(); itPV!=vertexResults.end(); itPV++){
00449             histo->Fill(itPV->first-resultsMap["SC"].first);
00450           }
00451         }
00452       }
00453       else if(itM->second == "Lumibased BeamSpotFit"){
00454         if(resultsMap.find("BF") != resultsMap.end()){
00455           theValuesContainer_->Fill(position  ,resultsMap["BF"].first);//Value
00456           theValuesContainer_->Fill(position+1,resultsMap["BF"].second);//Error
00457           theValuesContainer_->Fill(position+2,1);//ok
00458         }
00459       }
00460       else if(itM->second == "Lumibased PrimaryVertex"){
00461         if(resultsMap.find("PV") != resultsMap.end()){
00462           theValuesContainer_->Fill(position  ,resultsMap["PV"].first);//Value
00463           theValuesContainer_->Fill(position+1,resultsMap["PV"].second);//Error
00464           theValuesContainer_->Fill(position+2,1);//ok
00465         }
00466       }
00467       else if(itM->second == "Lumibased DataBase"){
00468         if(resultsMap.find("DB") != resultsMap.end()){
00469           theValuesContainer_->Fill(position  ,resultsMap["DB"].first);//Value
00470           theValuesContainer_->Fill(position+1,resultsMap["DB"].second);//Error           
00471           theValuesContainer_->Fill(position+2,1);//ok
00472         }
00473       }
00474       else if(itM->second == "Lumibased Scalers"){
00475         if(resultsMap.find("SC") != resultsMap.end()){
00476           theValuesContainer_->Fill(position  ,resultsMap["SC"].first);//Value
00477           theValuesContainer_->Fill(position+1,resultsMap["SC"].second);//Error           
00478           theValuesContainer_->Fill(position+2,1);//ok
00479         }
00480       }
00481       else if(itM->second == "Lumibased PrimaryVertex-DataBase fit"){
00482         if(resultsMap.find("PV") != resultsMap.end() && resultsMap.find("DB") != resultsMap.end()){
00483           theValuesContainer_->Fill(position  ,resultsMap["PV"].first-resultsMap["DB"].first);//Value
00484           theValuesContainer_->Fill(position+1,std::sqrt(std::pow(resultsMap["PV"].second,2)+std::pow(resultsMap["DB"].second,2)));//Error        
00485           theValuesContainer_->Fill(position+2,1);//ok
00486         }
00487       }
00488       else if(itM->second == "Lumibased PrimaryVertex-Scalers fit"){
00489         if(resultsMap.find("PV") != resultsMap.end() && resultsMap.find("SC") != resultsMap.end()){
00490           theValuesContainer_->Fill(position  ,resultsMap["PV"].first-resultsMap["SC"].first);//Value
00491           theValuesContainer_->Fill(position+1,std::sqrt(std::pow(resultsMap["PV"].second,2)+std::pow(resultsMap["SC"].second,2)));//Error        
00492           theValuesContainer_->Fill(position+2,1);//ok
00493         }
00494       }
00495       else if(itM->second == "Lumibased Scalers-DataBase fit"){
00496         if(resultsMap.find("SC") != resultsMap.end() && resultsMap.find("DB") != resultsMap.end()){
00497           theValuesContainer_->Fill(position  ,resultsMap["SC"].first-resultsMap["DB"].first);//Value
00498           theValuesContainer_->Fill(position+1,std::sqrt(std::pow(resultsMap["SC"].second,2)+std::pow(resultsMap["DB"].second,2)));//Error        
00499           theValuesContainer_->Fill(position+2,1);//ok
00500         }
00501       }
00502       else if(itM->second == "Lumibased PrimaryVertex-DataBase"){
00503         if(resultsMap.find("DB") != resultsMap.end() && vertexResults.size() != 0){
00504           for(vector<pair<double,double> >::iterator itPV=vertexResults.begin(); itPV!=vertexResults.end(); itPV++){
00505             theValuesContainer_->Fill(position  ,(*itPV).first-resultsMap["DB"].first);//Value
00506           }
00507 /*
00508           double error = 0;
00509           if(vertexResults.size() != 0){
00510             for(vector<pair<double,double> >::iterator itPV=vertexResults.begin(); itPV!=vertexResults.end(); itPV++){
00511               error += std::pow((*itPV).first-resultsMap["DB"].first-theValuesContainer_->getTProfile()->GetBinContent(position+1),2.);
00512             }
00513             error = std::sqrt(error)/vertexResults.size();
00514           }
00515 //          theValuesContainer_->Fill(position+1,std::sqrt(std::pow((*itPV).second,2)+std::pow(resultsMap["DB"].second,2)));//Error       
00516           theValuesContainer_->Fill(position+1,error);//Error     
00517 */
00518           theValuesContainer_->Fill(position+1,theValuesContainer_->getTProfile()->GetBinError(position+1));//Error       
00519           theValuesContainer_->Fill(position+2,1);//ok
00520         }
00521       }
00522       else if(itM->second == "Lumibased PrimaryVertex-Scalers"){
00523         if(resultsMap.find("SC") != resultsMap.end() && vertexResults.size() != 0){
00524           for(vector<pair<double,double> >::iterator itPV=vertexResults.begin(); itPV!=vertexResults.end(); itPV++){
00525             theValuesContainer_->Fill(position  ,(*itPV).first-resultsMap["SC"].first);//Value
00526           }
00527 /*
00528           double error = 0;
00529           if(vertexResults.size() != 0){
00530             for(vector<pair<double,double> >::iterator itPV=vertexResults.begin(); itPV!=vertexResults.end(); itPV++){
00531               error += std::pow((*itPV).first-resultsMap["SC"].first-theValuesContainer_->getTProfile()->GetBinContent(position+1),2.);
00532             }
00533             error = std::sqrt(error)/vertexResults.size();
00534           }
00535 //          theValuesContainer_->Fill(position+1,std::sqrt(std::pow((*itPV).second,2)+std::pow(resultsMap["SC"].second,2)));//Error       
00536           theValuesContainer_->Fill(position+1,error);//Error     
00537 */
00538           theValuesContainer_->Fill(position+1,theValuesContainer_->getTProfile()->GetBinError(position+1));//Error       
00539           theValuesContainer_->Fill(position+2,1);//ok
00540         }
00541       }
00542 //      else if(itM->second == "Lumibased Scalers-DataBase"){
00543 //      if(resultsMap.find("SC") != resultsMap.end() && resultsMap.find("DB") != resultsMap.end()){
00544 //        itHHH->second->Fill(bin,resultsMap["SC"].first-resultsMap["DB"].first);
00545 //      }
00546 //    }
00547       else{
00548         LogInfo("AlcaBeamMonitor")
00549           << "The histosMap_ have a histogram named " << itM->second << " that I can't recognize in this loop!";
00550         //assert(0);
00551 
00552       }
00553     }
00554   }
00555 }
00556 
00557 //----------------------------------------------------------------------------------------------------------------------
00558 void AlcaBeamMonitor::endRun(const Run& iRun, const EventSetup& context){
00559 }
00560 
00561 //----------------------------------------------------------------------------------------------------------------------
00562 void AlcaBeamMonitor::endJob(const LuminosityBlock& iLumi, const EventSetup& iSetup){
00563 }
00564 
00565 
00566 DEFINE_FWK_MODULE(AlcaBeamMonitor);