CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/RecoVertex/BeamSpotProducer/src/PVFitter.cc

Go to the documentation of this file.
00001 
00014 #include "RecoVertex/BeamSpotProducer/interface/PVFitter.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 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00029 
00030 #include "TFitterMinuit.h"
00031 #include "Minuit2/FCNBase.h"
00032 #include "RecoVertex/BeamSpotProducer/interface/FcnBeamSpotFitPV.h"
00033 #include "FWCore/Utilities/interface/isFinite.h"
00034 
00035 #include "TF1.h"
00036 
00037 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00038 
00039 // ----------------------------------------------------------------------
00040 // Useful function:
00041 // ----------------------------------------------------------------------
00042 
00043 // static char * formatTime(const std::time_t & t)  {
00044 //   struct std::tm * ptm;
00045 //   ptm = gmtime(&t);
00046 //   static char ts[32];
00047 //   strftime(ts,sizeof(ts),"%Y.%m.%d %H:%M:%S %Z",ptm);
00048 //   return ts;
00049 // }
00050 
00051 PVFitter::PVFitter(const edm::ParameterSet& iConfig): ftree_(0)
00052 {
00053 
00054   debug_             = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("Debug");
00055   vertexLabel_     = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<edm::InputTag>("VertexCollection", edm::InputTag("offlinePrimaryVertices"));
00056   do3DFit_           = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("Apply3DFit");
00057   //writeTxt_          = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("WriteAscii");
00058   //outputTxt_         = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<std::string>("AsciiFileName");
00059 
00060   maxNrVertices_     = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<unsigned int>("maxNrStoredVertices");
00061   minNrVertices_     = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<unsigned int>("minNrVerticesForFit");
00062   minVtxNdf_         = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("minVertexNdf");
00063   maxVtxNormChi2_    = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("maxVertexNormChi2");
00064   minVtxTracks_      = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<unsigned int>("minVertexNTracks");
00065   minVtxWgt_         = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("minVertexMeanWeight");
00066   maxVtxR_           = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("maxVertexR");
00067   maxVtxZ_           = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("maxVertexZ");
00068   errorScale_        = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("errorScale");
00069   sigmaCut_          = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("nSigmaCut");
00070   fFitPerBunchCrossing=iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("FitPerBunchCrossing");
00071 
00072   // preset quality cut to "infinite"
00073   dynamicQualityCut_ = 1.e30;
00074 
00075   hPVx = new TH2F("hPVx","PVx vs PVz distribution",200,-maxVtxR_, maxVtxR_, 200, -maxVtxZ_, maxVtxZ_);
00076   hPVy = new TH2F("hPVy","PVy vs PVz distribution",200,-maxVtxR_, maxVtxR_, 200, -maxVtxZ_, maxVtxZ_);
00077 }
00078 
00079 PVFitter::~PVFitter() {
00080 
00081 }
00082 
00083 
00084 void PVFitter::readEvent(const edm::Event& iEvent)
00085 {
00086 
00087 //   frun = iEvent.id().run();
00088 //   const edm::TimeValue_t ftimestamp = iEvent.time().value();
00089 //   const std::time_t ftmptime = ftimestamp >> 32;
00090 
00091 //   if (fbeginLumiOfFit == -1) freftime[0] = freftime[1] = ftmptime;
00092 //   if (freftime[0] == 0 || ftmptime < freftime[0]) freftime[0] = ftmptime;
00093 //   const char* fbeginTime = formatTime(freftime[0]);
00094 //   sprintf(fbeginTimeOfFit,"%s",fbeginTime);
00095 
00096 //   if (freftime[1] == 0 || ftmptime > freftime[1]) freftime[1] = ftmptime;
00097 //   const char* fendTime = formatTime(freftime[1]);
00098 //   sprintf(fendTimeOfFit,"%s",fendTime);
00099 
00100 //   flumi = iEvent.luminosityBlock();
00101 //   frunFit = frun;
00102 
00103 //   if (fbeginLumiOfFit == -1 || fbeginLumiOfFit > flumi) fbeginLumiOfFit = flumi;
00104 //   if (fendLumiOfFit == -1 || fendLumiOfFit < flumi) fendLumiOfFit = flumi;
00105 //   std::cout << "flumi = " <<flumi<<"; fbeginLumiOfFit = " << fbeginLumiOfFit <<"; fendLumiOfFit = "<<fendLumiOfFit<<std::endl;
00106 
00107   //------ Primary Vertices
00108   edm::Handle< reco::VertexCollection > PVCollection;
00109   bool hasPVs = false;
00110   //edm::View<reco::Vertex> vertices;
00111   //const reco::VertexCollection & vertices = 0;
00112 
00113   if ( iEvent.getByLabel(vertexLabel_, PVCollection ) ) {
00114       //pv = *PVCollection;
00115       //vertices = *PVCollection;
00116       hasPVs = true;
00117   }
00118   //------
00119 
00120   if ( hasPVs ) {
00121 
00122       for (reco::VertexCollection::const_iterator pv = PVCollection->begin(); pv != PVCollection->end(); ++pv ) {
00123 
00124 
00125            //for ( size_t ipv=0; ipv != pv.size(); ++ipv ) {
00126 
00127           //--- vertex selection
00128           if ( pv->isFake() || pv->tracksSize()==0 )  continue;
00129           if ( pv->ndof() < minVtxNdf_ || (pv->ndof()+3.)/pv->tracksSize()<2*minVtxWgt_ )  continue;
00130           //---
00131 
00132           hPVx->Fill( pv->x(), pv->z() );
00133           hPVy->Fill( pv->y(), pv->z() );
00134 
00135           //
00136           // 3D fit section
00137           //
00138           // apply additional quality cut
00139           if ( pvQuality(*pv)>dynamicQualityCut_ )  continue;
00140           // if store exceeds max. size: reduce size and apply new quality cut
00141           if ( pvStore_.size()>=maxNrVertices_ ) {
00142              compressStore();
00143              if ( pvQuality(*pv)>dynamicQualityCut_ )  continue;
00144           }
00145           //
00146           // copy PV to store
00147           //
00148           int bx = iEvent.bunchCrossing();
00149           BeamSpotFitPVData pvData;
00150           pvData.bunchCrossing = bx;
00151           pvData.position[0] = pv->x();
00152           pvData.position[1] = pv->y();
00153           pvData.position[2] = pv->z();
00154           pvData.posError[0] = pv->xError();
00155           pvData.posError[1] = pv->yError();
00156           pvData.posError[2] = pv->zError();
00157           pvData.posCorr[0] = pv->covariance(0,1)/pv->xError()/pv->yError();
00158           pvData.posCorr[1] = pv->covariance(0,2)/pv->xError()/pv->zError();
00159           pvData.posCorr[2] = pv->covariance(1,2)/pv->yError()/pv->zError();
00160           pvStore_.push_back(pvData);
00161 
00162           if(ftree_ != 0){
00163             theBeamSpotTreeData_.run(iEvent.id().run());
00164             theBeamSpotTreeData_.lumi(iEvent.luminosityBlock());
00165             theBeamSpotTreeData_.bunchCrossing(bx);
00166             theBeamSpotTreeData_.pvData(pvData);
00167             ftree_->Fill();
00168           }
00169 
00170           if (fFitPerBunchCrossing) bxMap_[bx].push_back(pvData);
00171 
00172       }
00173 
00174   }
00175 
00176 
00177 
00178 
00179 }
00180 
00181 void PVFitter::setTree(TTree* tree){
00182   ftree_ = tree;
00183   theBeamSpotTreeData_.branch(ftree_);
00184 }
00185 
00186 bool PVFitter::runBXFitter() {
00187 
00188   edm::LogInfo("PVFitter") << " Number of bunch crossings: " << bxMap_.size() << std::endl;
00189 
00190   bool fit_ok = true;
00191 
00192   for ( std::map<int,std::vector<BeamSpotFitPVData> >::const_iterator pvStore = bxMap_.begin();
00193         pvStore!=bxMap_.end(); ++pvStore) {
00194 
00195     // first set null beam spot in case
00196     // fit fails
00197     fbspotMap[pvStore->first] = reco::BeamSpot();
00198 
00199     edm::LogInfo("PVFitter") << " Number of PVs collected for PVFitter: " << (pvStore->second).size() << " in bx: " << pvStore->first << std::endl;
00200 
00201     if ( (pvStore->second).size() <= minNrVertices_ ) {
00202         edm::LogWarning("PVFitter") << " not enough PVs, continue" << std::endl;
00203         fit_ok = false;
00204       continue;
00205     }
00206 
00207     //bool fit_ok = false;
00208     edm::LogInfo("PVFitter") << "Calculating beam spot with PVs ..." << std::endl;
00209 
00210     //
00211     // LL function and fitter
00212     //
00213     FcnBeamSpotFitPV* fcn = new FcnBeamSpotFitPV(pvStore->second);
00214     TFitterMinuit minuitx;
00215     minuitx.SetMinuitFCN(fcn);
00216     //
00217     // fit parameters: positions, widths, x-y correlations, tilts in xz and yz
00218     //
00219     minuitx.SetParameter(0,"x",0.,0.02,-10.,10.);
00220     minuitx.SetParameter(1,"y",0.,0.02,-10.,10.);
00221     minuitx.SetParameter(2,"z",0.,0.20,-30.,30.);
00222     minuitx.SetParameter(3,"ex",0.015,0.01,0.,10.);
00223     minuitx.SetParameter(4,"corrxy",0.,0.02,-1.,1.);
00224     minuitx.SetParameter(5,"ey",0.015,0.01,0.,10.);
00225     minuitx.SetParameter(6,"dxdz",0.,0.0002,-0.1,0.1);
00226     minuitx.SetParameter(7,"dydz",0.,0.0002,-0.1,0.1);
00227     minuitx.SetParameter(8,"ez",1.,0.1,0.,30.);
00228     minuitx.SetParameter(9,"scale",errorScale_,errorScale_/10.,errorScale_/2.,errorScale_*2.);
00229     //
00230     // first iteration without correlations
00231     //
00232     int ierr(0);
00233     minuitx.FixParameter(4);
00234     minuitx.FixParameter(6);
00235     minuitx.FixParameter(7);
00236     minuitx.FixParameter(9);
00237     minuitx.SetMaxIterations(100);
00238     //       minuitx.SetPrintLevel(3);
00239     minuitx.SetPrintLevel(0);
00240     minuitx.CreateMinimizer();
00241     ierr = minuitx.Minimize();
00242     if ( ierr ) {
00243         edm::LogInfo("PVFitter") << "3D beam spot fit failed in 1st iteration" << std::endl;
00244         fit_ok = false;
00245       continue;
00246     }
00247     //
00248     // refit with harder selection on vertices
00249     //
00250     fcn->setLimits(minuitx.GetParameter(0)-sigmaCut_*minuitx.GetParameter(3),
00251                    minuitx.GetParameter(0)+sigmaCut_*minuitx.GetParameter(3),
00252                    minuitx.GetParameter(1)-sigmaCut_*minuitx.GetParameter(5),
00253                    minuitx.GetParameter(1)+sigmaCut_*minuitx.GetParameter(5),
00254                    minuitx.GetParameter(2)-sigmaCut_*minuitx.GetParameter(8),
00255                    minuitx.GetParameter(2)+sigmaCut_*minuitx.GetParameter(8));
00256     ierr = minuitx.Minimize();
00257     if ( ierr ) {
00258       edm::LogInfo("PVFitter") << "3D beam spot fit failed in 2nd iteration" << std::endl;
00259       fit_ok = false;
00260       continue;
00261     }
00262     //
00263     // refit with correlations
00264     //
00265     minuitx.ReleaseParameter(4);
00266     minuitx.ReleaseParameter(6);
00267     minuitx.ReleaseParameter(7);
00268 
00269     ierr = minuitx.Minimize();
00270     if ( ierr ) {
00271         edm::LogInfo("PVFitter") << "3D beam spot fit failed in 3rd iteration" << std::endl;
00272         fit_ok = false;
00273       continue;
00274     }
00275     // refit with floating scale factor
00276     //   minuitx.ReleaseParameter(9);
00277     //   minuitx.Minimize();
00278 
00279     //minuitx.PrintResults(0,0);
00280 
00281     fwidthX = minuitx.GetParameter(3);
00282     fwidthY = minuitx.GetParameter(5);
00283     fwidthZ = minuitx.GetParameter(8);
00284     fwidthXerr = minuitx.GetParError(3);
00285     fwidthYerr = minuitx.GetParError(5);
00286     fwidthZerr = minuitx.GetParError(8);
00287 
00288     reco::BeamSpot::CovarianceMatrix matrix;
00289     // need to get the full cov matrix
00290     matrix(0,0) = pow( minuitx.GetParError(0), 2);
00291     matrix(1,1) = pow( minuitx.GetParError(1), 2);
00292     matrix(2,2) = pow( minuitx.GetParError(2), 2);
00293     matrix(3,3) = fwidthZerr * fwidthZerr;
00294     matrix(4,4) = pow( minuitx.GetParError(6), 2);
00295     matrix(5,5) = pow( minuitx.GetParError(7), 2);
00296     matrix(6,6) = fwidthXerr * fwidthXerr;
00297 
00298     fbeamspot = reco::BeamSpot( reco::BeamSpot::Point(minuitx.GetParameter(0),
00299                                                       minuitx.GetParameter(1),
00300                                                       minuitx.GetParameter(2) ),
00301                                 fwidthZ,
00302                                 minuitx.GetParameter(6), minuitx.GetParameter(7),
00303                                 fwidthX,
00304                                 matrix );
00305     fbeamspot.setBeamWidthX( fwidthX );
00306     fbeamspot.setBeamWidthY( fwidthY );
00307     fbeamspot.setType( reco::BeamSpot::Tracker );
00308 
00309     fbspotMap[pvStore->first] = fbeamspot;
00310     edm::LogInfo("PVFitter") << "3D PV fit done for this bunch crossing."<<std::endl;
00311     minuitx.Clear();
00312     //delete fcn;
00313     fit_ok = fit_ok & true;
00314   }
00315 
00316   return fit_ok;
00317 }
00318 
00319 
00320 bool PVFitter::runFitter() {
00321 
00322     edm::LogInfo("PVFitter") << " Number of PVs collected for PVFitter: " << pvStore_.size() << std::endl;
00323 
00324     if ( pvStore_.size() <= minNrVertices_ ) return false;
00325 
00326     //bool fit_ok = false;
00327 
00328     if ( ! do3DFit_ ) {
00329       TH1F *h1PVx = (TH1F*) hPVx->ProjectionX("h1PVx", 0, -1, "e");
00330       TH1F *h1PVy = (TH1F*) hPVy->ProjectionX("h1PVy", 0, -1, "e");
00331       TH1F *h1PVz = (TH1F*) hPVx->ProjectionY("h1PVz", 0, -1, "e");
00332 
00333       h1PVx->Fit("gaus","QLM0");
00334       h1PVy->Fit("gaus","QLM0");
00335       h1PVz->Fit("gaus","QLM0");
00336 
00337       TF1 *gausx = h1PVx->GetFunction("gaus");
00338       TF1 *gausy = h1PVy->GetFunction("gaus");
00339       TF1 *gausz = h1PVz->GetFunction("gaus");
00340 
00341       fwidthX = gausx->GetParameter(2);
00342       fwidthY = gausy->GetParameter(2);
00343       fwidthZ = gausz->GetParameter(2);
00344       fwidthXerr = gausx->GetParError(2);
00345       fwidthYerr = gausy->GetParError(2);
00346       fwidthZerr = gausz->GetParError(2);
00347 
00348       reco::BeamSpot::CovarianceMatrix matrix;
00349       matrix(2,2) = gausz->GetParError(1) * gausz->GetParError(1);
00350       matrix(3,3) = fwidthZerr * fwidthZerr;
00351       matrix(6,6) = fwidthXerr * fwidthXerr;
00352 
00353       fbeamspot = reco::BeamSpot( reco::BeamSpot::Point(gausx->GetParameter(1),
00354                                                         gausy->GetParameter(1),
00355                                                         gausz->GetParameter(1) ),
00356                                   fwidthZ,
00357                                   0., 0.,
00358                                   fwidthX,
00359                                   matrix );
00360       fbeamspot.setBeamWidthX( fwidthX );
00361       fbeamspot.setBeamWidthY( fwidthY );
00362       fbeamspot.setType(reco::BeamSpot::Tracker);
00363 
00364     }
00365     else { // do 3D fit
00366       //
00367       // LL function and fitter
00368       //
00369       FcnBeamSpotFitPV* fcn = new FcnBeamSpotFitPV(pvStore_);
00370       TFitterMinuit minuitx;
00371       minuitx.SetMinuitFCN(fcn);
00372       //
00373       // fit parameters: positions, widths, x-y correlations, tilts in xz and yz
00374       //
00375       minuitx.SetParameter(0,"x",0.,0.02,-10.,10.);
00376       minuitx.SetParameter(1,"y",0.,0.02,-10.,10.);
00377       minuitx.SetParameter(2,"z",0.,0.20,-30.,30.);
00378       //minuitx.SetParameter(3,"ex",0.015,0.01,0.,10.);
00379       //minuitx.SetParameter(3,"ex",0.015,0.01,0.0001,10.);
00380       minuitx.SetParameter(3,"ex",0.005,0.0005,0.0001,0.05);
00381       minuitx.SetParameter(4,"corrxy",0.,0.02,-1.,1.);
00382       //minuitx.SetParameter(5,"ey",0.015,0.01,0.,10.);
00383       //minuitx.SetParameter(5,"ey",0.015,0.01,0.0001,10.);
00384       minuitx.SetParameter(5,"ey",0.005,0.0005,0.0001,0.05);
00385       minuitx.SetParameter(6,"dxdz",0.,0.0002,-0.1,0.1);
00386       minuitx.SetParameter(7,"dydz",0.,0.0002,-0.1,0.1);
00387       //minuitx.SetParameter(8,"ez",1.,0.1,0.,30.);
00388       minuitx.SetParameter(8,"ez",1.,0.1,1.0,30.);
00389       minuitx.SetParameter(9,"scale",errorScale_,errorScale_/10.,errorScale_/2.,errorScale_*2.);
00390       //
00391       // first iteration without correlations
00392       //
00393       int ierr(0);
00394       minuitx.FixParameter(4);
00395       minuitx.FixParameter(6);
00396       minuitx.FixParameter(7);
00397       minuitx.FixParameter(9);
00398       minuitx.SetMaxIterations(100);
00399 //       minuitx.SetPrintLevel(3);
00400       minuitx.SetPrintLevel(0);
00401       minuitx.CreateMinimizer();
00402       ierr = minuitx.Minimize();
00403       if ( ierr ) {
00404           edm::LogWarning("PVFitter") << "3D beam spot fit failed in 1st iteration" << std::endl;
00405           return false;
00406       }
00407       //
00408       // refit with harder selection on vertices
00409       //
00410       fcn->setLimits(minuitx.GetParameter(0)-sigmaCut_*minuitx.GetParameter(3),
00411                      minuitx.GetParameter(0)+sigmaCut_*minuitx.GetParameter(3),
00412                      minuitx.GetParameter(1)-sigmaCut_*minuitx.GetParameter(5),
00413                      minuitx.GetParameter(1)+sigmaCut_*minuitx.GetParameter(5),
00414                      minuitx.GetParameter(2)-sigmaCut_*minuitx.GetParameter(8),
00415                      minuitx.GetParameter(2)+sigmaCut_*minuitx.GetParameter(8));
00416       ierr = minuitx.Minimize();
00417       if ( ierr ) {
00418           edm::LogWarning("PVFitter") << "3D beam spot fit failed in 2nd iteration" << std::endl;
00419           return false;
00420       }
00421       //
00422       // refit with correlations
00423       //
00424       minuitx.ReleaseParameter(4);
00425       minuitx.ReleaseParameter(6);
00426       minuitx.ReleaseParameter(7);
00427       ierr = minuitx.Minimize();
00428       if ( ierr ) {
00429           edm::LogWarning("PVFitter") << "3D beam spot fit failed in 3rd iteration" << std::endl;
00430           return false;
00431       }
00432       // refit with floating scale factor
00433       //   minuitx.ReleaseParameter(9);
00434       //   minuitx.Minimize();
00435 
00436       //minuitx.PrintResults(0,0);
00437 
00438       fwidthX = minuitx.GetParameter(3);
00439       fwidthY = minuitx.GetParameter(5);
00440       fwidthZ = minuitx.GetParameter(8);
00441       fwidthXerr = minuitx.GetParError(3);
00442       fwidthYerr = minuitx.GetParError(5);
00443       fwidthZerr = minuitx.GetParError(8);
00444 
00445       // check errors on widths and sigmaZ for nan
00446       if ( edm::isNotFinite(fwidthXerr) || edm::isNotFinite(fwidthYerr) || edm::isNotFinite(fwidthZerr) ) {
00447           edm::LogWarning("PVFitter") << "3D beam spot fit returns nan in 3rd iteration" << std::endl;
00448           return false;
00449       }
00450 
00451       reco::BeamSpot::CovarianceMatrix matrix;
00452       // need to get the full cov matrix
00453       matrix(0,0) = pow( minuitx.GetParError(0), 2);
00454       matrix(1,1) = pow( minuitx.GetParError(1), 2);
00455       matrix(2,2) = pow( minuitx.GetParError(2), 2);
00456       matrix(3,3) = fwidthZerr * fwidthZerr;
00457       matrix(6,6) = fwidthXerr * fwidthXerr;
00458 
00459       fbeamspot = reco::BeamSpot( reco::BeamSpot::Point(minuitx.GetParameter(0),
00460                                                         minuitx.GetParameter(1),
00461                                                         minuitx.GetParameter(2) ),
00462                                   fwidthZ,
00463                                   minuitx.GetParameter(6), minuitx.GetParameter(7),
00464                                   fwidthX,
00465                                   matrix );
00466       fbeamspot.setBeamWidthX( fwidthX );
00467       fbeamspot.setBeamWidthY( fwidthY );
00468       fbeamspot.setType(reco::BeamSpot::Tracker);
00469     }
00470 
00471     return true; //FIXME: Need to add quality test for the fit results!
00472 }
00473 
00474 void PVFitter::dumpTxtFile(){
00475 /*
00476   fasciiFile << "Runnumber " << frun << std::endl;
00477   fasciiFile << "BeginTimeOfFit " << fbeginTimeOfFit << std::endl;
00478   fasciiFile << "EndTimeOfFit " << fendTimeOfFit << std::endl;
00479   fasciiFile << "LumiRange " << fbeginLumiOfFit << " - " << fendLumiOfFit << std::endl;
00480   fasciiFile << "Type " << fbeamspot.type() << std::endl;
00481   fasciiFile << "X0 " << fbeamspot.x0() << std::endl;
00482   fasciiFile << "Y0 " << fbeamspot.y0() << std::endl;
00483   fasciiFile << "Z0 " << fbeamspot.z0() << std::endl;
00484   fasciiFile << "sigmaZ0 " << fbeamspot.sigmaZ() << std::endl;
00485   fasciiFile << "dxdz " << fbeamspot.dxdz() << std::endl;
00486   fasciiFile << "dydz " << fbeamspot.dydz() << std::endl;
00487   if (inputBeamWidth_ > 0 ) {
00488     fasciiFile << "BeamWidthX " << inputBeamWidth_ << std::endl;
00489     fasciiFile << "BeamWidthY " << inputBeamWidth_ << std::endl;
00490   } else {
00491     fasciiFile << "BeamWidthX " << fbeamspot.BeamWidthX() << std::endl;
00492     fasciiFile << "BeamWidthY " << fbeamspot.BeamWidthY() << std::endl;
00493   }
00494 
00495   for (int i = 0; i<6; ++i) {
00496     fasciiFile << "Cov("<<i<<",j) ";
00497     for (int j=0; j<7; ++j) {
00498       fasciiFile << fbeamspot.covariance(i,j) << " ";
00499     }
00500     fasciiFile << std::endl;
00501   }
00502   // beam width error
00503   if (inputBeamWidth_ > 0 ) {
00504     fasciiFile << "Cov(6,j) 0 0 0 0 0 0 " << "1e-4" << std::endl;
00505   } else {
00506     fasciiFile << "Cov(6,j) 0 0 0 0 0 0 " << fbeamspot.covariance(6,6) << std::endl;
00507   }
00508   fasciiFile << "EmittanceX " << fbeamspot.emittanceX() << std::endl;
00509   fasciiFile << "EmittanceY " << fbeamspot.emittanceY() << std::endl;
00510   fasciiFile << "BetaStar " << fbeamspot.betaStar() << std::endl;
00511 
00512 */
00513 }
00514 
00515 
00516 void
00517 PVFitter::compressStore ()
00518 {
00519   //
00520   // fill vertex qualities
00521   //
00522   pvQualities_.resize(pvStore_.size());
00523   for ( unsigned int i=0; i<pvStore_.size(); ++i )  pvQualities_[i] = pvQuality(pvStore_[i]);
00524   sort(pvQualities_.begin(),pvQualities_.end());
00525   //
00526   // Set new quality cut to median. This cut will be used to reduce the
00527   // number of vertices in the store and also apply to all new vertices
00528   // until the next reset
00529   //
00530   dynamicQualityCut_ = pvQualities_[pvQualities_.size()/2];
00531   //
00532   // remove all vertices failing the cut from the store
00533   //   (to be moved to a more efficient memory management!)
00534   //
00535   unsigned int iwrite(0);
00536   for ( unsigned int i=0; i<pvStore_.size(); ++i ) {
00537     if ( pvQuality(pvStore_[i])>dynamicQualityCut_ )  continue;
00538     if ( i!=iwrite )  pvStore_[iwrite] = pvStore_[i];
00539     ++iwrite;
00540   }
00541   pvStore_.resize(iwrite);
00542   edm::LogInfo("PVFitter") << "Reduced primary vertex store size to "
00543                            << pvStore_.size() << " ; new dynamic quality cut = "
00544                            << dynamicQualityCut_ << std::endl;
00545 
00546 }
00547 
00548 double
00549 PVFitter::pvQuality (const reco::Vertex& pv) const
00550 {
00551   //
00552   // determinant of the transverse part of the PV covariance matrix
00553   //
00554   return
00555     pv.covariance(0,0)*pv.covariance(1,1)-
00556     pv.covariance(0,1)*pv.covariance(0,1);
00557 }
00558 
00559 double
00560 PVFitter::pvQuality (const BeamSpotFitPVData& pv) const
00561 {
00562   //
00563   // determinant of the transverse part of the PV covariance matrix
00564   //
00565   double ex = pv.posError[0];
00566   double ey = pv.posError[1];
00567   return ex*ex*ey*ey*(1-pv.posCorr[0]*pv.posCorr[0]);
00568 }
00569