CMS 3D CMS Logo

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