CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/RecoVertex/BeamSpotProducer/src/BeamFitter.cc

Go to the documentation of this file.
00001 
00014 #include "RecoVertex/BeamSpotProducer/interface/BeamFitter.h"
00015 #include "FWCore/ServiceRegistry/interface/Service.h"
00016 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
00017 #include "CondFormats/BeamSpotObjects/interface/BeamSpotObjects.h"
00018 
00019 #include "FWCore/Utilities/interface/InputTag.h"
00020 #include "DataFormats/Common/interface/View.h"
00021 
00022 #include "DataFormats/TrackCandidate/interface/TrackCandidate.h"
00023 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
00024 #include "DataFormats/TrackReco/interface/Track.h"
00025 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00026 #include "DataFormats/TrackReco/interface/HitPattern.h"
00027 #include "DataFormats/VertexReco/interface/Vertex.h"
00028 
00029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00030 
00031 // ----------------------------------------------------------------------
00032 // Useful function:
00033 // ----------------------------------------------------------------------
00034 const char * BeamFitter::formatBTime(const std::time_t & t)  {
00035   struct std::tm * ptm;
00036   ptm = gmtime(&t);
00037   static char ts[] = "yyyy.mn.dd hh:mm:ss zzz ";
00038   strftime(ts,sizeof(ts),"%Y.%m.%d %H:%M:%S %Z",ptm);
00039   return ts;
00040 }
00041 
00042 BeamFitter::BeamFitter(const edm::ParameterSet& iConfig): fPVTree_(0)
00043 {
00044 
00045   debug_             = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<bool>("Debug");
00046   tracksLabel_       = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<edm::InputTag>("TrackCollection");
00047   vertexLabel_       = iConfig.getUntrackedParameter<edm::InputTag>("primaryVertex", edm::InputTag("offlinePrimaryVertices"));
00048   writeTxt_          = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<bool>("WriteAscii");
00049   outputTxt_         = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<std::string>("AsciiFileName");
00050   appendRunTxt_      = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<bool>("AppendRunToFileName");
00051   writeDIPTxt_       = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<bool>("WriteDIPAscii");
00052   outputDIPTxt_      = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<std::string>("DIPFileName");
00053   saveNtuple_        = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<bool>("SaveNtuple");
00054   saveBeamFit_       = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<bool>("SaveFitResults");
00055   savePVVertices_    = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<bool>("SavePVVertices");
00056   isMuon_            = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<bool>("IsMuonCollection");
00057 
00058   trk_MinpT_         = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<double>("MinimumPt");
00059   trk_MaxEta_        = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<double>("MaximumEta");
00060   trk_MaxIP_         = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<double>("MaximumImpactParameter");
00061   trk_MaxZ_          = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<double>("MaximumZ");
00062   trk_MinNTotLayers_ = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<int>("MinimumTotalLayers");
00063   trk_MinNPixLayers_ = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<int>("MinimumPixelLayers");
00064   trk_MaxNormChi2_   = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<double>("MaximumNormChi2");
00065   trk_Algorithm_     = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<std::vector<std::string> >("TrackAlgorithm");
00066   trk_Quality_       = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<std::vector<std::string> >("TrackQuality");
00067   min_Ntrks_         = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<int>("MinimumInputTracks");
00068   convergence_       = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<double>("FractionOfFittedTrks");
00069   inputBeamWidth_    = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<double>("InputBeamWidth",-1.);
00070 
00071   for (unsigned int j=0;j<trk_Algorithm_.size();j++)
00072     algorithm_.push_back(reco::TrackBase::algoByName(trk_Algorithm_[j]));
00073   for (unsigned int j=0;j<trk_Quality_.size();j++)
00074     quality_.push_back(reco::TrackBase::qualityByName(trk_Quality_[j]));
00075 
00076   if (saveNtuple_ || saveBeamFit_ || savePVVertices_) {
00077     outputfilename_ = iConfig.getParameter<edm::ParameterSet>("BeamFitter").getUntrackedParameter<std::string>("OutputFileName");
00078     file_ = TFile::Open(outputfilename_.c_str(),"RECREATE");
00079   }
00080   if (saveNtuple_) {
00081     ftree_ = new TTree("mytree","mytree");
00082     ftree_->AutoSave();
00083 
00084     ftree_->Branch("pt",&fpt,"fpt/D");
00085     ftree_->Branch("d0",&fd0,"fd0/D");
00086     ftree_->Branch("d0bs",&fd0bs,"fd0bs/D");
00087     ftree_->Branch("sigmad0",&fsigmad0,"fsigmad0/D");
00088     ftree_->Branch("phi0",&fphi0,"fphi0/D");
00089     ftree_->Branch("z0",&fz0,"fz0/D");
00090     ftree_->Branch("sigmaz0",&fsigmaz0,"fsigmaz0/D");
00091     ftree_->Branch("theta",&ftheta,"ftheta/D");
00092     ftree_->Branch("eta",&feta,"feta/D");
00093     ftree_->Branch("charge",&fcharge,"fcharge/I");
00094     ftree_->Branch("normchi2",&fnormchi2,"fnormchi2/D");
00095     ftree_->Branch("nTotLayerMeas",&fnTotLayerMeas,"fnTotLayerMeas/i");
00096     ftree_->Branch("nStripLayerMeas",&fnStripLayerMeas,"fnStripLayerMeas/i");
00097     ftree_->Branch("nPixelLayerMeas",&fnPixelLayerMeas,"fnPixelLayerMeas/i");
00098     ftree_->Branch("nTIBLayerMeas",&fnTIBLayerMeas,"fnTIBLayerMeas/i");
00099     ftree_->Branch("nTOBLayerMeas",&fnTOBLayerMeas,"fnTOBLayerMeas/i");
00100     ftree_->Branch("nTIDLayerMeas",&fnTIDLayerMeas,"fnTIDLayerMeas/i");
00101     ftree_->Branch("nTECLayerMeas",&fnTECLayerMeas,"fnTECLayerMeas/i");
00102     ftree_->Branch("nPXBLayerMeas",&fnPXBLayerMeas,"fnPXBLayerMeas/i");
00103     ftree_->Branch("nPXFLayerMeas",&fnPXFLayerMeas,"fnPXFLayerMeas/i");
00104     ftree_->Branch("cov",&fcov,"fcov[7][7]/D");
00105     ftree_->Branch("vx",&fvx,"fvx/D");
00106     ftree_->Branch("vy",&fvy,"fvy/D");
00107     ftree_->Branch("quality",&fquality,"fquality/O");
00108     ftree_->Branch("algo",&falgo,"falgo/O");
00109     ftree_->Branch("run",&frun,"frun/i");
00110     ftree_->Branch("lumi",&flumi,"flumi/i");
00111     ftree_->Branch("pvValid",&fpvValid,"fpvValid/O");
00112     ftree_->Branch("pvx", &fpvx, "fpvx/D");
00113     ftree_->Branch("pvy", &fpvy, "fpvy/D");
00114     ftree_->Branch("pvz", &fpvz, "fpvz/D");
00115   }
00116   if (saveBeamFit_){
00117     ftreeFit_ = new TTree("fitResults","fitResults");
00118     ftreeFit_->AutoSave();
00119     ftreeFit_->Branch("run",&frunFit,"frunFit/i");
00120     ftreeFit_->Branch("beginLumi",&fbeginLumiOfFit,"fbeginLumiOfFit/i");
00121     ftreeFit_->Branch("endLumi",&fendLumiOfFit,"fendLumiOfFit/i");
00122     ftreeFit_->Branch("beginTime",fbeginTimeOfFit,"fbeginTimeOfFit/C");
00123     ftreeFit_->Branch("endTime",fendTimeOfFit,"fendTimeOfFit/C");
00124     ftreeFit_->Branch("x",&fx,"fx/D");
00125     ftreeFit_->Branch("y",&fy,"fy/D");
00126     ftreeFit_->Branch("z",&fz,"fz/D");
00127     ftreeFit_->Branch("sigmaZ",&fsigmaZ,"fsigmaZ/D");
00128     ftreeFit_->Branch("dxdz",&fdxdz,"fdxdz/D");
00129     ftreeFit_->Branch("dydz",&fdydz,"fdydz/D");
00130     ftreeFit_->Branch("xErr",&fxErr,"fxErr/D");
00131     ftreeFit_->Branch("yErr",&fyErr,"fyErr/D");
00132     ftreeFit_->Branch("zErr",&fzErr,"fzErr/D");
00133     ftreeFit_->Branch("sigmaZErr",&fsigmaZErr,"fsigmaZErr/D");
00134     ftreeFit_->Branch("dxdzErr",&fdxdzErr,"fdxdzErr/D");
00135     ftreeFit_->Branch("dydzErr",&fdydzErr,"fdydzErr/D");
00136     ftreeFit_->Branch("widthX",&fwidthX,"fwidthX/D");
00137     ftreeFit_->Branch("widthY",&fwidthY,"fwidthY/D");
00138     ftreeFit_->Branch("widthXErr",&fwidthXErr,"fwidthXErr/D");
00139     ftreeFit_->Branch("widthYErr",&fwidthYErr,"fwidthYErr/D");
00140   }
00141 
00142   fBSvector.clear();
00143   ftotal_tracks = 0;
00144   fnTotLayerMeas = fnPixelLayerMeas = fnStripLayerMeas = fnTIBLayerMeas = 0;
00145   fnTIDLayerMeas = fnTOBLayerMeas = fnTECLayerMeas = fnPXBLayerMeas = fnPXFLayerMeas = 0;
00146   frun = flumi = -1;
00147   frunFit = fbeginLumiOfFit = fendLumiOfFit = -1;
00148   fquality = falgo = true;
00149   fpvValid = true;
00150   fpvx = fpvy = fpvz = 0;
00151   fitted_ = false;
00152   resetRefTime();
00153 
00154   //debug histograms
00155   h1ntrks = new TH1F("h1ntrks","number of tracks per event",50,0,50);
00156   h1vz_event = new TH1F("h1vz_event","track Vz", 50, -30, 30 );
00157   h1cutFlow = new TH1F("h1cutFlow","Cut flow table of track selection", 9, 0, 9);
00158   h1cutFlow->GetXaxis()->SetBinLabel(1,"No cut");
00159   h1cutFlow->GetXaxis()->SetBinLabel(2,"Traker hits");
00160   h1cutFlow->GetXaxis()->SetBinLabel(3,"Pixel hits");
00161   h1cutFlow->GetXaxis()->SetBinLabel(4,"norm. #chi^{2}");
00162   h1cutFlow->GetXaxis()->SetBinLabel(5,"algo");
00163   h1cutFlow->GetXaxis()->SetBinLabel(6,"quality");
00164   h1cutFlow->GetXaxis()->SetBinLabel(7,"d_{0}");
00165   h1cutFlow->GetXaxis()->SetBinLabel(8,"z_{0}");
00166   h1cutFlow->GetXaxis()->SetBinLabel(9,"p_{T}");
00167   resetCutFlow();
00168 
00169   // Primary vertex fitter
00170   MyPVFitter = new PVFitter(iConfig);
00171   MyPVFitter->resetAll();
00172   if (savePVVertices_){
00173     fPVTree_ = new TTree("PrimaryVertices","PrimaryVertices");
00174     MyPVFitter->setTree(fPVTree_);
00175   }
00176 
00177   // check filename
00178   ffilename_changed = false;
00179   if (writeDIPTxt_)
00180     fasciiDIP.open(outputDIPTxt_.c_str());
00181 }
00182 
00183 BeamFitter::~BeamFitter() {
00184 
00185   if (saveNtuple_) {
00186     file_->cd();
00187     if (fitted_ && h1z) h1z->Write();
00188     h1ntrks->Write();
00189     h1vz_event->Write();
00190     if (h1cutFlow) h1cutFlow->Write();
00191     ftree_->Write();
00192   }
00193   if (saveBeamFit_){
00194     file_->cd();
00195     ftreeFit_->Write();
00196   }
00197   if (savePVVertices_){
00198     file_->cd();
00199     fPVTree_->Write();
00200   }
00201 
00202 
00203   if (saveNtuple_ || saveBeamFit_ || savePVVertices_){
00204     file_->Close();
00205     delete file_;
00206   }
00207   delete MyPVFitter;
00208 }
00209 
00210 
00211 void BeamFitter::readEvent(const edm::Event& iEvent)
00212 {
00213 
00214 
00215   frun = iEvent.id().run();
00216   const edm::TimeValue_t ftimestamp = iEvent.time().value();
00217   const std::time_t ftmptime = ftimestamp >> 32;
00218 
00219   if (fbeginLumiOfFit == -1) freftime[0] = freftime[1] = ftmptime;
00220   if (freftime[0] == 0 || ftmptime < freftime[0]) freftime[0] = ftmptime;
00221   if (freftime[1] == 0 || ftmptime > freftime[1]) freftime[1] = ftmptime;
00222   flumi = iEvent.luminosityBlock();
00223   frunFit = frun;
00224   if (fbeginLumiOfFit == -1 || fbeginLumiOfFit > flumi) fbeginLumiOfFit = flumi;
00225   if (fendLumiOfFit == -1 || fendLumiOfFit < flumi) fendLumiOfFit = flumi;
00226 
00227   edm::Handle<reco::TrackCollection> TrackCollection;
00228   iEvent.getByLabel(tracksLabel_, TrackCollection);
00229 
00230   //------ Primary Vertices
00231   edm::Handle< edm::View<reco::Vertex> > PVCollection;
00232   bool hasPVs = false;
00233   edm::View<reco::Vertex> pv;
00234   if ( iEvent.getByLabel(vertexLabel_, PVCollection ) ) {
00235       pv = *PVCollection;
00236       hasPVs = true;
00237   }
00238   //------
00239 
00240   //------ Beam Spot in current event
00241   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00242   const reco::BeamSpot *refBS =  0;
00243   if ( iEvent.getByLabel("offlineBeamSpot",recoBeamSpotHandle) )
00244       refBS = recoBeamSpotHandle.product();
00245   //-------
00246 
00247   const reco::TrackCollection *tracks = TrackCollection.product();
00248 
00249   double eventZ = 0;
00250   double averageZ = 0;
00251 
00252   for ( reco::TrackCollection::const_iterator track = tracks->begin();
00253         track != tracks->end();
00254         ++track ) {
00255 
00256     if ( ! isMuon_) {
00257 
00258       const reco::HitPattern& trkHP = track->hitPattern();
00259 
00260       fnPixelLayerMeas = trkHP.pixelLayersWithMeasurement();
00261       fnStripLayerMeas = trkHP.stripLayersWithMeasurement();
00262       fnTotLayerMeas = trkHP.trackerLayersWithMeasurement();
00263       fnPXBLayerMeas = trkHP.pixelBarrelLayersWithMeasurement();
00264       fnPXFLayerMeas = trkHP.pixelEndcapLayersWithMeasurement();
00265       fnTIBLayerMeas = trkHP.stripTIBLayersWithMeasurement();
00266       fnTIDLayerMeas = trkHP.stripTIDLayersWithMeasurement();
00267       fnTOBLayerMeas = trkHP.stripTOBLayersWithMeasurement();
00268       fnTECLayerMeas = trkHP.stripTECLayersWithMeasurement();
00269     } else {
00270 
00271       fnTotLayerMeas = track->numberOfValidHits();
00272 
00273     }
00274 
00275     fpt = track->pt();
00276     feta = track->eta();
00277     fphi0 = track->phi();
00278     fcharge = track->charge();
00279     fnormchi2 = track->normalizedChi2();
00280     fd0 = track->d0();
00281     if (refBS) fd0bs = -1*track->dxy(refBS->position());
00282     else fd0bs = 0.;
00283 
00284     fsigmad0 = track->d0Error();
00285     fz0 = track->dz();
00286     fsigmaz0 = track->dzError();
00287     ftheta = track->theta();
00288     fvx = track->vx();
00289     fvy = track->vy();
00290 
00291     for (int i=0; i<5; ++i) {
00292       for (int j=0; j<5; ++j) {
00293         fcov[i][j] = track->covariance(i,j);
00294       }
00295     }
00296 
00297     fquality = true;
00298     falgo = true;
00299 
00300     if (! isMuon_ ) {
00301       if (quality_.size()!=0) {
00302           fquality = false;
00303           for (unsigned int i = 0; i<quality_.size();++i) {
00304               if(debug_) edm::LogInfo("BeamFitter") << "quality_[" << i << "] = " << track->qualityName(quality_[i]) << std::endl;
00305               if (track->quality(quality_[i])) {
00306                   fquality = true;
00307                   break;
00308               }
00309           }
00310       }
00311 
00312 
00313       // Track algorithm
00314 
00315       if (algorithm_.size()!=0) {
00316         if (std::find(algorithm_.begin(),algorithm_.end(),track->algo())==algorithm_.end())
00317           falgo = false;
00318       }
00319 
00320     }
00321 
00322     // check if we have a valid PV
00323     fpvValid = false;
00324 
00325     if ( hasPVs ) {
00326 
00327         for ( size_t ipv=0; ipv != pv.size(); ++ipv ) {
00328 
00329             if (! pv[ipv].isFake() ) fpvValid = true;
00330 
00331             if ( ipv==0 && !pv[0].isFake() ) { fpvx = pv[0].x(); fpvy = pv[0].y(); fpvz = pv[0].z(); } // fix this later
00332 
00333 
00334         }
00335 
00336     }
00337 
00338 
00339     if (saveNtuple_) ftree_->Fill();
00340     ftotal_tracks++;
00341 
00342     countPass[0] = ftotal_tracks;
00343     // Track selection
00344     if (fnTotLayerMeas >= trk_MinNTotLayers_) { countPass[1] += 1;
00345       if (fnPixelLayerMeas >= trk_MinNPixLayers_) { countPass[2] += 1;
00346         if (fnormchi2 < trk_MaxNormChi2_) { countPass[3] += 1;
00347           if (falgo) {countPass[4] += 1;
00348             if (fquality) { countPass[5] += 1;
00349               if (std::abs( fd0 ) < trk_MaxIP_) { countPass[6] += 1;
00350                 if (std::abs( fz0 ) < trk_MaxZ_){ countPass[7] += 1;
00351                   if (fpt > trk_MinpT_) {
00352                     countPass[8] += 1;
00353                     if (std::abs( feta ) < trk_MaxEta_
00354                         //&& fpvValid
00355                         ) {
00356                       if (debug_) {
00357                   edm::LogInfo("BeamFitter") << "Selected track quality = " << track->qualityMask()
00358                                              << "; track algorithm = " << track->algoName() << "= TrackAlgorithm: " << track->algo() << std::endl;
00359                       }
00360                       BSTrkParameters BSTrk(fz0,fsigmaz0,fd0,fsigmad0,fphi0,fpt,0.,0.);
00361                       BSTrk.setVx(fvx);
00362                       BSTrk.setVy(fvy);
00363                       fBSvector.push_back(BSTrk);
00364                       averageZ += fz0;
00365                     }
00366                   }
00367                 }
00368               }
00369             }
00370           }
00371         }
00372       }
00373     }// track selection
00374 
00375   }// tracks
00376 
00377   averageZ = averageZ/(float)(fBSvector.size());
00378 
00379   for( std::vector<BSTrkParameters>::const_iterator iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam) {
00380 
00381     eventZ += fabs( iparam->z0() - averageZ );
00382 
00383   }
00384 
00385   h1ntrks->Fill( fBSvector.size() );
00386   h1vz_event->Fill( eventZ/(float)(fBSvector.size() ) ) ;
00387   for (unsigned int i=0; i < sizeof(countPass)/sizeof(countPass[0]); i++)
00388     h1cutFlow->SetBinContent(i+1,countPass[i]);
00389 
00390   MyPVFitter->readEvent(iEvent);
00391 
00392 }
00393 
00394 bool BeamFitter::runPVandTrkFitter() {
00395 // run both PV and track fitters
00396     bool fit_ok = false;
00397     bool pv_fit_ok = false;
00398     reco::BeamSpot bspotPV;
00399     reco::BeamSpot bspotTrk;
00400 
00401     // First run PV fitter
00402     if ( MyPVFitter->IsFitPerBunchCrossing() ){
00403       const char* fbeginTime = formatBTime(freftime[0]);
00404       sprintf(fbeginTimeOfFit,"%s",fbeginTime);
00405       const char* fendTime = formatBTime(freftime[1]);
00406       sprintf(fendTimeOfFit,"%s",fendTime);
00407 
00408       edm::LogInfo("BeamFitter") << " [BeamFitterBxDebugTime] freftime[0] = " << freftime[0]
00409                                  << "; address =  " << &freftime[0]
00410                                  << " = " << fbeginTimeOfFit << std::endl;
00411       edm::LogInfo("BeamFitter") << " [BeamFitterBxDebugTime] freftime[1] = " << freftime[1]
00412                                  << "; address =  " << &freftime[1]
00413                                  << " = " << fendTimeOfFit << std::endl;
00414 
00415       if ( MyPVFitter->runBXFitter() ) {
00416         fbspotPVMap = MyPVFitter->getBeamSpotMap();
00417         pv_fit_ok = true;
00418       }
00419       if(writeTxt_ ) dumpTxtFile(outputTxt_,true); // all reaults
00420       if(writeDIPTxt_) dumpTxtFile(outputDIPTxt_,false); // for DQM/DIP
00421       return pv_fit_ok;
00422     }
00423 
00424     if ( MyPVFitter->runFitter() ) {
00425 
00426         bspotPV = MyPVFitter->getBeamSpot();
00427 
00428         // take beam width from PV fit and pass it to track fitter
00429         // assume circular transverse beam width
00430         inputBeamWidth_ = sqrt( pow(bspotPV.BeamWidthX(),2) + pow(bspotPV.BeamWidthY(),2) )/sqrt(2);
00431         pv_fit_ok = true;
00432 
00433     } else {
00434         // problems with PV fit
00435         bspotPV.setType(reco::BeamSpot::Unknown);
00436         bspotTrk.setType(reco::BeamSpot::Unknown); //propagate error to trk beam spot
00437     }
00438 
00439     if ( runFitterNoTxt() ) {
00440 
00441         bspotTrk = fbeamspot;
00442         fit_ok = true;
00443     } else {
00444       // beamfit failed, flagged as empty beam spot
00445       bspotTrk.setType(reco::BeamSpot::Fake);
00446       fit_ok = false;
00447     }
00448 
00449     // combined results into one single beam spot
00450 
00451     // to pass errors I have to create another beam spot object
00452     reco::BeamSpot::CovarianceMatrix matrix;
00453     for (int j = 0 ; j < 7 ; ++j) {
00454         for(int k = j ; k < 7 ; ++k) {
00455             matrix(j,k) = bspotTrk.covariance(j,k);
00456         }
00457     }
00458     // change beam width error to one from PV
00459     if (pv_fit_ok && fit_ok ) {
00460       matrix(6,6) = MyPVFitter->getWidthXerr() * MyPVFitter->getWidthXerr();
00461 
00462       // get Z and sigmaZ from PV fit
00463       matrix(2,2) = bspotPV.covariance(2,2);
00464       matrix(3,3) = bspotPV.covariance(3,3);
00465       reco::BeamSpot tmpbs(reco::BeamSpot::Point(bspotTrk.x0(), bspotTrk.y0(),
00466                                                  bspotPV.z0() ),
00467                            bspotPV.sigmaZ() ,
00468                            bspotTrk.dxdz(),
00469                            bspotTrk.dydz(),
00470                            bspotPV.BeamWidthX(),
00471                            matrix,
00472                            bspotPV.type() );
00473       tmpbs.setBeamWidthY( bspotPV.BeamWidthY() );
00474       // overwrite beam spot result
00475       fbeamspot = tmpbs;
00476     }
00477     if (pv_fit_ok && fit_ok) {
00478       fbeamspot.setType(bspotPV.type());
00479     }
00480     else if(!pv_fit_ok && fit_ok){
00481       fbeamspot.setType(reco::BeamSpot::Unknown);
00482     }
00483     else if(pv_fit_ok && !fit_ok){
00484       fbeamspot.setType(reco::BeamSpot::Unknown);
00485     }
00486     else if(!pv_fit_ok && !fit_ok){
00487       fbeamspot.setType(reco::BeamSpot::Fake);
00488     }
00489 
00490     if(writeTxt_ ) dumpTxtFile(outputTxt_,true); // all reaults
00491     if(writeDIPTxt_) dumpTxtFile(outputDIPTxt_,false); // for DQM/DIP
00492 
00493     return fit_ok && pv_fit_ok;
00494 }
00495 
00496 bool BeamFitter::runFitterNoTxt() {
00497   const char* fbeginTime = formatBTime(freftime[0]);
00498   sprintf(fbeginTimeOfFit,"%s",fbeginTime);
00499   const char* fendTime = formatBTime(freftime[1]);
00500   sprintf(fendTimeOfFit,"%s",fendTime);
00501 
00502   edm::LogInfo("BeamFitter") << " [BeamFitterDebugTime] freftime[0] = " << freftime[0]
00503                              << "; address =  " << &freftime[0]
00504                              << " = " << fbeginTimeOfFit << std::endl;
00505   edm::LogInfo("BeamFitter") << " [BeamFitterDebugTime] freftime[1] = " << freftime[1]
00506                              << "; address =  " << &freftime[1]
00507                              << " = " << fendTimeOfFit << std::endl;
00508 
00509   if (fbeginLumiOfFit == -1 || fendLumiOfFit == -1) {
00510       edm::LogWarning("BeamFitter") << "No event read! No Fitting!" << std::endl;
00511     return false;
00512   }
00513 
00514   bool fit_ok = false;
00515   // default fit to extract beam spot info
00516   if(fBSvector.size() > 1 ){
00517 
00518       edm::LogInfo("BeamFitter") << "Calculating beam spot..." << std::endl
00519                                  << "We will use " << fBSvector.size() << " good tracks out of " << ftotal_tracks << std::endl;
00520 
00521     BSFitter *myalgo = new BSFitter( fBSvector );
00522     myalgo->SetMaximumZ( trk_MaxZ_ );
00523     myalgo->SetConvergence( convergence_ );
00524     myalgo->SetMinimumNTrks( min_Ntrks_ );
00525     if (inputBeamWidth_ > 0 ) myalgo->SetInputBeamWidth( inputBeamWidth_ );
00526 
00527     fbeamspot = myalgo->Fit();
00528 
00529 
00530     // retrieve histogram for Vz
00531     h1z = (TH1F*) myalgo->GetVzHisto();
00532 
00533     delete myalgo;
00534     if ( fbeamspot.type() != 0 ) {// save all results except for Fake (all 0.)
00535       fit_ok = true;
00536       if (saveBeamFit_){
00537         fx = fbeamspot.x0();
00538         fy = fbeamspot.y0();
00539         fz = fbeamspot.z0();
00540         fsigmaZ = fbeamspot.sigmaZ();
00541         fdxdz = fbeamspot.dxdz();
00542         fdydz = fbeamspot.dydz();
00543         fwidthX = fbeamspot.BeamWidthX();
00544         fwidthY = fbeamspot.BeamWidthY();
00545         fxErr = fbeamspot.x0Error();
00546         fyErr = fbeamspot.y0Error();
00547         fzErr = fbeamspot.z0Error();
00548         fsigmaZErr = fbeamspot.sigmaZ0Error();
00549         fdxdzErr = fbeamspot.dxdzError();
00550         fdydzErr = fbeamspot.dydzError();
00551         fwidthXErr = fbeamspot.BeamWidthXError();
00552         fwidthYErr = fbeamspot.BeamWidthYError();
00553         ftreeFit_->Fill();
00554       }
00555     }
00556   }
00557   else{ // tracks <= 1
00558     reco::BeamSpot tmpbs;
00559     fbeamspot = tmpbs;
00560     fbeamspot.setType(reco::BeamSpot::Fake);
00561     edm::LogInfo("BeamFitter") << "Not enough good tracks selected! No beam fit!" << std::endl;
00562 
00563   }
00564   fitted_ = true;
00565   return fit_ok;
00566 
00567 }
00568 
00569 bool BeamFitter::runFitter() {
00570 
00571     bool fit_ok = runFitterNoTxt();
00572 
00573     if(writeTxt_ ) dumpTxtFile(outputTxt_,true); // all reaults
00574     if(writeDIPTxt_) dumpTxtFile(outputDIPTxt_,false); // for DQM/DIP
00575 
00576     return fit_ok;
00577 }
00578 
00579 bool BeamFitter::runBeamWidthFitter() {
00580   bool widthfit_ok = false;
00581   // default fit to extract beam spot info
00582   if(fBSvector.size() > 1 ){
00583 
00584       edm::LogInfo("BeamFitter") << "Calculating beam spot positions('d0-phi0' method) and width using llh Fit"<< std::endl
00585                                  << "We will use " << fBSvector.size() << " good tracks out of " << ftotal_tracks << std::endl;
00586 
00587         BSFitter *myalgo = new BSFitter( fBSvector );
00588         myalgo->SetMaximumZ( trk_MaxZ_ );
00589         myalgo->SetConvergence( convergence_ );
00590         myalgo->SetMinimumNTrks(min_Ntrks_);
00591         if (inputBeamWidth_ > 0 ) myalgo->SetInputBeamWidth( inputBeamWidth_ );
00592 
00593 
00594    myalgo->SetFitVariable(std::string("d*z"));
00595    myalgo->SetFitType(std::string("likelihood"));
00596    fbeamWidthFit = myalgo->Fit();
00597 
00598    //Add to .txt file
00599    if(writeTxt_   ) dumpBWTxtFile(outputTxt_);
00600 
00601    delete myalgo;
00602 
00603    // not fake
00604    if ( fbeamspot.type() != 0 )
00605        widthfit_ok = true;
00606   }
00607   else{
00608     fbeamspot.setType(reco::BeamSpot::Fake);
00609     edm::LogWarning("BeamFitter") << "Not enough good tracks selected! No beam fit!" << std::endl;
00610   }
00611   return widthfit_ok;
00612 }
00613 
00614 void BeamFitter::dumpBWTxtFile(std::string & fileName ){
00615     std::ofstream outFile;
00616     outFile.open(fileName.c_str(),std::ios::app);
00617     outFile<<"-------------------------------------------------------------------------------------------------------------------------------------------------------------"<<std::endl;
00618     outFile<<"Beam width(in cm) from Log-likelihood fit (Here we assume a symmetric beam(SigmaX=SigmaY)!)"<<std::endl;
00619     outFile<<"   "<<std::endl;
00620     outFile << "BeamWidth =  " <<fbeamWidthFit.BeamWidthX() <<" +/- "<<fbeamWidthFit.BeamWidthXError() << std::endl;
00621     outFile.close();
00622 }
00623 
00624 void BeamFitter::dumpTxtFile(std::string & fileName, bool append){
00625   std::ofstream outFile;
00626 
00627   std::string tmpname = outputTxt_;
00628   char index[15];
00629   if (appendRunTxt_ && !ffilename_changed ) {
00630       sprintf(index,"%s%i","_Run", frun );
00631       tmpname.insert(outputTxt_.length()-4,index);
00632       fileName = tmpname;
00633       ffilename_changed = true;
00634   }
00635 
00636   if (!append)
00637     outFile.open(fileName.c_str());
00638   else
00639     outFile.open(fileName.c_str(),std::ios::app);
00640 
00641   if ( MyPVFitter->IsFitPerBunchCrossing() ) {
00642 
00643     for (std::map<int,reco::BeamSpot>::const_iterator abspot = fbspotPVMap.begin(); abspot!= fbspotPVMap.end(); ++abspot) {
00644       reco::BeamSpot beamspottmp = abspot->second;
00645       int bx = abspot->first;
00646 
00647       outFile << "Runnumber " << frun << " bx " << bx << std::endl;
00648       outFile << "BeginTimeOfFit " << fbeginTimeOfFit << " " << freftime[0] << std::endl;
00649       outFile << "EndTimeOfFit " << fendTimeOfFit << " " << freftime[1] << std::endl;
00650       outFile << "LumiRange " << fbeginLumiOfFit << " - " << fendLumiOfFit << std::endl;
00651       outFile << "Type " << beamspottmp.type() << std::endl;
00652       outFile << "X0 " << beamspottmp.x0() << std::endl;
00653       outFile << "Y0 " << beamspottmp.y0() << std::endl;
00654       outFile << "Z0 " << beamspottmp.z0() << std::endl;
00655       outFile << "sigmaZ0 " << beamspottmp.sigmaZ() << std::endl;
00656       outFile << "dxdz " << beamspottmp.dxdz() << std::endl;
00657       outFile << "dydz " << beamspottmp.dydz() << std::endl;
00658       outFile << "BeamWidthX " << beamspottmp.BeamWidthX() << std::endl;
00659       outFile << "BeamWidthY " << beamspottmp.BeamWidthY() << std::endl;
00660       for (int i = 0; i<6; ++i) {
00661         outFile << "Cov("<<i<<",j) ";
00662         for (int j=0; j<7; ++j) {
00663           outFile << beamspottmp.covariance(i,j) << " ";
00664         }
00665         outFile << std::endl;
00666       }
00667       outFile << "Cov(6,j) 0 0 0 0 0 0 " << beamspottmp.covariance(6,6) << std::endl;
00668       //}
00669       outFile << "EmittanceX " << beamspottmp.emittanceX() << std::endl;
00670       outFile << "EmittanceY " << beamspottmp.emittanceY() << std::endl;
00671       outFile << "BetaStar " << beamspottmp.betaStar() << std::endl;
00672 
00673     }
00674   }
00675   else {
00676     outFile << "Runnumber " << frun << std::endl;
00677     outFile << "BeginTimeOfFit " << fbeginTimeOfFit << " " << freftime[0] << std::endl;
00678     outFile << "EndTimeOfFit " << fendTimeOfFit << " " << freftime[1] << std::endl;
00679     outFile << "LumiRange " << fbeginLumiOfFit << " - " << fendLumiOfFit << std::endl;
00680     outFile << "Type " << fbeamspot.type() << std::endl;
00681     outFile << "X0 " << fbeamspot.x0() << std::endl;
00682     outFile << "Y0 " << fbeamspot.y0() << std::endl;
00683     outFile << "Z0 " << fbeamspot.z0() << std::endl;
00684     outFile << "sigmaZ0 " << fbeamspot.sigmaZ() << std::endl;
00685     outFile << "dxdz " << fbeamspot.dxdz() << std::endl;
00686     outFile << "dydz " << fbeamspot.dydz() << std::endl;
00687     //  if (inputBeamWidth_ > 0 ) {
00688     //    outFile << "BeamWidthX " << inputBeamWidth_ << std::endl;
00689     //    outFile << "BeamWidthY " << inputBeamWidth_ << std::endl;
00690     //  } else {
00691     outFile << "BeamWidthX " << fbeamspot.BeamWidthX() << std::endl;
00692     outFile << "BeamWidthY " << fbeamspot.BeamWidthY() << std::endl;
00693     //  }
00694 
00695     for (int i = 0; i<6; ++i) {
00696       outFile << "Cov("<<i<<",j) ";
00697       for (int j=0; j<7; ++j) {
00698         outFile << fbeamspot.covariance(i,j) << " ";
00699       }
00700       outFile << std::endl;
00701     }
00702 
00703     // beam width error
00704     //if (inputBeamWidth_ > 0 ) {
00705     //  outFile << "Cov(6,j) 0 0 0 0 0 0 " << "1e-4" << std::endl;
00706     //} else {
00707     outFile << "Cov(6,j) 0 0 0 0 0 0 " << fbeamspot.covariance(6,6) << std::endl;
00708     //}
00709     outFile << "EmittanceX " << fbeamspot.emittanceX() << std::endl;
00710     outFile << "EmittanceY " << fbeamspot.emittanceY() << std::endl;
00711     outFile << "BetaStar " << fbeamspot.betaStar() << std::endl;
00712   }
00713   outFile.close();
00714 }
00715 
00716 void BeamFitter::write2DB(){
00717   BeamSpotObjects *pBSObjects = new BeamSpotObjects();
00718 
00719   pBSObjects->SetPosition(fbeamspot.position().X(),fbeamspot.position().Y(),fbeamspot.position().Z());
00720   //std::cout << " wrote: x= " << fbeamspot.position().X() << " y= "<< fbeamspot.position().Y() << " z= " << fbeamspot.position().Z() << std::endl;
00721   pBSObjects->SetSigmaZ(fbeamspot.sigmaZ());
00722   pBSObjects->Setdxdz(fbeamspot.dxdz());
00723   pBSObjects->Setdydz(fbeamspot.dydz());
00724   //if (inputBeamWidth_ > 0 ) {
00725   //  std::cout << " beam width value forced to be " << inputBeamWidth_ << std::endl;
00726   //  pBSObjects->SetBeamWidthX(inputBeamWidth_);
00727   //  pBSObjects->SetBeamWidthY(inputBeamWidth_);
00728   //} else {
00729     // need to fix this
00730     //std::cout << " using default value, 15e-4, for beam width!!!"<<std::endl;
00731   pBSObjects->SetBeamWidthX(fbeamspot.BeamWidthX() );
00732   pBSObjects->SetBeamWidthY(fbeamspot.BeamWidthY() );
00733   //}
00734 
00735   for (int i = 0; i<7; ++i) {
00736     for (int j=0; j<7; ++j) {
00737       pBSObjects->SetCovariance(i,j,fbeamspot.covariance(i,j));
00738     }
00739   }
00740   edm::Service<cond::service::PoolDBOutputService> poolDbService;
00741   if( poolDbService.isAvailable() ) {
00742     std::cout << "poolDBService available"<<std::endl;
00743     if ( poolDbService->isNewTagRequest( "BeamSpotObjectsRcd" ) ) {
00744       std::cout << "new tag requested" << std::endl;
00745       poolDbService->createNewIOV<BeamSpotObjects>( pBSObjects, poolDbService->beginOfTime(),poolDbService->endOfTime(),
00746                                                     "BeamSpotObjectsRcd"  );
00747     }
00748     else {
00749       std::cout << "no new tag requested" << std::endl;
00750       poolDbService->appendSinceTime<BeamSpotObjects>( pBSObjects, poolDbService->currentTime(),
00751                                                        "BeamSpotObjectsRcd" );
00752     }
00753   }
00754 }
00755 
00756 void BeamFitter::runAllFitter() {
00757   if(fBSvector.size()!=0){
00758     BSFitter *myalgo = new BSFitter( fBSvector );
00759     fbeamspot = myalgo->Fit_d0phi();
00760 
00761     // iterative
00762     if(debug_)
00763       std::cout << " d0-phi Iterative:" << std::endl;
00764     BSFitter *myitealgo = new BSFitter( fBSvector );
00765     myitealgo->Setd0Cut_d0phi(4.0);
00766     reco::BeamSpot beam_ite = myitealgo->Fit_ited0phi();
00767     if (debug_){
00768       std::cout << beam_ite << std::endl;
00769       std::cout << "\n Now run tests of the different fits\n";
00770     }
00771     // from here are just tests
00772     std::string fit_type = "chi2";
00773     myalgo->SetFitVariable(std::string("z"));
00774     myalgo->SetFitType(std::string("chi2"));
00775     reco::BeamSpot beam_fit_z_chi2 = myalgo->Fit();
00776     if (debug_){
00777       std::cout << " z Chi2 Fit ONLY:" << std::endl;
00778       std::cout << beam_fit_z_chi2 << std::endl;
00779     }
00780 
00781     fit_type = "combined";
00782     myalgo->SetFitVariable("z");
00783     myalgo->SetFitType("combined");
00784     reco::BeamSpot beam_fit_z_lh = myalgo->Fit();
00785     if (debug_){
00786       std::cout << " z Log-Likelihood Fit ONLY:" << std::endl;
00787       std::cout << beam_fit_z_lh << std::endl;
00788     }
00789 
00790     myalgo->SetFitVariable("d");
00791     myalgo->SetFitType("d0phi");
00792     reco::BeamSpot beam_fit_dphi = myalgo->Fit();
00793     if (debug_){
00794       std::cout << " d0-phi0 Fit ONLY:" << std::endl;
00795       std::cout << beam_fit_dphi << std::endl;
00796     }
00797     /*
00798     myalgo->SetFitVariable(std::string("d*z"));
00799     myalgo->SetFitType(std::string("likelihood"));
00800     reco::BeamSpot beam_fit_dz_lh = myalgo->Fit();
00801     if (debug_){
00802       std::cout << " Log-Likelihood Fit:" << std::endl;
00803       std::cout << beam_fit_dz_lh << std::endl;
00804     }
00805 
00806     myalgo->SetFitVariable(std::string("d*z"));
00807     myalgo->SetFitType(std::string("resolution"));
00808     reco::BeamSpot beam_fit_dresz_lh = myalgo->Fit();
00809     if(debug_){
00810       std::cout << " IP Resolution Fit" << std::endl;
00811       std::cout << beam_fit_dresz_lh << std::endl;
00812 
00813       std::cout << "c0 = " << myalgo->GetResPar0() << " +- " << myalgo->GetResPar0Err() << std::endl;
00814       std::cout << "c1 = " << myalgo->GetResPar1() << " +- " << myalgo->GetResPar1Err() << std::endl;
00815     }
00816     */
00817   }
00818   else
00819     if (debug_) std::cout << "No good track selected! No beam fit!" << std::endl;
00820 }