CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/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   if ( hasPVs ) {
00120 
00121       for (reco::VertexCollection::const_iterator pv = PVCollection->begin(); pv != PVCollection->end(); ++pv ) {
00122 
00123 
00124            //for ( size_t ipv=0; ipv != pv.size(); ++ipv ) {
00125 
00126           //--- vertex selection
00127           if ( pv->isFake() || pv->tracksSize()==0 )  continue;
00128           if ( pv->ndof() < minVtxNdf_ || (pv->ndof()+3.)/pv->tracksSize()<2*minVtxWgt_ )  continue;
00129           //---
00130 
00131           hPVx->Fill( pv->x(), pv->z() );
00132           hPVy->Fill( pv->y(), pv->z() );
00133 
00134           //
00135           // 3D fit section
00136           //
00137           // apply additional quality cut
00138           if ( pvQuality(*pv)>dynamicQualityCut_ )  continue;
00139           // if store exceeds max. size: reduce size and apply new quality cut
00140           if ( pvStore_.size()>=maxNrVertices_ ) {
00141              compressStore();
00142              if ( pvQuality(*pv)>dynamicQualityCut_ )  continue;
00143           }
00144           //
00145           // copy PV to store
00146           //
00147           int bx = iEvent.bunchCrossing();
00148           BeamSpotFitPVData pvData;
00149           pvData.bunchCrossing = bx;
00150           pvData.position[0] = pv->x();
00151           pvData.position[1] = pv->y();
00152           pvData.position[2] = pv->z();
00153           pvData.posError[0] = pv->xError();
00154           pvData.posError[1] = pv->yError();
00155           pvData.posError[2] = pv->zError();
00156           pvData.posCorr[0] = pv->covariance(0,1)/pv->xError()/pv->yError();
00157           pvData.posCorr[1] = pv->covariance(0,2)/pv->xError()/pv->zError();
00158           pvData.posCorr[2] = pv->covariance(1,2)/pv->yError()/pv->zError();
00159           pvStore_.push_back(pvData);
00160 
00161           if(ftree_ != 0){
00162             theBeamSpotTreeData_.run(iEvent.id().run());
00163             theBeamSpotTreeData_.lumi(iEvent.luminosityBlock());
00164             theBeamSpotTreeData_.bunchCrossing(bx);
00165             theBeamSpotTreeData_.pvData(pvData);
00166             ftree_->Fill();
00167           }
00168 
00169           if (fFitPerBunchCrossing) bxMap_[bx].push_back(pvData);
00170 
00171       }
00172 
00173   }
00174 
00175 
00176 
00177 
00178 }
00179 
00180 void PVFitter::setTree(TTree* tree){
00181   ftree_ = tree;
00182   theBeamSpotTreeData_.branch(ftree_);
00183 }
00184 
00185 bool PVFitter::runBXFitter() {
00186 
00187   edm::LogInfo("PVFitter") << " Number of bunch crossings: " << bxMap_.size() << std::endl;
00188 
00189   bool fit_ok = true;
00190 
00191   for ( std::map<int,std::vector<BeamSpotFitPVData> >::const_iterator pvStore = bxMap_.begin();
00192         pvStore!=bxMap_.end(); ++pvStore) {
00193 
00194     // first set null beam spot in case
00195     // fit fails
00196     fbspotMap[pvStore->first] = reco::BeamSpot();
00197 
00198     edm::LogInfo("PVFitter") << " Number of PVs collected for PVFitter: " << (pvStore->second).size() << " in bx: " << pvStore->first << std::endl;
00199 
00200     if ( (pvStore->second).size() <= minNrVertices_ ) {
00201         edm::LogWarning("PVFitter") << " not enough PVs, continue" << std::endl;
00202         fit_ok = false;
00203       continue;
00204     }
00205 
00206     //bool fit_ok = false;
00207     edm::LogInfo("PVFitter") << "Calculating beam spot with PVs ..." << std::endl;
00208 
00209     //
00210     // LL function and fitter
00211     //
00212     FcnBeamSpotFitPV* fcn = new FcnBeamSpotFitPV(pvStore->second);
00213     TFitterMinuit minuitx;
00214     minuitx.SetMinuitFCN(fcn);
00215     //
00216     // fit parameters: positions, widths, x-y correlations, tilts in xz and yz
00217     //
00218     minuitx.SetParameter(0,"x",0.,0.02,-10.,10.);
00219     minuitx.SetParameter(1,"y",0.,0.02,-10.,10.);
00220     minuitx.SetParameter(2,"z",0.,0.20,-30.,30.);
00221     minuitx.SetParameter(3,"ex",0.015,0.01,0.,10.);
00222     minuitx.SetParameter(4,"corrxy",0.,0.02,-1.,1.);
00223     minuitx.SetParameter(5,"ey",0.015,0.01,0.,10.);
00224     minuitx.SetParameter(6,"dxdz",0.,0.0002,-0.1,0.1);
00225     minuitx.SetParameter(7,"dydz",0.,0.0002,-0.1,0.1);
00226     minuitx.SetParameter(8,"ez",1.,0.1,0.,30.);
00227     minuitx.SetParameter(9,"scale",errorScale_,errorScale_/10.,errorScale_/2.,errorScale_*2.);
00228     //
00229     // first iteration without correlations
00230     //
00231     int ierr(0);
00232     minuitx.FixParameter(4);
00233     minuitx.FixParameter(6);
00234     minuitx.FixParameter(7);
00235     minuitx.FixParameter(9);
00236     minuitx.SetMaxIterations(100);
00237     //       minuitx.SetPrintLevel(3);
00238     minuitx.SetPrintLevel(0);
00239     minuitx.CreateMinimizer();
00240     ierr = minuitx.Minimize();
00241     if ( ierr ) {
00242         edm::LogInfo("PVFitter") << "3D beam spot fit failed in 1st iteration" << std::endl;
00243         fit_ok = false;
00244       continue;
00245     }
00246     //
00247     // refit with harder selection on vertices
00248     //
00249     fcn->setLimits(minuitx.GetParameter(0)-sigmaCut_*minuitx.GetParameter(3),
00250                    minuitx.GetParameter(0)+sigmaCut_*minuitx.GetParameter(3),
00251                    minuitx.GetParameter(1)-sigmaCut_*minuitx.GetParameter(5),
00252                    minuitx.GetParameter(1)+sigmaCut_*minuitx.GetParameter(5),
00253                    minuitx.GetParameter(2)-sigmaCut_*minuitx.GetParameter(8),
00254                    minuitx.GetParameter(2)+sigmaCut_*minuitx.GetParameter(8));
00255     ierr = minuitx.Minimize();
00256     if ( ierr ) {
00257       edm::LogInfo("PVFitter") << "3D beam spot fit failed in 2nd iteration" << std::endl;
00258       fit_ok = false;
00259       continue;
00260     }
00261     //
00262     // refit with correlations
00263     //
00264     minuitx.ReleaseParameter(4);
00265     minuitx.ReleaseParameter(6);
00266     minuitx.ReleaseParameter(7);
00267 
00268     ierr = minuitx.Minimize();
00269     if ( ierr ) {
00270         edm::LogInfo("PVFitter") << "3D beam spot fit failed in 3rd iteration" << std::endl;
00271         fit_ok = false;
00272       continue;
00273     }
00274     // refit with floating scale factor
00275     //   minuitx.ReleaseParameter(9);
00276     //   minuitx.Minimize();
00277 
00278     //minuitx.PrintResults(0,0);
00279 
00280     fwidthX = minuitx.GetParameter(3);
00281     fwidthY = minuitx.GetParameter(5);
00282     fwidthZ = minuitx.GetParameter(8);
00283     fwidthXerr = minuitx.GetParError(3);
00284     fwidthYerr = minuitx.GetParError(5);
00285     fwidthZerr = minuitx.GetParError(8);
00286 
00287     reco::BeamSpot::CovarianceMatrix matrix;
00288     // need to get the full cov matrix
00289     matrix(0,0) = pow( minuitx.GetParError(0), 2);
00290     matrix(1,1) = pow( minuitx.GetParError(1), 2);
00291     matrix(2,2) = pow( minuitx.GetParError(2), 2);
00292     matrix(3,3) = fwidthZerr * fwidthZerr;
00293     matrix(4,4) = pow( minuitx.GetParError(6), 2);
00294     matrix(5,5) = pow( minuitx.GetParError(7), 2);
00295     matrix(6,6) = fwidthXerr * fwidthXerr;
00296 
00297     fbeamspot = reco::BeamSpot( reco::BeamSpot::Point(minuitx.GetParameter(0),
00298                                                       minuitx.GetParameter(1),
00299                                                       minuitx.GetParameter(2) ),
00300                                 fwidthZ,
00301                                 minuitx.GetParameter(6), minuitx.GetParameter(7),
00302                                 fwidthX,
00303                                 matrix );
00304     fbeamspot.setBeamWidthX( fwidthX );
00305     fbeamspot.setBeamWidthY( fwidthY );
00306     fbeamspot.setType( reco::BeamSpot::Tracker );
00307 
00308     fbspotMap[pvStore->first] = fbeamspot;
00309     edm::LogInfo("PVFitter") << "3D PV fit done for this bunch crossing."<<std::endl;
00310     minuitx.Clear();
00311     //delete fcn;
00312     fit_ok = fit_ok & true;
00313   }
00314 
00315   return fit_ok;
00316 }
00317 
00318 
00319 bool PVFitter::runFitter() {
00320 
00321     edm::LogInfo("PVFitter") << " Number of PVs collected for PVFitter: " << pvStore_.size() << std::endl;
00322 
00323     if ( pvStore_.size() <= minNrVertices_ ) return false;
00324 
00325     //bool fit_ok = false;
00326 
00327     if ( ! do3DFit_ ) {
00328       TH1F *h1PVx = (TH1F*) hPVx->ProjectionX("h1PVx", 0, -1, "e");
00329       TH1F *h1PVy = (TH1F*) hPVy->ProjectionX("h1PVy", 0, -1, "e");
00330       TH1F *h1PVz = (TH1F*) hPVx->ProjectionY("h1PVz", 0, -1, "e");
00331 
00332       h1PVx->Fit("gaus","QLM0");
00333       h1PVy->Fit("gaus","QLM0");
00334       h1PVz->Fit("gaus","QLM0");
00335 
00336       TF1 *gausx = h1PVx->GetFunction("gaus");
00337       TF1 *gausy = h1PVy->GetFunction("gaus");
00338       TF1 *gausz = h1PVz->GetFunction("gaus");
00339 
00340       fwidthX = gausx->GetParameter(2);
00341       fwidthY = gausy->GetParameter(2);
00342       fwidthZ = gausz->GetParameter(2);
00343       fwidthXerr = gausx->GetParError(2);
00344       fwidthYerr = gausy->GetParError(2);
00345       fwidthZerr = gausz->GetParError(2);
00346 
00347       reco::BeamSpot::CovarianceMatrix matrix;
00348       matrix(2,2) = gausz->GetParError(1) * gausz->GetParError(1);
00349       matrix(3,3) = fwidthZerr * fwidthZerr;
00350       matrix(6,6) = fwidthXerr * fwidthXerr;
00351 
00352       fbeamspot = reco::BeamSpot( reco::BeamSpot::Point(gausx->GetParameter(1),
00353                                                         gausy->GetParameter(1),
00354                                                         gausz->GetParameter(1) ),
00355                                   fwidthZ,
00356                                   0., 0.,
00357                                   fwidthX,
00358                                   matrix );
00359       fbeamspot.setBeamWidthX( fwidthX );
00360       fbeamspot.setBeamWidthY( fwidthY );
00361       fbeamspot.setType(reco::BeamSpot::Tracker);
00362 
00363     }
00364     else { // do 3D fit
00365       //
00366       // LL function and fitter
00367       //
00368       FcnBeamSpotFitPV* fcn = new FcnBeamSpotFitPV(pvStore_);
00369       TFitterMinuit minuitx;
00370       minuitx.SetMinuitFCN(fcn);
00371       //
00372       // fit parameters: positions, widths, x-y correlations, tilts in xz and yz
00373       //
00374       minuitx.SetParameter(0,"x",0.,0.02,-10.,10.);
00375       minuitx.SetParameter(1,"y",0.,0.02,-10.,10.);
00376       minuitx.SetParameter(2,"z",0.,0.20,-30.,30.);
00377       //minuitx.SetParameter(3,"ex",0.015,0.01,0.,10.);
00378       minuitx.SetParameter(3,"ex",0.015,0.01,0.0001,10.);
00379       minuitx.SetParameter(4,"corrxy",0.,0.02,-1.,1.);
00380       //minuitx.SetParameter(5,"ey",0.015,0.01,0.,10.);
00381       minuitx.SetParameter(5,"ey",0.015,0.01,0.0001,10.);
00382       minuitx.SetParameter(6,"dxdz",0.,0.0002,-0.1,0.1);
00383       minuitx.SetParameter(7,"dydz",0.,0.0002,-0.1,0.1);
00384       //minuitx.SetParameter(8,"ez",1.,0.1,0.,30.);
00385       minuitx.SetParameter(8,"ez",1.,0.1,1.0,30.);
00386       minuitx.SetParameter(9,"scale",errorScale_,errorScale_/10.,errorScale_/2.,errorScale_*2.);
00387       //
00388       // first iteration without correlations
00389       //
00390       int ierr(0);
00391       minuitx.FixParameter(4);
00392       minuitx.FixParameter(6);
00393       minuitx.FixParameter(7);
00394       minuitx.FixParameter(9);
00395       minuitx.SetMaxIterations(100);
00396 //       minuitx.SetPrintLevel(3);
00397       minuitx.SetPrintLevel(0);
00398       minuitx.CreateMinimizer();
00399       ierr = minuitx.Minimize();
00400       if ( ierr ) {
00401           edm::LogWarning("PVFitter") << "3D beam spot fit failed in 1st iteration" << std::endl;
00402           return false;
00403       }
00404       //
00405       // refit with harder selection on vertices
00406       //
00407       fcn->setLimits(minuitx.GetParameter(0)-sigmaCut_*minuitx.GetParameter(3),
00408                      minuitx.GetParameter(0)+sigmaCut_*minuitx.GetParameter(3),
00409                      minuitx.GetParameter(1)-sigmaCut_*minuitx.GetParameter(5),
00410                      minuitx.GetParameter(1)+sigmaCut_*minuitx.GetParameter(5),
00411                      minuitx.GetParameter(2)-sigmaCut_*minuitx.GetParameter(8),
00412                      minuitx.GetParameter(2)+sigmaCut_*minuitx.GetParameter(8));
00413       ierr = minuitx.Minimize();
00414       if ( ierr ) {
00415           edm::LogWarning("PVFitter") << "3D beam spot fit failed in 2nd iteration" << std::endl;
00416           return false;
00417       }
00418       //
00419       // refit with correlations
00420       //
00421       minuitx.ReleaseParameter(4);
00422       minuitx.ReleaseParameter(6);
00423       minuitx.ReleaseParameter(7);
00424       ierr = minuitx.Minimize();
00425       if ( ierr ) {
00426           edm::LogWarning("PVFitter") << "3D beam spot fit failed in 3rd iteration" << std::endl;
00427           return false;
00428       }
00429       // refit with floating scale factor
00430       //   minuitx.ReleaseParameter(9);
00431       //   minuitx.Minimize();
00432 
00433       //minuitx.PrintResults(0,0);
00434 
00435       fwidthX = minuitx.GetParameter(3);
00436       fwidthY = minuitx.GetParameter(5);
00437       fwidthZ = minuitx.GetParameter(8);
00438       fwidthXerr = minuitx.GetParError(3);
00439       fwidthYerr = minuitx.GetParError(5);
00440       fwidthZerr = minuitx.GetParError(8);
00441 
00442       // check errors on widths and sigmaZ for nan
00443       if ( isnan(fwidthXerr) || isnan(fwidthYerr) || isnan(fwidthZerr) ) {
00444           edm::LogWarning("PVFitter") << "3D beam spot fit returns nan in 3rd iteration" << std::endl;
00445           return false;
00446       }
00447 
00448       reco::BeamSpot::CovarianceMatrix matrix;
00449       // need to get the full cov matrix
00450       matrix(0,0) = pow( minuitx.GetParError(0), 2);
00451       matrix(1,1) = pow( minuitx.GetParError(1), 2);
00452       matrix(2,2) = pow( minuitx.GetParError(2), 2);
00453       matrix(3,3) = fwidthZerr * fwidthZerr;
00454       matrix(6,6) = fwidthXerr * fwidthXerr;
00455 
00456       fbeamspot = reco::BeamSpot( reco::BeamSpot::Point(minuitx.GetParameter(0),
00457                                                         minuitx.GetParameter(1),
00458                                                         minuitx.GetParameter(2) ),
00459                                   fwidthZ,
00460                                   minuitx.GetParameter(6), minuitx.GetParameter(7),
00461                                   fwidthX,
00462                                   matrix );
00463       fbeamspot.setBeamWidthX( fwidthX );
00464       fbeamspot.setBeamWidthY( fwidthY );
00465       fbeamspot.setType(reco::BeamSpot::Tracker);
00466     }
00467 
00468     return true; //FIXME: Need to add quality test for the fit results!
00469 }
00470 
00471 void PVFitter::dumpTxtFile(){
00472 /*
00473   fasciiFile << "Runnumber " << frun << std::endl;
00474   fasciiFile << "BeginTimeOfFit " << fbeginTimeOfFit << std::endl;
00475   fasciiFile << "EndTimeOfFit " << fendTimeOfFit << std::endl;
00476   fasciiFile << "LumiRange " << fbeginLumiOfFit << " - " << fendLumiOfFit << std::endl;
00477   fasciiFile << "Type " << fbeamspot.type() << std::endl;
00478   fasciiFile << "X0 " << fbeamspot.x0() << std::endl;
00479   fasciiFile << "Y0 " << fbeamspot.y0() << std::endl;
00480   fasciiFile << "Z0 " << fbeamspot.z0() << std::endl;
00481   fasciiFile << "sigmaZ0 " << fbeamspot.sigmaZ() << std::endl;
00482   fasciiFile << "dxdz " << fbeamspot.dxdz() << std::endl;
00483   fasciiFile << "dydz " << fbeamspot.dydz() << std::endl;
00484   if (inputBeamWidth_ > 0 ) {
00485     fasciiFile << "BeamWidthX " << inputBeamWidth_ << std::endl;
00486     fasciiFile << "BeamWidthY " << inputBeamWidth_ << std::endl;
00487   } else {
00488     fasciiFile << "BeamWidthX " << fbeamspot.BeamWidthX() << std::endl;
00489     fasciiFile << "BeamWidthY " << fbeamspot.BeamWidthY() << std::endl;
00490   }
00491 
00492   for (int i = 0; i<6; ++i) {
00493     fasciiFile << "Cov("<<i<<",j) ";
00494     for (int j=0; j<7; ++j) {
00495       fasciiFile << fbeamspot.covariance(i,j) << " ";
00496     }
00497     fasciiFile << std::endl;
00498   }
00499   // beam width error
00500   if (inputBeamWidth_ > 0 ) {
00501     fasciiFile << "Cov(6,j) 0 0 0 0 0 0 " << "1e-4" << std::endl;
00502   } else {
00503     fasciiFile << "Cov(6,j) 0 0 0 0 0 0 " << fbeamspot.covariance(6,6) << std::endl;
00504   }
00505   fasciiFile << "EmittanceX " << fbeamspot.emittanceX() << std::endl;
00506   fasciiFile << "EmittanceY " << fbeamspot.emittanceY() << std::endl;
00507   fasciiFile << "BetaStar " << fbeamspot.betaStar() << std::endl;
00508 
00509 */
00510 }
00511 
00512 
00513 void
00514 PVFitter::compressStore ()
00515 {
00516   //
00517   // fill vertex qualities
00518   //
00519   pvQualities_.resize(pvStore_.size());
00520   for ( unsigned int i=0; i<pvStore_.size(); ++i )  pvQualities_[i] = pvQuality(pvStore_[i]);
00521   sort(pvQualities_.begin(),pvQualities_.end());
00522   //
00523   // Set new quality cut to median. This cut will be used to reduce the
00524   // number of vertices in the store and also apply to all new vertices
00525   // until the next reset
00526   //
00527   dynamicQualityCut_ = pvQualities_[pvQualities_.size()/2];
00528   //
00529   // remove all vertices failing the cut from the store
00530   //   (to be moved to a more efficient memory management!)
00531   //
00532   unsigned int iwrite(0);
00533   for ( unsigned int i=0; i<pvStore_.size(); ++i ) {
00534     if ( pvQuality(pvStore_[i])>dynamicQualityCut_ )  continue;
00535     if ( i!=iwrite )  pvStore_[iwrite] = pvStore_[i];
00536     ++iwrite;
00537   }
00538   pvStore_.resize(iwrite);
00539   edm::LogInfo("PVFitter") << "Reduced primary vertex store size to "
00540                            << pvStore_.size() << " ; new dynamic quality cut = "
00541                            << dynamicQualityCut_ << std::endl;
00542 
00543 }
00544 
00545 double
00546 PVFitter::pvQuality (const reco::Vertex& pv) const
00547 {
00548   //
00549   // determinant of the transverse part of the PV covariance matrix
00550   //
00551   return
00552     pv.covariance(0,0)*pv.covariance(1,1)-
00553     pv.covariance(0,1)*pv.covariance(0,1);
00554 }
00555 
00556 double
00557 PVFitter::pvQuality (const BeamSpotFitPVData& pv) const
00558 {
00559   //
00560   // determinant of the transverse part of the PV covariance matrix
00561   //
00562   double ex = pv.posError[0];
00563   double ey = pv.posError[1];
00564   return ex*ex*ey*ey*(1-pv.posCorr[0]*pv.posCorr[0]);
00565 }
00566